Home | History | Annotate | Download | only in CholmodSupport
      1 // This file is part of Eigen, a lightweight C++ template library
      2 // for linear algebra.
      3 //
      4 // Copyright (C) 2008-2010 Gael Guennebaud <gael.guennebaud (at) inria.fr>
      5 //
      6 // This Source Code Form is subject to the terms of the Mozilla
      7 // Public License v. 2.0. If a copy of the MPL was not distributed
      8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
      9 
     10 #ifndef EIGEN_CHOLMODSUPPORT_H
     11 #define EIGEN_CHOLMODSUPPORT_H
     12 
     13 namespace Eigen {
     14 
     15 namespace internal {
     16 
     17 template<typename Scalar> struct cholmod_configure_matrix;
     18 
     19 template<> struct cholmod_configure_matrix<double> {
     20   template<typename CholmodType>
     21   static void run(CholmodType& mat) {
     22     mat.xtype = CHOLMOD_REAL;
     23     mat.dtype = CHOLMOD_DOUBLE;
     24   }
     25 };
     26 
     27 template<> struct cholmod_configure_matrix<std::complex<double> > {
     28   template<typename CholmodType>
     29   static void run(CholmodType& mat) {
     30     mat.xtype = CHOLMOD_COMPLEX;
     31     mat.dtype = CHOLMOD_DOUBLE;
     32   }
     33 };
     34 
     35 // Other scalar types are not yet suppotred by Cholmod
     36 // template<> struct cholmod_configure_matrix<float> {
     37 //   template<typename CholmodType>
     38 //   static void run(CholmodType& mat) {
     39 //     mat.xtype = CHOLMOD_REAL;
     40 //     mat.dtype = CHOLMOD_SINGLE;
     41 //   }
     42 // };
     43 //
     44 // template<> struct cholmod_configure_matrix<std::complex<float> > {
     45 //   template<typename CholmodType>
     46 //   static void run(CholmodType& mat) {
     47 //     mat.xtype = CHOLMOD_COMPLEX;
     48 //     mat.dtype = CHOLMOD_SINGLE;
     49 //   }
     50 // };
     51 
     52 } // namespace internal
     53 
     54 /** Wraps the Eigen sparse matrix \a mat into a Cholmod sparse matrix object.
     55   * Note that the data are shared.
     56   */
     57 template<typename _Scalar, int _Options, typename _StorageIndex>
     58 cholmod_sparse viewAsCholmod(Ref<SparseMatrix<_Scalar,_Options,_StorageIndex> > mat)
     59 {
     60   cholmod_sparse res;
     61   res.nzmax   = mat.nonZeros();
     62   res.nrow    = mat.rows();
     63   res.ncol    = mat.cols();
     64   res.p       = mat.outerIndexPtr();
     65   res.i       = mat.innerIndexPtr();
     66   res.x       = mat.valuePtr();
     67   res.z       = 0;
     68   res.sorted  = 1;
     69   if(mat.isCompressed())
     70   {
     71     res.packed  = 1;
     72     res.nz = 0;
     73   }
     74   else
     75   {
     76     res.packed  = 0;
     77     res.nz = mat.innerNonZeroPtr();
     78   }
     79 
     80   res.dtype   = 0;
     81   res.stype   = -1;
     82 
     83   if (internal::is_same<_StorageIndex,int>::value)
     84   {
     85     res.itype = CHOLMOD_INT;
     86   }
     87   else if (internal::is_same<_StorageIndex,long>::value)
     88   {
     89     res.itype = CHOLMOD_LONG;
     90   }
     91   else
     92   {
     93     eigen_assert(false && "Index type not supported yet");
     94   }
     95 
     96   // setup res.xtype
     97   internal::cholmod_configure_matrix<_Scalar>::run(res);
     98 
     99   res.stype = 0;
    100 
    101   return res;
    102 }
    103 
    104 template<typename _Scalar, int _Options, typename _Index>
    105 const cholmod_sparse viewAsCholmod(const SparseMatrix<_Scalar,_Options,_Index>& mat)
    106 {
    107   cholmod_sparse res = viewAsCholmod(Ref<SparseMatrix<_Scalar,_Options,_Index> >(mat.const_cast_derived()));
    108   return res;
    109 }
    110 
    111 template<typename _Scalar, int _Options, typename _Index>
    112 const cholmod_sparse viewAsCholmod(const SparseVector<_Scalar,_Options,_Index>& mat)
    113 {
    114   cholmod_sparse res = viewAsCholmod(Ref<SparseMatrix<_Scalar,_Options,_Index> >(mat.const_cast_derived()));
    115   return res;
    116 }
    117 
    118 /** Returns a view of the Eigen sparse matrix \a mat as Cholmod sparse matrix.
    119   * The data are not copied but shared. */
    120 template<typename _Scalar, int _Options, typename _Index, unsigned int UpLo>
    121 cholmod_sparse viewAsCholmod(const SparseSelfAdjointView<const SparseMatrix<_Scalar,_Options,_Index>, UpLo>& mat)
    122 {
    123   cholmod_sparse res = viewAsCholmod(Ref<SparseMatrix<_Scalar,_Options,_Index> >(mat.matrix().const_cast_derived()));
    124 
    125   if(UpLo==Upper) res.stype =  1;
    126   if(UpLo==Lower) res.stype = -1;
    127 
    128   return res;
    129 }
    130 
    131 /** Returns a view of the Eigen \b dense matrix \a mat as Cholmod dense matrix.
    132   * The data are not copied but shared. */
    133 template<typename Derived>
    134 cholmod_dense viewAsCholmod(MatrixBase<Derived>& mat)
    135 {
    136   EIGEN_STATIC_ASSERT((internal::traits<Derived>::Flags&RowMajorBit)==0,THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
    137   typedef typename Derived::Scalar Scalar;
    138 
    139   cholmod_dense res;
    140   res.nrow   = mat.rows();
    141   res.ncol   = mat.cols();
    142   res.nzmax  = res.nrow * res.ncol;
    143   res.d      = Derived::IsVectorAtCompileTime ? mat.derived().size() : mat.derived().outerStride();
    144   res.x      = (void*)(mat.derived().data());
    145   res.z      = 0;
    146 
    147   internal::cholmod_configure_matrix<Scalar>::run(res);
    148 
    149   return res;
    150 }
    151 
    152 /** Returns a view of the Cholmod sparse matrix \a cm as an Eigen sparse matrix.
    153   * The data are not copied but shared. */
    154 template<typename Scalar, int Flags, typename StorageIndex>
    155 MappedSparseMatrix<Scalar,Flags,StorageIndex> viewAsEigen(cholmod_sparse& cm)
    156 {
    157   return MappedSparseMatrix<Scalar,Flags,StorageIndex>
    158          (cm.nrow, cm.ncol, static_cast<StorageIndex*>(cm.p)[cm.ncol],
    159           static_cast<StorageIndex*>(cm.p), static_cast<StorageIndex*>(cm.i),static_cast<Scalar*>(cm.x) );
    160 }
    161 
    162 enum CholmodMode {
    163   CholmodAuto, CholmodSimplicialLLt, CholmodSupernodalLLt, CholmodLDLt
    164 };
    165 
    166 
    167 /** \ingroup CholmodSupport_Module
    168   * \class CholmodBase
    169   * \brief The base class for the direct Cholesky factorization of Cholmod
    170   * \sa class CholmodSupernodalLLT, class CholmodSimplicialLDLT, class CholmodSimplicialLLT
    171   */
    172 template<typename _MatrixType, int _UpLo, typename Derived>
    173 class CholmodBase : public SparseSolverBase<Derived>
    174 {
    175   protected:
    176     typedef SparseSolverBase<Derived> Base;
    177     using Base::derived;
    178     using Base::m_isInitialized;
    179   public:
    180     typedef _MatrixType MatrixType;
    181     enum { UpLo = _UpLo };
    182     typedef typename MatrixType::Scalar Scalar;
    183     typedef typename MatrixType::RealScalar RealScalar;
    184     typedef MatrixType CholMatrixType;
    185     typedef typename MatrixType::StorageIndex StorageIndex;
    186     enum {
    187       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
    188       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
    189     };
    190 
    191   public:
    192 
    193     CholmodBase()
    194       : m_cholmodFactor(0), m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false)
    195     {
    196       EIGEN_STATIC_ASSERT((internal::is_same<double,RealScalar>::value), CHOLMOD_SUPPORTS_DOUBLE_PRECISION_ONLY);
    197       m_shiftOffset[0] = m_shiftOffset[1] = 0.0;
    198       cholmod_start(&m_cholmod);
    199     }
    200 
    201     explicit CholmodBase(const MatrixType& matrix)
    202       : m_cholmodFactor(0), m_info(Success), m_factorizationIsOk(false), m_analysisIsOk(false)
    203     {
    204       EIGEN_STATIC_ASSERT((internal::is_same<double,RealScalar>::value), CHOLMOD_SUPPORTS_DOUBLE_PRECISION_ONLY);
    205       m_shiftOffset[0] = m_shiftOffset[1] = 0.0;
    206       cholmod_start(&m_cholmod);
    207       compute(matrix);
    208     }
    209 
    210     ~CholmodBase()
    211     {
    212       if(m_cholmodFactor)
    213         cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
    214       cholmod_finish(&m_cholmod);
    215     }
    216 
    217     inline StorageIndex cols() const { return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); }
    218     inline StorageIndex rows() const { return internal::convert_index<StorageIndex, Index>(m_cholmodFactor->n); }
    219 
    220     /** \brief Reports whether previous computation was successful.
    221       *
    222       * \returns \c Success if computation was succesful,
    223       *          \c NumericalIssue if the matrix.appears to be negative.
    224       */
    225     ComputationInfo info() const
    226     {
    227       eigen_assert(m_isInitialized && "Decomposition is not initialized.");
    228       return m_info;
    229     }
    230 
    231     /** Computes the sparse Cholesky decomposition of \a matrix */
    232     Derived& compute(const MatrixType& matrix)
    233     {
    234       analyzePattern(matrix);
    235       factorize(matrix);
    236       return derived();
    237     }
    238 
    239     /** Performs a symbolic decomposition on the sparsity pattern of \a matrix.
    240       *
    241       * This function is particularly useful when solving for several problems having the same structure.
    242       *
    243       * \sa factorize()
    244       */
    245     void analyzePattern(const MatrixType& matrix)
    246     {
    247       if(m_cholmodFactor)
    248       {
    249         cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
    250         m_cholmodFactor = 0;
    251       }
    252       cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
    253       m_cholmodFactor = cholmod_analyze(&A, &m_cholmod);
    254 
    255       this->m_isInitialized = true;
    256       this->m_info = Success;
    257       m_analysisIsOk = true;
    258       m_factorizationIsOk = false;
    259     }
    260 
    261     /** Performs a numeric decomposition of \a matrix
    262       *
    263       * The given matrix must have the same sparsity pattern as the matrix on which the symbolic decomposition has been performed.
    264       *
    265       * \sa analyzePattern()
    266       */
    267     void factorize(const MatrixType& matrix)
    268     {
    269       eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
    270       cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
    271       cholmod_factorize_p(&A, m_shiftOffset, 0, 0, m_cholmodFactor, &m_cholmod);
    272 
    273       // If the factorization failed, minor is the column at which it did. On success minor == n.
    274       this->m_info = (m_cholmodFactor->minor == m_cholmodFactor->n ? Success : NumericalIssue);
    275       m_factorizationIsOk = true;
    276     }
    277 
    278     /** Returns a reference to the Cholmod's configuration structure to get a full control over the performed operations.
    279      *  See the Cholmod user guide for details. */
    280     cholmod_common& cholmod() { return m_cholmod; }
    281 
    282     #ifndef EIGEN_PARSED_BY_DOXYGEN
    283     /** \internal */
    284     template<typename Rhs,typename Dest>
    285     void _solve_impl(const MatrixBase<Rhs> &b, MatrixBase<Dest> &dest) const
    286     {
    287       eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
    288       const Index size = m_cholmodFactor->n;
    289       EIGEN_UNUSED_VARIABLE(size);
    290       eigen_assert(size==b.rows());
    291 
    292       // Cholmod needs column-major stoarge without inner-stride, which corresponds to the default behavior of Ref.
    293       Ref<const Matrix<typename Rhs::Scalar,Dynamic,Dynamic,ColMajor> > b_ref(b.derived());
    294 
    295       cholmod_dense b_cd = viewAsCholmod(b_ref);
    296       cholmod_dense* x_cd = cholmod_solve(CHOLMOD_A, m_cholmodFactor, &b_cd, &m_cholmod);
    297       if(!x_cd)
    298       {
    299         this->m_info = NumericalIssue;
    300         return;
    301       }
    302       // TODO optimize this copy by swapping when possible (be careful with alignment, etc.)
    303       dest = Matrix<Scalar,Dest::RowsAtCompileTime,Dest::ColsAtCompileTime>::Map(reinterpret_cast<Scalar*>(x_cd->x),b.rows(),b.cols());
    304       cholmod_free_dense(&x_cd, &m_cholmod);
    305     }
    306 
    307     /** \internal */
    308     template<typename RhsDerived, typename DestDerived>
    309     void _solve_impl(const SparseMatrixBase<RhsDerived> &b, SparseMatrixBase<DestDerived> &dest) const
    310     {
    311       eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
    312       const Index size = m_cholmodFactor->n;
    313       EIGEN_UNUSED_VARIABLE(size);
    314       eigen_assert(size==b.rows());
    315 
    316       // note: cs stands for Cholmod Sparse
    317       Ref<SparseMatrix<typename RhsDerived::Scalar,ColMajor,typename RhsDerived::StorageIndex> > b_ref(b.const_cast_derived());
    318       cholmod_sparse b_cs = viewAsCholmod(b_ref);
    319       cholmod_sparse* x_cs = cholmod_spsolve(CHOLMOD_A, m_cholmodFactor, &b_cs, &m_cholmod);
    320       if(!x_cs)
    321       {
    322         this->m_info = NumericalIssue;
    323         return;
    324       }
    325       // TODO optimize this copy by swapping when possible (be careful with alignment, etc.)
    326       dest.derived() = viewAsEigen<typename DestDerived::Scalar,ColMajor,typename DestDerived::StorageIndex>(*x_cs);
    327       cholmod_free_sparse(&x_cs, &m_cholmod);
    328     }
    329     #endif // EIGEN_PARSED_BY_DOXYGEN
    330 
    331 
    332     /** Sets the shift parameter that will be used to adjust the diagonal coefficients during the numerical factorization.
    333       *
    334       * During the numerical factorization, an offset term is added to the diagonal coefficients:\n
    335       * \c d_ii = \a offset + \c d_ii
    336       *
    337       * The default is \a offset=0.
    338       *
    339       * \returns a reference to \c *this.
    340       */
    341     Derived& setShift(const RealScalar& offset)
    342     {
    343       m_shiftOffset[0] = double(offset);
    344       return derived();
    345     }
    346 
    347     /** \returns the determinant of the underlying matrix from the current factorization */
    348     Scalar determinant() const
    349     {
    350       using std::exp;
    351       return exp(logDeterminant());
    352     }
    353 
    354     /** \returns the log determinant of the underlying matrix from the current factorization */
    355     Scalar logDeterminant() const
    356     {
    357       using std::log;
    358       using numext::real;
    359       eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
    360 
    361       RealScalar logDet = 0;
    362       Scalar *x = static_cast<Scalar*>(m_cholmodFactor->x);
    363       if (m_cholmodFactor->is_super)
    364       {
    365         // Supernodal factorization stored as a packed list of dense column-major blocs,
    366         // as described by the following structure:
    367 
    368         // super[k] == index of the first column of the j-th super node
    369         StorageIndex *super = static_cast<StorageIndex*>(m_cholmodFactor->super);
    370         // pi[k] == offset to the description of row indices
    371         StorageIndex *pi = static_cast<StorageIndex*>(m_cholmodFactor->pi);
    372         // px[k] == offset to the respective dense block
    373         StorageIndex *px = static_cast<StorageIndex*>(m_cholmodFactor->px);
    374 
    375         Index nb_super_nodes = m_cholmodFactor->nsuper;
    376         for (Index k=0; k < nb_super_nodes; ++k)
    377         {
    378           StorageIndex ncols = super[k + 1] - super[k];
    379           StorageIndex nrows = pi[k + 1] - pi[k];
    380 
    381           Map<const Array<Scalar,1,Dynamic>, 0, InnerStride<> > sk(x + px[k], ncols, InnerStride<>(nrows+1));
    382           logDet += sk.real().log().sum();
    383         }
    384       }
    385       else
    386       {
    387         // Simplicial factorization stored as standard CSC matrix.
    388         StorageIndex *p = static_cast<StorageIndex*>(m_cholmodFactor->p);
    389         Index size = m_cholmodFactor->n;
    390         for (Index k=0; k<size; ++k)
    391           logDet += log(real( x[p[k]] ));
    392       }
    393       if (m_cholmodFactor->is_ll)
    394         logDet *= 2.0;
    395       return logDet;
    396     };
    397 
    398     template<typename Stream>
    399     void dumpMemory(Stream& /*s*/)
    400     {}
    401 
    402   protected:
    403     mutable cholmod_common m_cholmod;
    404     cholmod_factor* m_cholmodFactor;
    405     double m_shiftOffset[2];
    406     mutable ComputationInfo m_info;
    407     int m_factorizationIsOk;
    408     int m_analysisIsOk;
    409 };
    410 
    411 /** \ingroup CholmodSupport_Module
    412   * \class CholmodSimplicialLLT
    413   * \brief A simplicial direct Cholesky (LLT) factorization and solver based on Cholmod
    414   *
    415   * This class allows to solve for A.X = B sparse linear problems via a simplicial LL^T Cholesky factorization
    416   * using the Cholmod library.
    417   * This simplicial variant is equivalent to Eigen's built-in SimplicialLLT class. Therefore, it has little practical interest.
    418   * The sparse matrix A must be selfadjoint and positive definite. The vectors or matrices
    419   * X and B can be either dense or sparse.
    420   *
    421   * \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
    422   * \tparam _UpLo the triangular part that will be used for the computations. It can be Lower
    423   *               or Upper. Default is Lower.
    424   *
    425   * \implsparsesolverconcept
    426   *
    427   * This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
    428   *
    429   * \warning Only double precision real and complex scalar types are supported by Cholmod.
    430   *
    431   * \sa \ref TutorialSparseSolverConcept, class CholmodSupernodalLLT, class SimplicialLLT
    432   */
    433 template<typename _MatrixType, int _UpLo = Lower>
    434 class CholmodSimplicialLLT : public CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLLT<_MatrixType, _UpLo> >
    435 {
    436     typedef CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLLT> Base;
    437     using Base::m_cholmod;
    438 
    439   public:
    440 
    441     typedef _MatrixType MatrixType;
    442 
    443     CholmodSimplicialLLT() : Base() { init(); }
    444 
    445     CholmodSimplicialLLT(const MatrixType& matrix) : Base()
    446     {
    447       init();
    448       this->compute(matrix);
    449     }
    450 
    451     ~CholmodSimplicialLLT() {}
    452   protected:
    453     void init()
    454     {
    455       m_cholmod.final_asis = 0;
    456       m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
    457       m_cholmod.final_ll = 1;
    458     }
    459 };
    460 
    461 
    462 /** \ingroup CholmodSupport_Module
    463   * \class CholmodSimplicialLDLT
    464   * \brief A simplicial direct Cholesky (LDLT) factorization and solver based on Cholmod
    465   *
    466   * This class allows to solve for A.X = B sparse linear problems via a simplicial LDL^T Cholesky factorization
    467   * using the Cholmod library.
    468   * This simplicial variant is equivalent to Eigen's built-in SimplicialLDLT class. Therefore, it has little practical interest.
    469   * The sparse matrix A must be selfadjoint and positive definite. The vectors or matrices
    470   * X and B can be either dense or sparse.
    471   *
    472   * \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
    473   * \tparam _UpLo the triangular part that will be used for the computations. It can be Lower
    474   *               or Upper. Default is Lower.
    475   *
    476   * \implsparsesolverconcept
    477   *
    478   * This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
    479   *
    480   * \warning Only double precision real and complex scalar types are supported by Cholmod.
    481   *
    482   * \sa \ref TutorialSparseSolverConcept, class CholmodSupernodalLLT, class SimplicialLDLT
    483   */
    484 template<typename _MatrixType, int _UpLo = Lower>
    485 class CholmodSimplicialLDLT : public CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLDLT<_MatrixType, _UpLo> >
    486 {
    487     typedef CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLDLT> Base;
    488     using Base::m_cholmod;
    489 
    490   public:
    491 
    492     typedef _MatrixType MatrixType;
    493 
    494     CholmodSimplicialLDLT() : Base() { init(); }
    495 
    496     CholmodSimplicialLDLT(const MatrixType& matrix) : Base()
    497     {
    498       init();
    499       this->compute(matrix);
    500     }
    501 
    502     ~CholmodSimplicialLDLT() {}
    503   protected:
    504     void init()
    505     {
    506       m_cholmod.final_asis = 1;
    507       m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
    508     }
    509 };
    510 
    511 /** \ingroup CholmodSupport_Module
    512   * \class CholmodSupernodalLLT
    513   * \brief A supernodal Cholesky (LLT) factorization and solver based on Cholmod
    514   *
    515   * This class allows to solve for A.X = B sparse linear problems via a supernodal LL^T Cholesky factorization
    516   * using the Cholmod library.
    517   * This supernodal variant performs best on dense enough problems, e.g., 3D FEM, or very high order 2D FEM.
    518   * The sparse matrix A must be selfadjoint and positive definite. The vectors or matrices
    519   * X and B can be either dense or sparse.
    520   *
    521   * \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
    522   * \tparam _UpLo the triangular part that will be used for the computations. It can be Lower
    523   *               or Upper. Default is Lower.
    524   *
    525   * \implsparsesolverconcept
    526   *
    527   * This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
    528   *
    529   * \warning Only double precision real and complex scalar types are supported by Cholmod.
    530   *
    531   * \sa \ref TutorialSparseSolverConcept
    532   */
    533 template<typename _MatrixType, int _UpLo = Lower>
    534 class CholmodSupernodalLLT : public CholmodBase<_MatrixType, _UpLo, CholmodSupernodalLLT<_MatrixType, _UpLo> >
    535 {
    536     typedef CholmodBase<_MatrixType, _UpLo, CholmodSupernodalLLT> Base;
    537     using Base::m_cholmod;
    538 
    539   public:
    540 
    541     typedef _MatrixType MatrixType;
    542 
    543     CholmodSupernodalLLT() : Base() { init(); }
    544 
    545     CholmodSupernodalLLT(const MatrixType& matrix) : Base()
    546     {
    547       init();
    548       this->compute(matrix);
    549     }
    550 
    551     ~CholmodSupernodalLLT() {}
    552   protected:
    553     void init()
    554     {
    555       m_cholmod.final_asis = 1;
    556       m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
    557     }
    558 };
    559 
    560 /** \ingroup CholmodSupport_Module
    561   * \class CholmodDecomposition
    562   * \brief A general Cholesky factorization and solver based on Cholmod
    563   *
    564   * This class allows to solve for A.X = B sparse linear problems via a LL^T or LDL^T Cholesky factorization
    565   * using the Cholmod library. The sparse matrix A must be selfadjoint and positive definite. The vectors or matrices
    566   * X and B can be either dense or sparse.
    567   *
    568   * This variant permits to change the underlying Cholesky method at runtime.
    569   * On the other hand, it does not provide access to the result of the factorization.
    570   * The default is to let Cholmod automatically choose between a simplicial and supernodal factorization.
    571   *
    572   * \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
    573   * \tparam _UpLo the triangular part that will be used for the computations. It can be Lower
    574   *               or Upper. Default is Lower.
    575   *
    576   * \implsparsesolverconcept
    577   *
    578   * This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
    579   *
    580   * \warning Only double precision real and complex scalar types are supported by Cholmod.
    581   *
    582   * \sa \ref TutorialSparseSolverConcept
    583   */
    584 template<typename _MatrixType, int _UpLo = Lower>
    585 class CholmodDecomposition : public CholmodBase<_MatrixType, _UpLo, CholmodDecomposition<_MatrixType, _UpLo> >
    586 {
    587     typedef CholmodBase<_MatrixType, _UpLo, CholmodDecomposition> Base;
    588     using Base::m_cholmod;
    589 
    590   public:
    591 
    592     typedef _MatrixType MatrixType;
    593 
    594     CholmodDecomposition() : Base() { init(); }
    595 
    596     CholmodDecomposition(const MatrixType& matrix) : Base()
    597     {
    598       init();
    599       this->compute(matrix);
    600     }
    601 
    602     ~CholmodDecomposition() {}
    603 
    604     void setMode(CholmodMode mode)
    605     {
    606       switch(mode)
    607       {
    608         case CholmodAuto:
    609           m_cholmod.final_asis = 1;
    610           m_cholmod.supernodal = CHOLMOD_AUTO;
    611           break;
    612         case CholmodSimplicialLLt:
    613           m_cholmod.final_asis = 0;
    614           m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
    615           m_cholmod.final_ll = 1;
    616           break;
    617         case CholmodSupernodalLLt:
    618           m_cholmod.final_asis = 1;
    619           m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
    620           break;
    621         case CholmodLDLt:
    622           m_cholmod.final_asis = 1;
    623           m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
    624           break;
    625         default:
    626           break;
    627       }
    628     }
    629   protected:
    630     void init()
    631     {
    632       m_cholmod.final_asis = 1;
    633       m_cholmod.supernodal = CHOLMOD_AUTO;
    634     }
    635 };
    636 
    637 } // end namespace Eigen
    638 
    639 #endif // EIGEN_CHOLMODSUPPORT_H
    640