// This file is part of Eigen, a lightweight C++ template library // for linear algebra. Eigen itself is part of the KDE project. // // Copyright (C) 2008 Gael Guennebaud // // Eigen is free software; you can redistribute it and/or // modify it under the terms of the GNU Lesser General Public // License as published by the Free Software Foundation; either // version 3 of the License, or (at your option) any later version. // // Alternatively, you can redistribute it and/or // modify it under the terms of the GNU General Public License as // published by the Free Software Foundation; either version 2 of // the License, or (at your option) any later version. // // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the // GNU General Public License for more details. // // You should have received a copy of the GNU Lesser General Public // License and a copy of the GNU General Public License along with // Eigen. If not, see . #ifndef EIGEN_CHOLMODSUPPORT_H #define EIGEN_CHOLMODSUPPORT_H template cholmod_sparse SparseMatrix::asCholmodMatrix() { cholmod_sparse res; res.nzmax = nonZeros(); res.nrow = rows();; res.ncol = cols(); res.p = _outerIndexPtr(); res.i = _innerIndexPtr(); res.x = _valuePtr(); res.xtype = CHOLMOD_REAL; res.itype = CHOLMOD_INT; res.sorted = 1; res.packed = 1; res.dtype = 0; res.stype = -1; if (ei_is_same_type::ret) { res.xtype = CHOLMOD_REAL; res.dtype = 1; } else if (ei_is_same_type::ret) { res.xtype = CHOLMOD_REAL; res.dtype = 0; } else if (ei_is_same_type >::ret) { res.xtype = CHOLMOD_COMPLEX; res.dtype = 1; } else if (ei_is_same_type >::ret) { res.xtype = CHOLMOD_COMPLEX; res.dtype = 0; } else { ei_assert(false && "Scalar type not supported by CHOLMOD"); } if (Flags & SelfAdjoint) { if (Flags & Upper) res.stype = 1; else if (Flags & Lower) res.stype = -1; else res.stype = 0; } else res.stype = 0; return res; } template SparseMatrix SparseMatrix::Map(cholmod_sparse& cm) { SparseMatrix res; res.m_innerSize = cm.nrow; res.m_outerSize = cm.ncol; res.m_outerIndex = reinterpret_cast(cm.p); SparseArray data = SparseArray::Map( reinterpret_cast(cm.i), reinterpret_cast(cm.x), res.m_outerIndex[cm.ncol]); res.m_data.swap(data); res.markAsRValue(); return res; } template class SparseCholesky : public SparseCholesky { protected: typedef SparseCholesky Base; using Base::Scalar; using Base::RealScalar; using Base::MatrixLIsDirty; using Base::SupernodalFactorIsDirty; using Base::m_flags; using Base::m_matrix; using Base::m_status; public: SparseCholesky(const MatrixType& matrix, int flags = 0) : Base(matrix, flags), m_cholmodFactor(0) { cholmod_start(&m_cholmod); compute(matrix); } ~SparseCholesky() { if (m_cholmodFactor) cholmod_free_factor(&m_cholmodFactor, &m_cholmod); cholmod_finish(&m_cholmod); } inline const typename Base::CholMatrixType& matrixL(void) const; template void solveInPlace(MatrixBase &b) const; void compute(const MatrixType& matrix); protected: mutable cholmod_common m_cholmod; cholmod_factor* m_cholmodFactor; }; template void SparseCholesky::compute(const MatrixType& a) { if (m_cholmodFactor) { cholmod_free_factor(&m_cholmodFactor, &m_cholmod); m_cholmodFactor = 0; } cholmod_sparse A = const_cast(a).asCholmodMatrix(); if (m_flags&IncompleteFactorization) { m_cholmod.nmethods = 1; m_cholmod.method [0].ordering = CHOLMOD_NATURAL; m_cholmod.postorder = 0; } else { m_cholmod.nmethods = 1; m_cholmod.method[0].ordering = CHOLMOD_NATURAL; m_cholmod.postorder = 0; } m_cholmod.final_ll = 1; m_cholmodFactor = cholmod_analyze(&A, &m_cholmod); cholmod_factorize(&A, m_cholmodFactor, &m_cholmod); m_status = (m_status & ~SupernodalFactorIsDirty) | MatrixLIsDirty; } template inline const typename SparseCholesky::CholMatrixType& SparseCholesky::matrixL() const { if (m_status & MatrixLIsDirty) { ei_assert(!(m_status & SupernodalFactorIsDirty)); cholmod_sparse* cmRes = cholmod_factor_to_sparse(m_cholmodFactor, &m_cholmod); const_cast(m_matrix) = Base::CholMatrixType::Map(*cmRes); free(cmRes); m_status = (m_status & ~MatrixLIsDirty); } return m_matrix; } template template void SparseCholesky::solveInPlace(MatrixBase &b) const { const int size = m_matrix.rows(); ei_assert(size==b.rows()); if (m_status & MatrixLIsDirty) { // ei_assert(!(m_status & SupernodalFactorIsDirty)); // taucs_supernodal_solve_llt(m_taucsSupernodalFactor,double* b); matrixL(); } // else { Base::solveInPlace(b); } } #endif // EIGEN_CHOLMODSUPPORT_H