Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-04-04 08:48:40

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2017 Kyle Macfarlan <kyle.macfarlan@gmail.com>
0005 //
0006 // This Source Code Form is subject to the terms of the Mozilla
0007 // Public License v. 2.0. If a copy of the MPL was not distributed
0008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0009 
0010 #ifndef EIGEN_KLUSUPPORT_H
0011 #define EIGEN_KLUSUPPORT_H
0012 
0013 namespace Eigen {
0014 
0015 /* TODO extract L, extract U, compute det, etc... */
0016 
0017 /** \ingroup KLUSupport_Module
0018   * \brief A sparse LU factorization and solver based on KLU
0019   *
0020   * This class allows to solve for A.X = B sparse linear problems via a LU factorization
0021   * using the KLU library. The sparse matrix A must be squared and full rank.
0022   * The vectors or matrices X and B can be either dense or sparse.
0023   *
0024   * \warning The input matrix A should be in a \b compressed and \b column-major form.
0025   * Otherwise an expensive copy will be made. You can call the inexpensive makeCompressed() to get a compressed matrix.
0026   * \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
0027   *
0028   * \implsparsesolverconcept
0029   *
0030   * \sa \ref TutorialSparseSolverConcept, class UmfPackLU, class SparseLU
0031   */
0032 
0033 
0034 inline int klu_solve(klu_symbolic *Symbolic, klu_numeric *Numeric, Index ldim, Index nrhs, double B [ ], klu_common *Common, double) {
0035    return klu_solve(Symbolic, Numeric, internal::convert_index<int>(ldim), internal::convert_index<int>(nrhs), B, Common);
0036 }
0037 
0038 inline int klu_solve(klu_symbolic *Symbolic, klu_numeric *Numeric, Index ldim, Index nrhs, std::complex<double>B[], klu_common *Common, std::complex<double>) {
0039    return klu_z_solve(Symbolic, Numeric, internal::convert_index<int>(ldim), internal::convert_index<int>(nrhs), &numext::real_ref(B[0]), Common);
0040 }
0041 
0042 inline int klu_tsolve(klu_symbolic *Symbolic, klu_numeric *Numeric, Index ldim, Index nrhs, double B[], klu_common *Common, double) {
0043    return klu_tsolve(Symbolic, Numeric, internal::convert_index<int>(ldim), internal::convert_index<int>(nrhs), B, Common);
0044 }
0045 
0046 inline int klu_tsolve(klu_symbolic *Symbolic, klu_numeric *Numeric, Index ldim, Index nrhs, std::complex<double>B[], klu_common *Common, std::complex<double>) {
0047    return klu_z_tsolve(Symbolic, Numeric, internal::convert_index<int>(ldim), internal::convert_index<int>(nrhs), &numext::real_ref(B[0]), 0, Common);
0048 }
0049 
0050 inline klu_numeric* klu_factor(int Ap [ ], int Ai [ ], double Ax [ ], klu_symbolic *Symbolic, klu_common *Common, double) {
0051    return klu_factor(Ap, Ai, Ax, Symbolic, Common);
0052 }
0053 
0054 inline klu_numeric* klu_factor(int Ap[], int Ai[], std::complex<double> Ax[], klu_symbolic *Symbolic, klu_common *Common, std::complex<double>) {
0055    return klu_z_factor(Ap, Ai, &numext::real_ref(Ax[0]), Symbolic, Common);
0056 }
0057 
0058 
0059 template<typename _MatrixType>
0060 class KLU : public SparseSolverBase<KLU<_MatrixType> >
0061 {
0062   protected:
0063     typedef SparseSolverBase<KLU<_MatrixType> > Base;
0064     using Base::m_isInitialized;
0065   public:
0066     using Base::_solve_impl;
0067     typedef _MatrixType MatrixType;
0068     typedef typename MatrixType::Scalar Scalar;
0069     typedef typename MatrixType::RealScalar RealScalar;
0070     typedef typename MatrixType::StorageIndex StorageIndex;
0071     typedef Matrix<Scalar,Dynamic,1> Vector;
0072     typedef Matrix<int, 1, MatrixType::ColsAtCompileTime> IntRowVectorType;
0073     typedef Matrix<int, MatrixType::RowsAtCompileTime, 1> IntColVectorType;
0074     typedef SparseMatrix<Scalar> LUMatrixType;
0075     typedef SparseMatrix<Scalar,ColMajor,int> KLUMatrixType;
0076     typedef Ref<const KLUMatrixType, StandardCompressedFormat> KLUMatrixRef;
0077     enum {
0078       ColsAtCompileTime = MatrixType::ColsAtCompileTime,
0079       MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
0080     };
0081 
0082   public:
0083 
0084     KLU()
0085       : m_dummy(0,0), mp_matrix(m_dummy)
0086     {
0087       init();
0088     }
0089 
0090     template<typename InputMatrixType>
0091     explicit KLU(const InputMatrixType& matrix)
0092       : mp_matrix(matrix)
0093     {
0094       init();
0095       compute(matrix);
0096     }
0097 
0098     ~KLU()
0099     {
0100       if(m_symbolic) klu_free_symbolic(&m_symbolic,&m_common);
0101       if(m_numeric)  klu_free_numeric(&m_numeric,&m_common);
0102     }
0103 
0104     EIGEN_CONSTEXPR inline Index rows() const EIGEN_NOEXCEPT { return mp_matrix.rows(); }
0105     EIGEN_CONSTEXPR inline Index cols() const EIGEN_NOEXCEPT { return mp_matrix.cols(); }
0106 
0107     /** \brief Reports whether previous computation was successful.
0108       *
0109       * \returns \c Success if computation was successful,
0110       *          \c NumericalIssue if the matrix.appears to be negative.
0111       */
0112     ComputationInfo info() const
0113     {
0114       eigen_assert(m_isInitialized && "Decomposition is not initialized.");
0115       return m_info;
0116     }
0117 #if 0 // not implemented yet
0118     inline const LUMatrixType& matrixL() const
0119     {
0120       if (m_extractedDataAreDirty) extractData();
0121       return m_l;
0122     }
0123 
0124     inline const LUMatrixType& matrixU() const
0125     {
0126       if (m_extractedDataAreDirty) extractData();
0127       return m_u;
0128     }
0129 
0130     inline const IntColVectorType& permutationP() const
0131     {
0132       if (m_extractedDataAreDirty) extractData();
0133       return m_p;
0134     }
0135 
0136     inline const IntRowVectorType& permutationQ() const
0137     {
0138       if (m_extractedDataAreDirty) extractData();
0139       return m_q;
0140     }
0141 #endif
0142     /** Computes the sparse Cholesky decomposition of \a matrix
0143      *  Note that the matrix should be column-major, and in compressed format for best performance.
0144      *  \sa SparseMatrix::makeCompressed().
0145      */
0146     template<typename InputMatrixType>
0147     void compute(const InputMatrixType& matrix)
0148     {
0149       if(m_symbolic) klu_free_symbolic(&m_symbolic, &m_common);
0150       if(m_numeric)  klu_free_numeric(&m_numeric, &m_common);
0151       grab(matrix.derived());
0152       analyzePattern_impl();
0153       factorize_impl();
0154     }
0155 
0156     /** Performs a symbolic decomposition on the sparcity of \a matrix.
0157       *
0158       * This function is particularly useful when solving for several problems having the same structure.
0159       *
0160       * \sa factorize(), compute()
0161       */
0162     template<typename InputMatrixType>
0163     void analyzePattern(const InputMatrixType& matrix)
0164     {
0165       if(m_symbolic) klu_free_symbolic(&m_symbolic, &m_common);
0166       if(m_numeric)  klu_free_numeric(&m_numeric, &m_common);
0167 
0168       grab(matrix.derived());
0169 
0170       analyzePattern_impl();
0171     }
0172 
0173 
0174     /** Provides access to the control settings array used by KLU.
0175       *
0176       * See KLU documentation for details.
0177       */
0178     inline const klu_common& kluCommon() const
0179     {
0180       return m_common;
0181     }
0182 
0183     /** Provides access to the control settings array used by UmfPack.
0184       *
0185       * If this array contains NaN's, the default values are used.
0186       *
0187       * See KLU documentation for details.
0188       */
0189     inline klu_common& kluCommon()
0190     {
0191       return m_common;
0192     }
0193 
0194     /** Performs a numeric decomposition of \a matrix
0195       *
0196       * The given matrix must has the same sparcity than the matrix on which the pattern anylysis has been performed.
0197       *
0198       * \sa analyzePattern(), compute()
0199       */
0200     template<typename InputMatrixType>
0201     void factorize(const InputMatrixType& matrix)
0202     {
0203       eigen_assert(m_analysisIsOk && "KLU: you must first call analyzePattern()");
0204       if(m_numeric)
0205         klu_free_numeric(&m_numeric,&m_common);
0206 
0207       grab(matrix.derived());
0208 
0209       factorize_impl();
0210     }
0211 
0212     /** \internal */
0213     template<typename BDerived,typename XDerived>
0214     bool _solve_impl(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const;
0215 
0216 #if 0 // not implemented yet
0217     Scalar determinant() const;
0218 
0219     void extractData() const;
0220 #endif
0221 
0222   protected:
0223 
0224     void init()
0225     {
0226       m_info                  = InvalidInput;
0227       m_isInitialized         = false;
0228       m_numeric               = 0;
0229       m_symbolic              = 0;
0230       m_extractedDataAreDirty = true;
0231 
0232       klu_defaults(&m_common);
0233     }
0234 
0235     void analyzePattern_impl()
0236     {
0237       m_info = InvalidInput;
0238       m_analysisIsOk = false;
0239       m_factorizationIsOk = false;
0240       m_symbolic = klu_analyze(internal::convert_index<int>(mp_matrix.rows()),
0241                                      const_cast<StorageIndex*>(mp_matrix.outerIndexPtr()), const_cast<StorageIndex*>(mp_matrix.innerIndexPtr()),
0242                                      &m_common);
0243       if (m_symbolic) {
0244          m_isInitialized = true;
0245          m_info = Success;
0246          m_analysisIsOk = true;
0247          m_extractedDataAreDirty = true;
0248       }
0249     }
0250 
0251     void factorize_impl()
0252     {
0253 
0254       m_numeric = klu_factor(const_cast<StorageIndex*>(mp_matrix.outerIndexPtr()), const_cast<StorageIndex*>(mp_matrix.innerIndexPtr()), const_cast<Scalar*>(mp_matrix.valuePtr()),
0255                                     m_symbolic, &m_common, Scalar());
0256 
0257 
0258       m_info = m_numeric ? Success : NumericalIssue;
0259       m_factorizationIsOk = m_numeric ? 1 : 0;
0260       m_extractedDataAreDirty = true;
0261     }
0262 
0263     template<typename MatrixDerived>
0264     void grab(const EigenBase<MatrixDerived> &A)
0265     {
0266       mp_matrix.~KLUMatrixRef();
0267       ::new (&mp_matrix) KLUMatrixRef(A.derived());
0268     }
0269 
0270     void grab(const KLUMatrixRef &A)
0271     {
0272       if(&(A.derived()) != &mp_matrix)
0273       {
0274         mp_matrix.~KLUMatrixRef();
0275         ::new (&mp_matrix) KLUMatrixRef(A);
0276       }
0277     }
0278 
0279     // cached data to reduce reallocation, etc.
0280 #if 0 // not implemented yet
0281     mutable LUMatrixType m_l;
0282     mutable LUMatrixType m_u;
0283     mutable IntColVectorType m_p;
0284     mutable IntRowVectorType m_q;
0285 #endif
0286 
0287     KLUMatrixType m_dummy;
0288     KLUMatrixRef mp_matrix;
0289 
0290     klu_numeric* m_numeric;
0291     klu_symbolic* m_symbolic;
0292     klu_common m_common;
0293     mutable ComputationInfo m_info;
0294     int m_factorizationIsOk;
0295     int m_analysisIsOk;
0296     mutable bool m_extractedDataAreDirty;
0297 
0298   private:
0299     KLU(const KLU& ) { }
0300 };
0301 
0302 #if 0 // not implemented yet
0303 template<typename MatrixType>
0304 void KLU<MatrixType>::extractData() const
0305 {
0306   if (m_extractedDataAreDirty)
0307   {
0308      eigen_assert(false && "KLU: extractData Not Yet Implemented");
0309 
0310     // get size of the data
0311     int lnz, unz, rows, cols, nz_udiag;
0312     umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());
0313 
0314     // allocate data
0315     m_l.resize(rows,(std::min)(rows,cols));
0316     m_l.resizeNonZeros(lnz);
0317 
0318     m_u.resize((std::min)(rows,cols),cols);
0319     m_u.resizeNonZeros(unz);
0320 
0321     m_p.resize(rows);
0322     m_q.resize(cols);
0323 
0324     // extract
0325     umfpack_get_numeric(m_l.outerIndexPtr(), m_l.innerIndexPtr(), m_l.valuePtr(),
0326                         m_u.outerIndexPtr(), m_u.innerIndexPtr(), m_u.valuePtr(),
0327                         m_p.data(), m_q.data(), 0, 0, 0, m_numeric);
0328 
0329     m_extractedDataAreDirty = false;
0330   }
0331 }
0332 
0333 template<typename MatrixType>
0334 typename KLU<MatrixType>::Scalar KLU<MatrixType>::determinant() const
0335 {
0336   eigen_assert(false && "KLU: extractData Not Yet Implemented");
0337   return Scalar();
0338 }
0339 #endif
0340 
0341 template<typename MatrixType>
0342 template<typename BDerived,typename XDerived>
0343 bool KLU<MatrixType>::_solve_impl(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const
0344 {
0345   Index rhsCols = b.cols();
0346   EIGEN_STATIC_ASSERT((XDerived::Flags&RowMajorBit)==0, THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
0347   eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or analyzePattern()/factorize()");
0348 
0349   x = b;
0350   int info = klu_solve(m_symbolic, m_numeric, b.rows(), rhsCols, x.const_cast_derived().data(), const_cast<klu_common*>(&m_common), Scalar());
0351 
0352   m_info = info!=0 ? Success : NumericalIssue;
0353   return true;
0354 }
0355 
0356 } // end namespace Eigen
0357 
0358 #endif // EIGEN_KLUSUPPORT_H