00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025 #ifndef EIGEN_CHOLMODSUPPORT_LEGACY_H
00026 #define EIGEN_CHOLMODSUPPORT_LEGACY_H
00027
00028 namespace internal {
00029
00030 template<typename Scalar, typename CholmodType>
00031 void cholmod_configure_matrix_legacy(CholmodType& mat)
00032 {
00033 if (internal::is_same<Scalar,float>::value)
00034 {
00035 mat.xtype = CHOLMOD_REAL;
00036 mat.dtype = CHOLMOD_SINGLE;
00037 }
00038 else if (internal::is_same<Scalar,double>::value)
00039 {
00040 mat.xtype = CHOLMOD_REAL;
00041 mat.dtype = CHOLMOD_DOUBLE;
00042 }
00043 else if (internal::is_same<Scalar,std::complex<float> >::value)
00044 {
00045 mat.xtype = CHOLMOD_COMPLEX;
00046 mat.dtype = CHOLMOD_SINGLE;
00047 }
00048 else if (internal::is_same<Scalar,std::complex<double> >::value)
00049 {
00050 mat.xtype = CHOLMOD_COMPLEX;
00051 mat.dtype = CHOLMOD_DOUBLE;
00052 }
00053 else
00054 {
00055 eigen_assert(false && "Scalar type not supported by CHOLMOD");
00056 }
00057 }
00058
00059 template<typename _MatrixType>
00060 cholmod_sparse cholmod_map_eigen_to_sparse(_MatrixType& mat)
00061 {
00062 typedef typename _MatrixType::Scalar Scalar;
00063 cholmod_sparse res;
00064 res.nzmax = mat.nonZeros();
00065 res.nrow = mat.rows();;
00066 res.ncol = mat.cols();
00067 res.p = mat._outerIndexPtr();
00068 res.i = mat._innerIndexPtr();
00069 res.x = mat._valuePtr();
00070 res.xtype = CHOLMOD_REAL;
00071 res.itype = CHOLMOD_INT;
00072 res.sorted = 1;
00073 res.packed = 1;
00074 res.dtype = 0;
00075 res.stype = -1;
00076
00077 internal::cholmod_configure_matrix_legacy<Scalar>(res);
00078
00079
00080 if (_MatrixType::Flags & SelfAdjoint)
00081 {
00082 if (_MatrixType::Flags & Upper)
00083 res.stype = 1;
00084 else if (_MatrixType::Flags & Lower)
00085 res.stype = -1;
00086 else
00087 res.stype = 0;
00088 }
00089 else
00090 res.stype = -1;
00091
00092 return res;
00093 }
00094
00095 template<typename Derived>
00096 cholmod_dense cholmod_map_eigen_to_dense(MatrixBase<Derived>& mat)
00097 {
00098 EIGEN_STATIC_ASSERT((internal::traits<Derived>::Flags&RowMajorBit)==0,THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
00099 typedef typename Derived::Scalar Scalar;
00100
00101 cholmod_dense res;
00102 res.nrow = mat.rows();
00103 res.ncol = mat.cols();
00104 res.nzmax = res.nrow * res.ncol;
00105 res.d = Derived::IsVectorAtCompileTime ? mat.derived().size() : mat.derived().outerStride();
00106 res.x = mat.derived().data();
00107 res.z = 0;
00108
00109 internal::cholmod_configure_matrix_legacy<Scalar>(res);
00110
00111 return res;
00112 }
00113
00114 template<typename Scalar, int Flags, typename Index>
00115 MappedSparseMatrix<Scalar,Flags,Index> map_cholmod_sparse_to_eigen(cholmod_sparse& cm)
00116 {
00117 return MappedSparseMatrix<Scalar,Flags,Index>
00118 (cm.nrow, cm.ncol, reinterpret_cast<Index*>(cm.p)[cm.ncol],
00119 reinterpret_cast<Index*>(cm.p), reinterpret_cast<Index*>(cm.i),reinterpret_cast<Scalar*>(cm.x) );
00120 }
00121
00122 }
00123
00124 template<typename _MatrixType>
00125 class SparseLLT<_MatrixType, Cholmod> : public SparseLLT<_MatrixType>
00126 {
00127 protected:
00128 typedef SparseLLT<_MatrixType> Base;
00129 typedef typename Base::Scalar Scalar;
00130 typedef typename Base::RealScalar RealScalar;
00131 typedef typename Base::CholMatrixType CholMatrixType;
00132 using Base::MatrixLIsDirty;
00133 using Base::SupernodalFactorIsDirty;
00134 using Base::m_flags;
00135 using Base::m_matrix;
00136 using Base::m_status;
00137
00138 public:
00139 typedef _MatrixType MatrixType;
00140 typedef typename MatrixType::Index Index;
00141
00142 SparseLLT(int flags = 0)
00143 : Base(flags), m_cholmodFactor(0)
00144 {
00145 cholmod_start(&m_cholmod);
00146 }
00147
00148 SparseLLT(const MatrixType& matrix, int flags = 0)
00149 : Base(flags), m_cholmodFactor(0)
00150 {
00151 cholmod_start(&m_cholmod);
00152 compute(matrix);
00153 }
00154
00155 ~SparseLLT()
00156 {
00157 if (m_cholmodFactor)
00158 cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
00159 cholmod_finish(&m_cholmod);
00160 }
00161
00162 inline const CholMatrixType& matrixL() const;
00163
00164 template<typename Derived>
00165 bool solveInPlace(MatrixBase<Derived> &b) const;
00166
00167 template<typename Rhs>
00168 inline const internal::solve_retval<SparseLLT<MatrixType, Cholmod>, Rhs>
00169 solve(const MatrixBase<Rhs>& b) const
00170 {
00171 eigen_assert(true && "SparseLLT is not initialized.");
00172 return internal::solve_retval<SparseLLT<MatrixType, Cholmod>, Rhs>(*this, b.derived());
00173 }
00174
00175 void compute(const MatrixType& matrix);
00176
00177 inline Index cols() const { return m_matrix.cols(); }
00178 inline Index rows() const { return m_matrix.rows(); }
00179
00180 inline const cholmod_factor* cholmodFactor() const
00181 { return m_cholmodFactor; }
00182
00183 inline cholmod_common* cholmodCommon() const
00184 { return &m_cholmod; }
00185
00186 bool succeeded() const;
00187
00188 protected:
00189 mutable cholmod_common m_cholmod;
00190 cholmod_factor* m_cholmodFactor;
00191 };
00192
00193
00194 namespace internal {
00195
00196 template<typename _MatrixType, typename Rhs>
00197 struct solve_retval<SparseLLT<_MatrixType, Cholmod>, Rhs>
00198 : solve_retval_base<SparseLLT<_MatrixType, Cholmod>, Rhs>
00199 {
00200 typedef SparseLLT<_MatrixType, Cholmod> SpLLTDecType;
00201 EIGEN_MAKE_SOLVE_HELPERS(SpLLTDecType,Rhs)
00202
00203 template<typename Dest> void evalTo(Dest& dst) const
00204 {
00205
00206 eigen_assert((Index)dec().cholmodFactor()->n==rhs().rows());
00207
00208 cholmod_factor* cholmodFactor = const_cast<cholmod_factor*>(dec().cholmodFactor());
00209 cholmod_common* cholmodCommon = const_cast<cholmod_common*>(dec().cholmodCommon());
00210
00211
00212
00213
00214
00215
00216 cholmod_dense cdb = internal::cholmod_map_eigen_to_dense(rhs().const_cast_derived());
00217 cholmod_dense* x = cholmod_solve(CHOLMOD_A, cholmodFactor, &cdb, cholmodCommon);
00218
00219 dst = Matrix<typename Base::Scalar,Dynamic,1>::Map(reinterpret_cast<typename Base::Scalar*>(x->x), rhs().rows());
00220
00221 cholmod_free_dense(&x, cholmodCommon);
00222
00223 }
00224
00225 };
00226
00227 }
00228
00229
00230
00231 template<typename _MatrixType>
00232 void SparseLLT<_MatrixType,Cholmod>::compute(const _MatrixType& a)
00233 {
00234 if (m_cholmodFactor)
00235 {
00236 cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
00237 m_cholmodFactor = 0;
00238 }
00239
00240 cholmod_sparse A = internal::cholmod_map_eigen_to_sparse(const_cast<_MatrixType&>(a));
00241
00242
00243
00244
00245
00246
00247
00248
00249
00250
00251
00252
00253
00254
00255
00256 m_cholmodFactor = cholmod_analyze(&A, &m_cholmod);
00257 cholmod_factorize(&A, m_cholmodFactor, &m_cholmod);
00258
00259 this->m_status = (this->m_status & ~Base::SupernodalFactorIsDirty) | Base::MatrixLIsDirty;
00260 }
00261
00262
00263
00264 template<typename _MatrixType>
00265 bool SparseLLT<_MatrixType,Cholmod>::succeeded() const
00266 { return true; }
00267
00268
00269
00270 template<typename _MatrixType>
00271 inline const typename SparseLLT<_MatrixType,Cholmod>::CholMatrixType&
00272 SparseLLT<_MatrixType,Cholmod>::matrixL() const
00273 {
00274 if (this->m_status & Base::MatrixLIsDirty)
00275 {
00276 eigen_assert(!(this->m_status & Base::SupernodalFactorIsDirty));
00277
00278 cholmod_sparse* cmRes = cholmod_factor_to_sparse(m_cholmodFactor, &m_cholmod);
00279 const_cast<typename Base::CholMatrixType&>(this->m_matrix) =
00280 internal::map_cholmod_sparse_to_eigen<Scalar,ColMajor,Index>(*cmRes);
00281 free(cmRes);
00282
00283 this->m_status = (this->m_status & ~Base::MatrixLIsDirty);
00284 }
00285 return this->m_matrix;
00286 }
00287
00288
00289
00290
00291 template<typename _MatrixType>
00292 template<typename Derived>
00293 bool SparseLLT<_MatrixType,Cholmod>::solveInPlace(MatrixBase<Derived> &b) const
00294 {
00295
00296 eigen_assert((Index)m_cholmodFactor->n==b.rows());
00297
00298
00299
00300
00301
00302
00303
00304 cholmod_dense cdb = internal::cholmod_map_eigen_to_dense(b);
00305
00306 cholmod_dense* x = cholmod_solve(CHOLMOD_A, m_cholmodFactor, &cdb, &m_cholmod);
00307 eigen_assert(x && "Eigen: cholmod_solve failed.");
00308
00309 b = Matrix<typename Base::Scalar,Dynamic,1>::Map(reinterpret_cast<typename Base::Scalar*>(x->x),b.rows());
00310 cholmod_free_dense(&x, &m_cholmod);
00311 return true;
00312 }
00313
00314
00315
00316
00317
00318
00319
00320
00321
00322
00323
00324 template<typename _MatrixType>
00325 class SparseLDLT<_MatrixType,Cholmod> : public SparseLDLT<_MatrixType>
00326 {
00327 protected:
00328 typedef SparseLDLT<_MatrixType> Base;
00329 typedef typename Base::Scalar Scalar;
00330 typedef typename Base::RealScalar RealScalar;
00331 using Base::MatrixLIsDirty;
00332 using Base::SupernodalFactorIsDirty;
00333 using Base::m_flags;
00334 using Base::m_matrix;
00335 using Base::m_status;
00336
00337 public:
00338 typedef _MatrixType MatrixType;
00339 typedef typename MatrixType::Index Index;
00340
00341 SparseLDLT(int flags = 0)
00342 : Base(flags), m_cholmodFactor(0)
00343 {
00344 cholmod_start(&m_cholmod);
00345 }
00346
00347 SparseLDLT(const _MatrixType& matrix, int flags = 0)
00348 : Base(flags), m_cholmodFactor(0)
00349 {
00350 cholmod_start(&m_cholmod);
00351 compute(matrix);
00352 }
00353
00354 ~SparseLDLT()
00355 {
00356 if (m_cholmodFactor)
00357 cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
00358 cholmod_finish(&m_cholmod);
00359 }
00360
00361 inline const typename Base::CholMatrixType& matrixL(void) const;
00362
00363 template<typename Derived>
00364 void solveInPlace(MatrixBase<Derived> &b) const;
00365
00366 template<typename Rhs>
00367 inline const internal::solve_retval<SparseLDLT<MatrixType, Cholmod>, Rhs>
00368 solve(const MatrixBase<Rhs>& b) const
00369 {
00370 eigen_assert(true && "SparseLDLT is not initialized.");
00371 return internal::solve_retval<SparseLDLT<MatrixType, Cholmod>, Rhs>(*this, b.derived());
00372 }
00373
00374 void compute(const _MatrixType& matrix);
00375
00376 inline Index cols() const { return m_matrix.cols(); }
00377 inline Index rows() const { return m_matrix.rows(); }
00378
00379 inline const cholmod_factor* cholmodFactor() const
00380 { return m_cholmodFactor; }
00381
00382 inline cholmod_common* cholmodCommon() const
00383 { return &m_cholmod; }
00384
00385 bool succeeded() const;
00386
00387 protected:
00388 mutable cholmod_common m_cholmod;
00389 cholmod_factor* m_cholmodFactor;
00390 };
00391
00392
00393
00394 namespace internal {
00395
00396 template<typename _MatrixType, typename Rhs>
00397 struct solve_retval<SparseLDLT<_MatrixType, Cholmod>, Rhs>
00398 : solve_retval_base<SparseLDLT<_MatrixType, Cholmod>, Rhs>
00399 {
00400 typedef SparseLDLT<_MatrixType, Cholmod> SpLDLTDecType;
00401 EIGEN_MAKE_SOLVE_HELPERS(SpLDLTDecType,Rhs)
00402
00403 template<typename Dest> void evalTo(Dest& dst) const
00404 {
00405
00406 eigen_assert((Index)dec().cholmodFactor()->n==rhs().rows());
00407
00408 cholmod_factor* cholmodFactor = const_cast<cholmod_factor*>(dec().cholmodFactor());
00409 cholmod_common* cholmodCommon = const_cast<cholmod_common*>(dec().cholmodCommon());
00410
00411
00412
00413
00414
00415
00416 cholmod_dense cdb = internal::cholmod_map_eigen_to_dense(rhs().const_cast_derived());
00417 cholmod_dense* x = cholmod_solve(CHOLMOD_LDLt, cholmodFactor, &cdb, cholmodCommon);
00418
00419 dst = Matrix<typename Base::Scalar,Dynamic,1>::Map(reinterpret_cast<typename Base::Scalar*>(x->x), rhs().rows());
00420 cholmod_free_dense(&x, cholmodCommon);
00421
00422 }
00423
00424 };
00425
00426
00427 }
00428
00429 template<typename _MatrixType>
00430 void SparseLDLT<_MatrixType,Cholmod>::compute(const _MatrixType& a)
00431 {
00432 if (m_cholmodFactor)
00433 {
00434 cholmod_free_factor(&m_cholmodFactor, &m_cholmod);
00435 m_cholmodFactor = 0;
00436 }
00437
00438 cholmod_sparse A = internal::cholmod_map_eigen_to_sparse(const_cast<_MatrixType&>(a));
00439
00440
00441 m_cholmod.supernodal = CHOLMOD_SIMPLICIAL;
00442
00443
00444 if (this->m_flags & IncompleteFactorization)
00445 {
00446 m_cholmod.nmethods = 1;
00447
00448 m_cholmod.method[0].ordering = CHOLMOD_COLAMD;
00449 m_cholmod.postorder = 1;
00450 }
00451 else
00452 {
00453 m_cholmod.nmethods = 1;
00454 m_cholmod.method[0].ordering = CHOLMOD_NATURAL;
00455 m_cholmod.postorder = 0;
00456 }
00457 m_cholmod.final_ll = 0;
00458 m_cholmodFactor = cholmod_analyze(&A, &m_cholmod);
00459 cholmod_factorize(&A, m_cholmodFactor, &m_cholmod);
00460
00461 this->m_status = (this->m_status & ~Base::SupernodalFactorIsDirty) | Base::MatrixLIsDirty;
00462 }
00463
00464
00465
00466 template<typename _MatrixType>
00467 bool SparseLDLT<_MatrixType,Cholmod>::succeeded() const
00468 { return true; }
00469
00470
00471 template<typename _MatrixType>
00472 inline const typename SparseLDLT<_MatrixType>::CholMatrixType&
00473 SparseLDLT<_MatrixType,Cholmod>::matrixL() const
00474 {
00475 if (this->m_status & Base::MatrixLIsDirty)
00476 {
00477 eigen_assert(!(this->m_status & Base::SupernodalFactorIsDirty));
00478
00479 cholmod_sparse* cmRes = cholmod_factor_to_sparse(m_cholmodFactor, &m_cholmod);
00480 const_cast<typename Base::CholMatrixType&>(this->m_matrix) = MappedSparseMatrix<Scalar>(*cmRes);
00481 free(cmRes);
00482
00483 this->m_status = (this->m_status & ~Base::MatrixLIsDirty);
00484 }
00485 return this->m_matrix;
00486 }
00487
00488
00489
00490
00491
00492
00493 template<typename _MatrixType>
00494 template<typename Derived>
00495 void SparseLDLT<_MatrixType,Cholmod>::solveInPlace(MatrixBase<Derived> &b) const
00496 {
00497
00498 eigen_assert((Index)m_cholmodFactor->n == b.rows());
00499
00500
00501
00502
00503
00504
00505
00506 cholmod_dense cdb = internal::cholmod_map_eigen_to_dense(b);
00507 cholmod_dense* x = cholmod_solve(CHOLMOD_A, m_cholmodFactor, &cdb, &m_cholmod);
00508 b = Matrix<typename Base::Scalar,Dynamic,1>::Map(reinterpret_cast<typename Base::Scalar*>(x->x),b.rows());
00509 cholmod_free_dense(&x, &m_cholmod);
00510 }
00511
00512
00513
00514
00515
00516
00517 #endif // EIGEN_CHOLMODSUPPORT_LEGACY_H