• Home
  • Line#
  • Scopes#
  • Navigate#
  • Raw
  • Download
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@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, typename CholmodType>
cholmod_configure_matrix(CholmodType & mat)18 void cholmod_configure_matrix(CholmodType& mat)
19 {
20   if (internal::is_same<Scalar,float>::value)
21   {
22     mat.xtype = CHOLMOD_REAL;
23     mat.dtype = CHOLMOD_SINGLE;
24   }
25   else if (internal::is_same<Scalar,double>::value)
26   {
27     mat.xtype = CHOLMOD_REAL;
28     mat.dtype = CHOLMOD_DOUBLE;
29   }
30   else if (internal::is_same<Scalar,std::complex<float> >::value)
31   {
32     mat.xtype = CHOLMOD_COMPLEX;
33     mat.dtype = CHOLMOD_SINGLE;
34   }
35   else if (internal::is_same<Scalar,std::complex<double> >::value)
36   {
37     mat.xtype = CHOLMOD_COMPLEX;
38     mat.dtype = CHOLMOD_DOUBLE;
39   }
40   else
41   {
42     eigen_assert(false && "Scalar type not supported by CHOLMOD");
43   }
44 }
45 
46 } // namespace internal
47 
48 /** Wraps the Eigen sparse matrix \a mat into a Cholmod sparse matrix object.
49   * Note that the data are shared.
50   */
51 template<typename _Scalar, int _Options, typename _Index>
viewAsCholmod(SparseMatrix<_Scalar,_Options,_Index> & mat)52 cholmod_sparse viewAsCholmod(SparseMatrix<_Scalar,_Options,_Index>& mat)
53 {
54   typedef SparseMatrix<_Scalar,_Options,_Index> MatrixType;
55   cholmod_sparse res;
56   res.nzmax   = mat.nonZeros();
57   res.nrow    = mat.rows();;
58   res.ncol    = mat.cols();
59   res.p       = mat.outerIndexPtr();
60   res.i       = mat.innerIndexPtr();
61   res.x       = mat.valuePtr();
62   res.sorted  = 1;
63   if(mat.isCompressed())
64   {
65     res.packed  = 1;
66   }
67   else
68   {
69     res.packed  = 0;
70     res.nz = mat.innerNonZeroPtr();
71   }
72 
73   res.dtype   = 0;
74   res.stype   = -1;
75 
76   if (internal::is_same<_Index,int>::value)
77   {
78     res.itype = CHOLMOD_INT;
79   }
80   else
81   {
82     eigen_assert(false && "Index type different than int is not supported yet");
83   }
84 
85   // setup res.xtype
86   internal::cholmod_configure_matrix<_Scalar>(res);
87 
88   res.stype = 0;
89 
90   return res;
91 }
92 
93 template<typename _Scalar, int _Options, typename _Index>
viewAsCholmod(const SparseMatrix<_Scalar,_Options,_Index> & mat)94 const cholmod_sparse viewAsCholmod(const SparseMatrix<_Scalar,_Options,_Index>& mat)
95 {
96   cholmod_sparse res = viewAsCholmod(mat.const_cast_derived());
97   return res;
98 }
99 
100 /** Returns a view of the Eigen sparse matrix \a mat as Cholmod sparse matrix.
101   * The data are not copied but shared. */
102 template<typename _Scalar, int _Options, typename _Index, unsigned int UpLo>
viewAsCholmod(const SparseSelfAdjointView<SparseMatrix<_Scalar,_Options,_Index>,UpLo> & mat)103 cholmod_sparse viewAsCholmod(const SparseSelfAdjointView<SparseMatrix<_Scalar,_Options,_Index>, UpLo>& mat)
104 {
105   cholmod_sparse res = viewAsCholmod(mat.matrix().const_cast_derived());
106 
107   if(UpLo==Upper) res.stype =  1;
108   if(UpLo==Lower) res.stype = -1;
109 
110   return res;
111 }
112 
113 /** Returns a view of the Eigen \b dense matrix \a mat as Cholmod dense matrix.
114   * The data are not copied but shared. */
115 template<typename Derived>
viewAsCholmod(MatrixBase<Derived> & mat)116 cholmod_dense viewAsCholmod(MatrixBase<Derived>& mat)
117 {
118   EIGEN_STATIC_ASSERT((internal::traits<Derived>::Flags&RowMajorBit)==0,THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
119   typedef typename Derived::Scalar Scalar;
120 
121   cholmod_dense res;
122   res.nrow   = mat.rows();
123   res.ncol   = mat.cols();
124   res.nzmax  = res.nrow * res.ncol;
125   res.d      = Derived::IsVectorAtCompileTime ? mat.derived().size() : mat.derived().outerStride();
126   res.x      = mat.derived().data();
127   res.z      = 0;
128 
129   internal::cholmod_configure_matrix<Scalar>(res);
130 
131   return res;
132 }
133 
134 /** Returns a view of the Cholmod sparse matrix \a cm as an Eigen sparse matrix.
135   * The data are not copied but shared. */
136 template<typename Scalar, int Flags, typename Index>
viewAsEigen(cholmod_sparse & cm)137 MappedSparseMatrix<Scalar,Flags,Index> viewAsEigen(cholmod_sparse& cm)
138 {
139   return MappedSparseMatrix<Scalar,Flags,Index>
140          (cm.nrow, cm.ncol, reinterpret_cast<Index*>(cm.p)[cm.ncol],
141           reinterpret_cast<Index*>(cm.p), reinterpret_cast<Index*>(cm.i),reinterpret_cast<Scalar*>(cm.x) );
142 }
143 
144 enum CholmodMode {
145   CholmodAuto, CholmodSimplicialLLt, CholmodSupernodalLLt, CholmodLDLt
146 };
147 
148 
149 /** \ingroup CholmodSupport_Module
150   * \class CholmodBase
151   * \brief The base class for the direct Cholesky factorization of Cholmod
152   * \sa class CholmodSupernodalLLT, class CholmodSimplicialLDLT, class CholmodSimplicialLLT
153   */
154 template<typename _MatrixType, int _UpLo, typename Derived>
155 class CholmodBase : internal::noncopyable
156 {
157   public:
158     typedef _MatrixType MatrixType;
159     enum { UpLo = _UpLo };
160     typedef typename MatrixType::Scalar Scalar;
161     typedef typename MatrixType::RealScalar RealScalar;
162     typedef MatrixType CholMatrixType;
163     typedef typename MatrixType::Index Index;
164 
165   public:
166 
CholmodBase()167     CholmodBase()
168       : m_cholmodFactor(0), m_info(Success), m_isInitialized(false)
169     {
170       cholmod_start(&m_cholmod);
171     }
172 
CholmodBase(const MatrixType & matrix)173     CholmodBase(const MatrixType& matrix)
174       : m_cholmodFactor(0), m_info(Success), m_isInitialized(false)
175     {
176       cholmod_start(&m_cholmod);
177       compute(matrix);
178     }
179 
~CholmodBase()180     ~CholmodBase()
181     {
182       if(m_cholmodFactor)
183         cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
184       cholmod_finish(&m_cholmod);
185     }
186 
cols()187     inline Index cols() const { return m_cholmodFactor->n; }
rows()188     inline Index rows() const { return m_cholmodFactor->n; }
189 
derived()190     Derived& derived() { return *static_cast<Derived*>(this); }
derived()191     const Derived& derived() const { return *static_cast<const Derived*>(this); }
192 
193     /** \brief Reports whether previous computation was successful.
194       *
195       * \returns \c Success if computation was succesful,
196       *          \c NumericalIssue if the matrix.appears to be negative.
197       */
info()198     ComputationInfo info() const
199     {
200       eigen_assert(m_isInitialized && "Decomposition is not initialized.");
201       return m_info;
202     }
203 
204     /** Computes the sparse Cholesky decomposition of \a matrix */
compute(const MatrixType & matrix)205     Derived& compute(const MatrixType& matrix)
206     {
207       analyzePattern(matrix);
208       factorize(matrix);
209       return derived();
210     }
211 
212     /** \returns the solution x of \f$ A x = b \f$ using the current decomposition of A.
213       *
214       * \sa compute()
215       */
216     template<typename Rhs>
217     inline const internal::solve_retval<CholmodBase, Rhs>
solve(const MatrixBase<Rhs> & b)218     solve(const MatrixBase<Rhs>& b) const
219     {
220       eigen_assert(m_isInitialized && "LLT is not initialized.");
221       eigen_assert(rows()==b.rows()
222                 && "CholmodDecomposition::solve(): invalid number of rows of the right hand side matrix b");
223       return internal::solve_retval<CholmodBase, Rhs>(*this, b.derived());
224     }
225 
226     /** \returns the solution x of \f$ A x = b \f$ using the current decomposition of A.
227       *
228       * \sa compute()
229       */
230     template<typename Rhs>
231     inline const internal::sparse_solve_retval<CholmodBase, Rhs>
solve(const SparseMatrixBase<Rhs> & b)232     solve(const SparseMatrixBase<Rhs>& b) const
233     {
234       eigen_assert(m_isInitialized && "LLT is not initialized.");
235       eigen_assert(rows()==b.rows()
236                 && "CholmodDecomposition::solve(): invalid number of rows of the right hand side matrix b");
237       return internal::sparse_solve_retval<CholmodBase, Rhs>(*this, b.derived());
238     }
239 
240     /** Performs a symbolic decomposition on the sparcity of \a matrix.
241       *
242       * This function is particularly useful when solving for several problems having the same structure.
243       *
244       * \sa factorize()
245       */
analyzePattern(const MatrixType & matrix)246     void analyzePattern(const MatrixType& matrix)
247     {
248       if(m_cholmodFactor)
249       {
250         cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
251         m_cholmodFactor = 0;
252       }
253       cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
254       m_cholmodFactor = cholmod_analyze(&A, &m_cholmod);
255 
256       this->m_isInitialized = true;
257       this->m_info = Success;
258       m_analysisIsOk = true;
259       m_factorizationIsOk = false;
260     }
261 
262     /** Performs a numeric decomposition of \a matrix
263       *
264       * The given matrix must has the same sparcity than the matrix on which the symbolic decomposition has been performed.
265       *
266       * \sa analyzePattern()
267       */
factorize(const MatrixType & matrix)268     void factorize(const MatrixType& matrix)
269     {
270       eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
271       cholmod_sparse A = viewAsCholmod(matrix.template selfadjointView<UpLo>());
272       cholmod_factorize(&A, m_cholmodFactor, &m_cholmod);
273 
274       this->m_info = Success;
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. */
cholmod()280     cholmod_common& cholmod() { return m_cholmod; }
281 
282     #ifndef EIGEN_PARSED_BY_DOXYGEN
283     /** \internal */
284     template<typename Rhs,typename Dest>
_solve(const MatrixBase<Rhs> & b,MatrixBase<Dest> & dest)285     void _solve(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_assert(size==b.rows());
290 
291       // note: cd stands for Cholmod Dense
292       cholmod_dense b_cd = viewAsCholmod(b.const_cast_derived());
293       cholmod_dense* x_cd = cholmod_solve(CHOLMOD_A, m_cholmodFactor, &b_cd, &m_cholmod);
294       if(!x_cd)
295       {
296         this->m_info = NumericalIssue;
297       }
298       // TODO optimize this copy by swapping when possible (be carreful with alignment, etc.)
299       dest = Matrix<Scalar,Dest::RowsAtCompileTime,Dest::ColsAtCompileTime>::Map(reinterpret_cast<Scalar*>(x_cd->x),b.rows(),b.cols());
300       cholmod_free_dense(&x_cd, &m_cholmod);
301     }
302 
303     /** \internal */
304     template<typename RhsScalar, int RhsOptions, typename RhsIndex, typename DestScalar, int DestOptions, typename DestIndex>
_solve(const SparseMatrix<RhsScalar,RhsOptions,RhsIndex> & b,SparseMatrix<DestScalar,DestOptions,DestIndex> & dest)305     void _solve(const SparseMatrix<RhsScalar,RhsOptions,RhsIndex> &b, SparseMatrix<DestScalar,DestOptions,DestIndex> &dest) const
306     {
307       eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or symbolic()/numeric()");
308       const Index size = m_cholmodFactor->n;
309       eigen_assert(size==b.rows());
310 
311       // note: cs stands for Cholmod Sparse
312       cholmod_sparse b_cs = viewAsCholmod(b);
313       cholmod_sparse* x_cs = cholmod_spsolve(CHOLMOD_A, m_cholmodFactor, &b_cs, &m_cholmod);
314       if(!x_cs)
315       {
316         this->m_info = NumericalIssue;
317       }
318       // TODO optimize this copy by swapping when possible (be carreful with alignment, etc.)
319       dest = viewAsEigen<DestScalar,DestOptions,DestIndex>(*x_cs);
320       cholmod_free_sparse(&x_cs, &m_cholmod);
321     }
322     #endif // EIGEN_PARSED_BY_DOXYGEN
323 
324     template<typename Stream>
dumpMemory(Stream & s)325     void dumpMemory(Stream& s)
326     {}
327 
328   protected:
329     mutable cholmod_common m_cholmod;
330     cholmod_factor* m_cholmodFactor;
331     mutable ComputationInfo m_info;
332     bool m_isInitialized;
333     int m_factorizationIsOk;
334     int m_analysisIsOk;
335 };
336 
337 /** \ingroup CholmodSupport_Module
338   * \class CholmodSimplicialLLT
339   * \brief A simplicial direct Cholesky (LLT) factorization and solver based on Cholmod
340   *
341   * This class allows to solve for A.X = B sparse linear problems via a simplicial LL^T Cholesky factorization
342   * using the Cholmod library.
343   * This simplicial variant is equivalent to Eigen's built-in SimplicialLLT class. Thefore, it has little practical interest.
344   * The sparse matrix A must be selfajoint and positive definite. The vectors or matrices
345   * X and B can be either dense or sparse.
346   *
347   * \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
348   * \tparam _UpLo the triangular part that will be used for the computations. It can be Lower
349   *               or Upper. Default is Lower.
350   *
351   * This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
352   *
353   * \sa \ref TutorialSparseDirectSolvers, class CholmodSupernodalLLT, class SimplicialLLT
354   */
355 template<typename _MatrixType, int _UpLo = Lower>
356 class CholmodSimplicialLLT : public CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLLT<_MatrixType, _UpLo> >
357 {
358     typedef CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLLT> Base;
359     using Base::m_cholmod;
360 
361   public:
362 
363     typedef _MatrixType MatrixType;
364 
CholmodSimplicialLLT()365     CholmodSimplicialLLT() : Base() { init(); }
366 
CholmodSimplicialLLT(const MatrixType & matrix)367     CholmodSimplicialLLT(const MatrixType& matrix) : Base()
368     {
369       init();
370       compute(matrix);
371     }
372 
~CholmodSimplicialLLT()373     ~CholmodSimplicialLLT() {}
374   protected:
init()375     void init()
376     {
377       m_cholmod.final_asis = 0;
378       m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
379       m_cholmod.final_ll = 1;
380     }
381 };
382 
383 
384 /** \ingroup CholmodSupport_Module
385   * \class CholmodSimplicialLDLT
386   * \brief A simplicial direct Cholesky (LDLT) factorization and solver based on Cholmod
387   *
388   * This class allows to solve for A.X = B sparse linear problems via a simplicial LDL^T Cholesky factorization
389   * using the Cholmod library.
390   * This simplicial variant is equivalent to Eigen's built-in SimplicialLDLT class. Thefore, it has little practical interest.
391   * The sparse matrix A must be selfajoint and positive definite. The vectors or matrices
392   * X and B can be either dense or sparse.
393   *
394   * \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
395   * \tparam _UpLo the triangular part that will be used for the computations. It can be Lower
396   *               or Upper. Default is Lower.
397   *
398   * This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
399   *
400   * \sa \ref TutorialSparseDirectSolvers, class CholmodSupernodalLLT, class SimplicialLDLT
401   */
402 template<typename _MatrixType, int _UpLo = Lower>
403 class CholmodSimplicialLDLT : public CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLDLT<_MatrixType, _UpLo> >
404 {
405     typedef CholmodBase<_MatrixType, _UpLo, CholmodSimplicialLDLT> Base;
406     using Base::m_cholmod;
407 
408   public:
409 
410     typedef _MatrixType MatrixType;
411 
CholmodSimplicialLDLT()412     CholmodSimplicialLDLT() : Base() { init(); }
413 
CholmodSimplicialLDLT(const MatrixType & matrix)414     CholmodSimplicialLDLT(const MatrixType& matrix) : Base()
415     {
416       init();
417       compute(matrix);
418     }
419 
~CholmodSimplicialLDLT()420     ~CholmodSimplicialLDLT() {}
421   protected:
init()422     void init()
423     {
424       m_cholmod.final_asis = 1;
425       m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
426     }
427 };
428 
429 /** \ingroup CholmodSupport_Module
430   * \class CholmodSupernodalLLT
431   * \brief A supernodal Cholesky (LLT) factorization and solver based on Cholmod
432   *
433   * This class allows to solve for A.X = B sparse linear problems via a supernodal LL^T Cholesky factorization
434   * using the Cholmod library.
435   * This supernodal variant performs best on dense enough problems, e.g., 3D FEM, or very high order 2D FEM.
436   * The sparse matrix A must be selfajoint and positive definite. The vectors or matrices
437   * X and B can be either dense or sparse.
438   *
439   * \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
440   * \tparam _UpLo the triangular part that will be used for the computations. It can be Lower
441   *               or Upper. Default is Lower.
442   *
443   * This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
444   *
445   * \sa \ref TutorialSparseDirectSolvers
446   */
447 template<typename _MatrixType, int _UpLo = Lower>
448 class CholmodSupernodalLLT : public CholmodBase<_MatrixType, _UpLo, CholmodSupernodalLLT<_MatrixType, _UpLo> >
449 {
450     typedef CholmodBase<_MatrixType, _UpLo, CholmodSupernodalLLT> Base;
451     using Base::m_cholmod;
452 
453   public:
454 
455     typedef _MatrixType MatrixType;
456 
CholmodSupernodalLLT()457     CholmodSupernodalLLT() : Base() { init(); }
458 
CholmodSupernodalLLT(const MatrixType & matrix)459     CholmodSupernodalLLT(const MatrixType& matrix) : Base()
460     {
461       init();
462       compute(matrix);
463     }
464 
~CholmodSupernodalLLT()465     ~CholmodSupernodalLLT() {}
466   protected:
init()467     void init()
468     {
469       m_cholmod.final_asis = 1;
470       m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
471     }
472 };
473 
474 /** \ingroup CholmodSupport_Module
475   * \class CholmodDecomposition
476   * \brief A general Cholesky factorization and solver based on Cholmod
477   *
478   * This class allows to solve for A.X = B sparse linear problems via a LL^T or LDL^T Cholesky factorization
479   * using the Cholmod library. The sparse matrix A must be selfajoint and positive definite. The vectors or matrices
480   * X and B can be either dense or sparse.
481   *
482   * This variant permits to change the underlying Cholesky method at runtime.
483   * On the other hand, it does not provide access to the result of the factorization.
484   * The default is to let Cholmod automatically choose between a simplicial and supernodal factorization.
485   *
486   * \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
487   * \tparam _UpLo the triangular part that will be used for the computations. It can be Lower
488   *               or Upper. Default is Lower.
489   *
490   * This class supports all kind of SparseMatrix<>: row or column major; upper, lower, or both; compressed or non compressed.
491   *
492   * \sa \ref TutorialSparseDirectSolvers
493   */
494 template<typename _MatrixType, int _UpLo = Lower>
495 class CholmodDecomposition : public CholmodBase<_MatrixType, _UpLo, CholmodDecomposition<_MatrixType, _UpLo> >
496 {
497     typedef CholmodBase<_MatrixType, _UpLo, CholmodDecomposition> Base;
498     using Base::m_cholmod;
499 
500   public:
501 
502     typedef _MatrixType MatrixType;
503 
CholmodDecomposition()504     CholmodDecomposition() : Base() { init(); }
505 
CholmodDecomposition(const MatrixType & matrix)506     CholmodDecomposition(const MatrixType& matrix) : Base()
507     {
508       init();
509       compute(matrix);
510     }
511 
~CholmodDecomposition()512     ~CholmodDecomposition() {}
513 
setMode(CholmodMode mode)514     void setMode(CholmodMode mode)
515     {
516       switch(mode)
517       {
518         case CholmodAuto:
519           m_cholmod.final_asis = 1;
520           m_cholmod.supernodal = CHOLMOD_AUTO;
521           break;
522         case CholmodSimplicialLLt:
523           m_cholmod.final_asis = 0;
524           m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
525           m_cholmod.final_ll = 1;
526           break;
527         case CholmodSupernodalLLt:
528           m_cholmod.final_asis = 1;
529           m_cholmod.supernodal = CHOLMOD_SUPERNODAL;
530           break;
531         case CholmodLDLt:
532           m_cholmod.final_asis = 1;
533           m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
534           break;
535         default:
536           break;
537       }
538     }
539   protected:
init()540     void init()
541     {
542       m_cholmod.final_asis = 1;
543       m_cholmod.supernodal = CHOLMOD_AUTO;
544     }
545 };
546 
547 namespace internal {
548 
549 template<typename _MatrixType, int _UpLo, typename Derived, typename Rhs>
550 struct solve_retval<CholmodBase<_MatrixType,_UpLo,Derived>, Rhs>
551   : solve_retval_base<CholmodBase<_MatrixType,_UpLo,Derived>, Rhs>
552 {
553   typedef CholmodBase<_MatrixType,_UpLo,Derived> Dec;
554   EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
555 
556   template<typename Dest> void evalTo(Dest& dst) const
557   {
558     dec()._solve(rhs(),dst);
559   }
560 };
561 
562 template<typename _MatrixType, int _UpLo, typename Derived, typename Rhs>
563 struct sparse_solve_retval<CholmodBase<_MatrixType,_UpLo,Derived>, Rhs>
564   : sparse_solve_retval_base<CholmodBase<_MatrixType,_UpLo,Derived>, Rhs>
565 {
566   typedef CholmodBase<_MatrixType,_UpLo,Derived> Dec;
567   EIGEN_MAKE_SPARSE_SOLVE_HELPERS(Dec,Rhs)
568 
569   template<typename Dest> void evalTo(Dest& dst) const
570   {
571     dec()._solve(rhs(),dst);
572   }
573 };
574 
575 } // end namespace internal
576 
577 } // end namespace Eigen
578 
579 #endif // EIGEN_CHOLMODSUPPORT_H
580