xref: /aosp_15_r20/external/eigen/Eigen/src/OrderingMethods/Ordering.h (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
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 <[email protected]>
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 
20 /** \internal
21   * \ingroup OrderingMethods_Module
22   * \param[in] A the input non-symmetric matrix
23   * \param[out] symmat the symmetric pattern A^T+A from the input matrix \a A.
24   * FIXME: The values should not be considered here
25   */
26 template<typename MatrixType>
ordering_helper_at_plus_a(const MatrixType & A,MatrixType & symmat)27 void ordering_helper_at_plus_a(const MatrixType& A, MatrixType& symmat)
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() = typename MatrixType::Scalar(0);
35   }
36   symmat = C + A;
37 }
38 
39 }
40 
41 /** \ingroup OrderingMethods_Module
42   * \class AMDOrdering
43   *
44   * Functor computing the \em approximate \em minimum \em degree ordering
45   * If the matrix is not structurally symmetric, an ordering of A^T+A is computed
46   * \tparam  StorageIndex The type of indices of the matrix
47   * \sa COLAMDOrdering
48   */
49 template <typename StorageIndex>
50 class AMDOrdering
51 {
52   public:
53     typedef PermutationMatrix<Dynamic, Dynamic, StorageIndex> PermutationType;
54 
55     /** Compute the permutation vector from a sparse matrix
56      * This routine is much faster if the input matrix is column-major
57      */
58     template <typename MatrixType>
operator()59     void operator()(const MatrixType& mat, PermutationType& perm)
60     {
61       // Compute the symmetric pattern
62       SparseMatrix<typename MatrixType::Scalar, ColMajor, StorageIndex> symm;
63       internal::ordering_helper_at_plus_a(mat,symm);
64 
65       // Call the AMD routine
66       //m_mat.prune(keep_diag());
67       internal::minimum_degree_ordering(symm, perm);
68     }
69 
70     /** Compute the permutation with a selfadjoint matrix */
71     template <typename SrcType, unsigned int SrcUpLo>
operator()72     void operator()(const SparseSelfAdjointView<SrcType, SrcUpLo>& mat, PermutationType& perm)
73     {
74       SparseMatrix<typename SrcType::Scalar, ColMajor, StorageIndex> C; C = mat;
75 
76       // Call the AMD routine
77       // m_mat.prune(keep_diag()); //Remove the diagonal elements
78       internal::minimum_degree_ordering(C, perm);
79     }
80 };
81 
82 /** \ingroup OrderingMethods_Module
83   * \class NaturalOrdering
84   *
85   * Functor computing the natural ordering (identity)
86   *
87   * \note Returns an empty permutation matrix
88   * \tparam  StorageIndex The type of indices of the matrix
89   */
90 template <typename StorageIndex>
91 class NaturalOrdering
92 {
93   public:
94     typedef PermutationMatrix<Dynamic, Dynamic, StorageIndex> PermutationType;
95 
96     /** Compute the permutation vector from a column-major sparse matrix */
97     template <typename MatrixType>
operator()98     void operator()(const MatrixType& /*mat*/, PermutationType& perm)
99     {
100       perm.resize(0);
101     }
102 
103 };
104 
105 /** \ingroup OrderingMethods_Module
106   * \class COLAMDOrdering
107   *
108   * \tparam  StorageIndex The type of indices of the matrix
109   *
110   * Functor computing the \em column \em approximate \em minimum \em degree ordering
111   * The matrix should be in column-major and \b compressed format (see SparseMatrix::makeCompressed()).
112   */
113 template<typename StorageIndex>
114 class COLAMDOrdering
115 {
116   public:
117     typedef PermutationMatrix<Dynamic, Dynamic, StorageIndex> PermutationType;
118     typedef Matrix<StorageIndex, Dynamic, 1> IndexVector;
119 
120     /** Compute the permutation vector \a perm form the sparse matrix \a mat
121       * \warning The input sparse matrix \a mat must be in compressed mode (see SparseMatrix::makeCompressed()).
122       */
123     template <typename MatrixType>
operator()124     void operator() (const MatrixType& mat, PermutationType& perm)
125     {
126       eigen_assert(mat.isCompressed() && "COLAMDOrdering requires a sparse matrix in compressed mode. Call .makeCompressed() before passing it to COLAMDOrdering");
127 
128       StorageIndex m = StorageIndex(mat.rows());
129       StorageIndex n = StorageIndex(mat.cols());
130       StorageIndex nnz = StorageIndex(mat.nonZeros());
131       // Get the recommended value of Alen to be used by colamd
132       StorageIndex Alen = internal::Colamd::recommended(nnz, m, n);
133       // Set the default parameters
134       double knobs [internal::Colamd::NKnobs];
135       StorageIndex stats [internal::Colamd::NStats];
136       internal::Colamd::set_defaults(knobs);
137 
138       IndexVector p(n+1), A(Alen);
139       for(StorageIndex i=0; i <= n; i++)   p(i) = mat.outerIndexPtr()[i];
140       for(StorageIndex i=0; i < nnz; i++)  A(i) = mat.innerIndexPtr()[i];
141       // Call Colamd routine to compute the ordering
142       StorageIndex info = internal::Colamd::compute_ordering(m, n, Alen, A.data(), p.data(), knobs, stats);
143       EIGEN_UNUSED_VARIABLE(info);
144       eigen_assert( info && "COLAMD failed " );
145 
146       perm.resize(n);
147       for (StorageIndex i = 0; i < n; i++) perm.indices()(p(i)) = i;
148     }
149 };
150 
151 } // end namespace Eigen
152 
153 #endif
154