Eigen::AMDOrdering< StorageIndex > Class Template Reference

Public Types

typedef PermutationMatrix< Dynamic, Dynamic, StorageIndex > PermutationType
 

Public Member Functions

template<typename MatrixType >
void operator() (const MatrixType &mat, PermutationType &perm)
 
template<typename SrcType , unsigned int SrcUpLo>
void operator() (const SparseSelfAdjointView< SrcType, SrcUpLo > &mat, PermutationType &perm)
 

Detailed Description

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

Functor computing the approximate minimum degree ordering If the matrix is not structurally symmetric, an ordering of A^T+A is computed

Template Parameters
StorageIndexThe type of indices of the matrix
See also
COLAMDOrdering

Definition at line 52 of file Ordering.h.

Member Typedef Documentation

◆ PermutationType

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

Definition at line 55 of file Ordering.h.

Member Function Documentation

◆ operator()() [1/2]

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

Compute the permutation vector from a sparse matrix This routine is much faster if the input matrix is column-major

Definition at line 61 of file Ordering.h.

62  {
63  // Compute the symmetric pattern
64  SparseMatrix<typename MatrixType::Scalar, ColMajor, StorageIndex> symm;
66 
67  // Call the AMD routine
68  //m_mat.prune(keep_diag());
70  }
void ordering_helper_at_plus_a(const MatrixType &A, MatrixType &symmat)
Definition: Ordering.h:29
void minimum_degree_ordering(SparseMatrix< Scalar, ColMajor, StorageIndex > &C, PermutationMatrix< Dynamic, Dynamic, StorageIndex > &perm)
Definition: Amd.h:86

◆ operator()() [2/2]

template<typename StorageIndex >
template<typename SrcType , unsigned int SrcUpLo>
void Eigen::AMDOrdering< StorageIndex >::operator() ( const SparseSelfAdjointView< SrcType, SrcUpLo > &  mat,
PermutationType perm 
)
inline

Compute the permutation with a selfadjoint matrix

Definition at line 74 of file Ordering.h.

75  {
76  SparseMatrix<typename SrcType::Scalar, ColMajor, StorageIndex> C; C = mat;
77 
78  // Call the AMD routine
79  // m_mat.prune(keep_diag()); //Remove the diagonal elements
81  }

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