Eigen::COLAMDOrdering< StorageIndex > Class Template Reference

Public Types

typedef Matrix< StorageIndex, Dynamic, 1 > IndexVector
 
typedef PermutationMatrix< Dynamic, Dynamic, StorageIndex > PermutationType
 

Public Member Functions

template<typename MatrixType >
void operator() (const MatrixType &mat, PermutationType &perm)
 

Detailed Description

template<typename StorageIndex>
class Eigen::COLAMDOrdering< StorageIndex >

Template Parameters
StorageIndexThe type of indices of the matrix

Functor computing the column approximate minimum degree ordering The matrix should be in column-major and compressed format (see SparseMatrix::makeCompressed()).

Definition at line 116 of file Ordering.h.

Member Typedef Documentation

◆ IndexVector

template<typename StorageIndex >
typedef Matrix<StorageIndex, Dynamic, 1> Eigen::COLAMDOrdering< StorageIndex >::IndexVector

Definition at line 120 of file Ordering.h.

◆ PermutationType

template<typename StorageIndex >
typedef PermutationMatrix<Dynamic, Dynamic, StorageIndex> Eigen::COLAMDOrdering< StorageIndex >::PermutationType

Definition at line 119 of file Ordering.h.

Member Function Documentation

◆ operator()()

template<typename StorageIndex >
template<typename MatrixType >
void Eigen::COLAMDOrdering< StorageIndex >::operator() ( const MatrixType mat,
PermutationType perm 
)
inline

Compute the permutation vector perm form the sparse matrix mat

Warning
The input sparse matrix mat must be in compressed mode (see SparseMatrix::makeCompressed()).

Definition at line 126 of file Ordering.h.

127  {
128  eigen_assert(mat.isCompressed() && "COLAMDOrdering requires a sparse matrix in compressed mode. Call .makeCompressed() before passing it to COLAMDOrdering");
129 
130  StorageIndex m = StorageIndex(mat.rows());
131  StorageIndex n = StorageIndex(mat.cols());
132  StorageIndex nnz = StorageIndex(mat.nonZeros());
133  // Get the recommended value of Alen to be used by colamd
134  StorageIndex Alen = internal::Colamd::recommended(nnz, m, n);
135  // Set the default parameters
136  double knobs [internal::Colamd::NKnobs];
137  StorageIndex stats [internal::Colamd::NStats];
138  internal::Colamd::set_defaults(knobs);
139 
140  IndexVector p(n+1), A(Alen);
141  for(StorageIndex i=0; i <= n; i++) p(i) = mat.outerIndexPtr()[i];
142  for(StorageIndex i=0; i < nnz; i++) A(i) = mat.innerIndexPtr()[i];
143  // Call Colamd routine to compute the ordering
144  StorageIndex info = internal::Colamd::compute_ordering(m, n, Alen, A.data(), p.data(), knobs, stats);
145  EIGEN_UNUSED_VARIABLE(info);
146  eigen_assert( info && "COLAMD failed " );
147 
148  perm.resize(n);
149  for (StorageIndex i = 0; i < n; i++) perm.indices()(p(i)) = i;
150  }
Matrix3f m
int n
MatrixXcf A
#define EIGEN_UNUSED_VARIABLE(var)
Definition: Macros.h:957
#define eigen_assert(x)
Definition: Macros.h:902
float * p
Matrix< StorageIndex, Dynamic, 1 > IndexVector
Definition: Ordering.h:120

The documentation for this class was generated from the following file: