Ordering.h
Go to the documentation of this file.
1 
2 // This file is part of Eigen, a lightweight C++ template library
3 // for linear algebra.
4 //
5 // Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_ORDERING_H
12 #define EIGEN_ORDERING_H
13 
14 namespace Eigen {
15 
16 #include "Eigen_Colamd.h"
17 
18 namespace internal {
19 
26 template<typename MatrixType>
28 {
29  MatrixType C;
30  C = A.transpose(); // NOTE: Could be costly
31  for (int i = 0; i < C.rows(); i++)
32  {
33  for (typename MatrixType::InnerIterator it(C, i); it; ++it)
34  it.valueRef() = 0.0;
35  }
36  symmat = C + A;
37 }
38 
39 }
40 
41 #ifndef EIGEN_MPL2_ONLY
42 
51 template <typename StorageIndex>
53 {
54  public:
56 
60  template <typename MatrixType>
62  {
63  // Compute the symmetric pattern
66 
67  // Call the AMD routine
68  //m_mat.prune(keep_diag());
70  }
71 
73  template <typename SrcType, unsigned int SrcUpLo>
75  {
77 
78  // Call the AMD routine
79  // m_mat.prune(keep_diag()); //Remove the diagonal elements
81  }
82 };
83 
84 #endif // EIGEN_MPL2_ONLY
85 
94 template <typename StorageIndex>
96 {
97  public:
99 
101  template <typename MatrixType>
102  void operator()(const MatrixType& /*mat*/, PermutationType& perm)
103  {
104  perm.resize(0);
105  }
106 
107 };
108 
117 template<typename StorageIndex>
119 {
120  public:
123 
127  template <typename MatrixType>
129  {
130  eigen_assert(mat.isCompressed() && "COLAMDOrdering requires a sparse matrix in compressed mode. Call .makeCompressed() before passing it to COLAMDOrdering");
131 
132  StorageIndex m = StorageIndex(mat.rows());
133  StorageIndex n = StorageIndex(mat.cols());
134  StorageIndex nnz = StorageIndex(mat.nonZeros());
135  // Get the recommended value of Alen to be used by colamd
136  StorageIndex Alen = internal::colamd_recommended(nnz, m, n);
137  // Set the default parameters
138  double knobs [COLAMD_KNOBS];
139  StorageIndex stats [COLAMD_STATS];
141 
142  IndexVector p(n+1), A(Alen);
143  for(StorageIndex i=0; i <= n; i++) p(i) = mat.outerIndexPtr()[i];
144  for(StorageIndex i=0; i < nnz; i++) A(i) = mat.innerIndexPtr()[i];
145  // Call Colamd routine to compute the ordering
146  StorageIndex info = internal::colamd(m, n, Alen, A.data(), p.data(), knobs, stats);
148  eigen_assert( info && "COLAMD failed " );
149 
150  perm.resize(n);
151  for (StorageIndex i = 0; i < n; i++) perm.indices()(p(i)) = i;
152  }
153 };
154 
155 } // end namespace Eigen
156 
157 #endif
COLAMD_STATS
#define COLAMD_STATS
Definition: Eigen_Colamd.h:63
Eigen
Definition: common.h:73
Eigen::SparseMatrix
A versatible sparse matrix representation.
Definition: SparseMatrix.h:96
MatrixType
Map< Matrix< Scalar, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > MatrixType
Definition: common.h:95
Eigen::PermutationMatrix::indices
const IndicesType & indices() const
Definition: PermutationMatrix.h:388
Eigen::COLAMDOrdering::PermutationType
PermutationMatrix< Dynamic, Dynamic, StorageIndex > PermutationType
Definition: Ordering.h:121
symm
int EIGEN_BLAS_FUNC() symm(const char *side, const char *uplo, const int *m, const int *n, const RealScalar *palpha, const RealScalar *pa, const int *lda, const RealScalar *pb, const int *ldb, const RealScalar *pbeta, RealScalar *pc, const int *ldc)
Definition: level3_impl.h:287
Eigen::COLAMDOrdering::operator()
void operator()(const MatrixType &mat, PermutationType &perm)
Definition: Ordering.h:128
eigen_assert
#define eigen_assert(x)
Definition: Macros.h:579
Eigen::AMDOrdering::PermutationType
PermutationMatrix< Dynamic, Dynamic, StorageIndex > PermutationType
Definition: Ordering.h:55
Eigen::PermutationBase::resize
void resize(Index newSize)
Definition: PermutationMatrix.h:136
Eigen::NaturalOrdering::operator()
void operator()(const MatrixType &, PermutationType &perm)
Definition: Ordering.h:102
Eigen::AMDOrdering::operator()
void operator()(const SparseSelfAdjointView< SrcType, SrcUpLo > &mat, PermutationType &perm)
Definition: Ordering.h:74
Eigen::AMDOrdering::operator()
void operator()(const MatrixType &mat, PermutationType &perm)
Definition: Ordering.h:61
info
else if n * info
Definition: cholesky.cpp:18
Eigen::internal::ordering_helper_at_plus_a
void ordering_helper_at_plus_a(const MatrixType &A, MatrixType &symmat)
Definition: Ordering.h:27
EIGEN_UNUSED_VARIABLE
#define EIGEN_UNUSED_VARIABLE(var)
Definition: Macros.h:618
Eigen::AMDOrdering
Definition: Ordering.h:52
COLAMD_KNOBS
#define COLAMD_KNOBS
Definition: Eigen_Colamd.h:60
Eigen::Map< Matrix< Scalar, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> >
Eigen::COLAMDOrdering
Definition: Ordering.h:118
Eigen::COLAMDOrdering::IndexVector
Matrix< StorageIndex, Dynamic, 1 > IndexVector
Definition: Ordering.h:122
Eigen::PermutationMatrix< Dynamic, Dynamic, StorageIndex >
internal::colamd
static bool colamd(IndexType n_row, IndexType n_col, IndexType Alen, IndexType *A, IndexType *p, double knobs[COLAMD_KNOBS], IndexType stats[COLAMD_STATS])
Computes a column ordering using the column approximate minimum degree ordering.
Definition: Eigen_Colamd.h:322
Eigen::internal::minimum_degree_ordering
void minimum_degree_ordering(SparseMatrix< Scalar, ColMajor, StorageIndex > &C, PermutationMatrix< Dynamic, Dynamic, StorageIndex > &perm)
Definition: Amd.h:94
Eigen_Colamd.h
internal::colamd_recommended
IndexType colamd_recommended(IndexType nnz, IndexType n_row, IndexType n_col)
Returns the recommended value of Alen.
Definition: Eigen_Colamd.h:257
internal::colamd_set_defaults
static void colamd_set_defaults(double knobs[COLAMD_KNOBS])
set default parameters The use of this routine is optional.
Definition: Eigen_Colamd.h:286
Eigen::SparseSelfAdjointView
Pseudo expression to manipulate a triangular sparse matrix as a selfadjoint matrix.
Definition: SparseSelfAdjointView.h:43
Eigen::Matrix< StorageIndex, Dynamic, 1 >
internal
Definition: BandTriangularSolver.h:13
n
PlainMatrixType mat * n
Definition: eigenvalues.cpp:41
A
MatrixType A(a, *n, *n, *lda)
Eigen::NaturalOrdering::PermutationType
PermutationMatrix< Dynamic, Dynamic, StorageIndex > PermutationType
Definition: Ordering.h:98
Eigen::NaturalOrdering
Definition: Ordering.h:95
mat
else mat
Definition: eigenvalues.cpp:43


control_box_rst
Author(s): Christoph Rösmann
autogenerated on Wed Mar 2 2022 00:05:58