Eigen::COLAMDOrdering< StorageIndex > Class Template Reference

#include <Ordering.h>

Public Types

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

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()).

Member Typedef Documentation

◆ IndexVector

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

◆ PermutationType

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

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()).
118  {
120  "COLAMDOrdering requires a sparse matrix in compressed mode. Call .makeCompressed() before passing it "
121  "to COLAMDOrdering");
122 
123  StorageIndex m = StorageIndex(mat.rows());
124  StorageIndex n = StorageIndex(mat.cols());
125  StorageIndex nnz = StorageIndex(mat.nonZeros());
126  // Get the recommended value of Alen to be used by colamd
127  StorageIndex Alen = internal::Colamd::recommended(nnz, m, n);
128  // Set the default parameters
129  double knobs[internal::Colamd::NKnobs];
130  StorageIndex stats[internal::Colamd::NStats];
132 
133  IndexVector p(n + 1), A(Alen);
134  for (StorageIndex i = 0; i <= n; i++) p(i) = mat.outerIndexPtr()[i];
135  for (StorageIndex i = 0; i < nnz; i++) A(i) = mat.innerIndexPtr()[i];
136  // Call Colamd routine to compute the ordering
137  StorageIndex info = internal::Colamd::compute_ordering(m, n, Alen, A.data(), p.data(), knobs, stats);
139  eigen_assert(info && "COLAMD failed ");
140 
141  perm.resize(n);
142  for (StorageIndex i = 0; i < n; i++) perm.indices()(p(i)) = i;
143  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9
const unsigned n
Definition: CG3DPackingUnitTest.cpp:11
#define EIGEN_UNUSED_VARIABLE(var)
Definition: Macros.h:966
#define eigen_assert(x)
Definition: Macros.h:910
float * p
Definition: Tutorial_Map_using.cpp:9
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
Matrix< StorageIndex, Dynamic, 1 > IndexVector
Definition: Ordering.h:112
Index nonZeros() const
Definition: SparseCompressedBase.h:64
Index cols() const
Definition: SparseMatrix.h:161
bool isCompressed() const
Definition: SparseCompressedBase.h:114
Index rows() const
Definition: SparseMatrix.h:159
const StorageIndex * outerIndexPtr() const
Definition: SparseMatrix.h:189
const StorageIndex * innerIndexPtr() const
Definition: SparseMatrix.h:180
int * m
Definition: level2_cplx_impl.h:294
int info
Definition: level2_cplx_impl.h:39
const int NKnobs
Definition: Eigen_Colamd.h:63
IndexType recommended(IndexType nnz, IndexType n_row, IndexType n_col)
Returns the recommended value of Alen.
Definition: Eigen_Colamd.h:267
static bool compute_ordering(IndexType n_row, IndexType n_col, IndexType Alen, IndexType *A, IndexType *p, double knobs[NKnobs], IndexType stats[NStats])
Computes a column ordering using the column approximate minimum degree ordering.
Definition: Eigen_Colamd.h:328
const int NStats
Definition: Eigen_Colamd.h:66
static void set_defaults(double knobs[NKnobs])
set default parameters The use of this routine is optional.
Definition: Eigen_Colamd.h:295

References Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::cols(), internal::Colamd::compute_ordering(), Eigen::PlainObjectBase< Derived >::data(), eigen_assert, EIGEN_UNUSED_VARIABLE, i, Eigen::PermutationMatrix< SizeAtCompileTime, MaxSizeAtCompileTime, StorageIndex_ >::indices(), info, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::innerIndexPtr(), Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::isCompressed(), m, n, internal::Colamd::NKnobs, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::nonZeros(), internal::Colamd::NStats, Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::outerIndexPtr(), p, internal::Colamd::recommended(), Eigen::PermutationBase< Derived >::resize(), Eigen::SparseMatrix< Scalar_, Options_, StorageIndex_ >::rows(), and internal::Colamd::set_defaults().


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