UmfPackSupport.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 //
00006 // This Source Code Form is subject to the terms of the Mozilla
00007 // Public License v. 2.0. If a copy of the MPL was not distributed
00008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00009 
00010 #ifndef EIGEN_UMFPACKSUPPORT_H
00011 #define EIGEN_UMFPACKSUPPORT_H
00012 
00013 namespace Eigen { 
00014 
00015 /* TODO extract L, extract U, compute det, etc... */
00016 
00017 // generic double/complex<double> wrapper functions:
00018 
00019 inline void umfpack_free_numeric(void **Numeric, double)
00020 { umfpack_di_free_numeric(Numeric); *Numeric = 0; }
00021 
00022 inline void umfpack_free_numeric(void **Numeric, std::complex<double>)
00023 { umfpack_zi_free_numeric(Numeric); *Numeric = 0; }
00024 
00025 inline void umfpack_free_symbolic(void **Symbolic, double)
00026 { umfpack_di_free_symbolic(Symbolic); *Symbolic = 0; }
00027 
00028 inline void umfpack_free_symbolic(void **Symbolic, std::complex<double>)
00029 { umfpack_zi_free_symbolic(Symbolic); *Symbolic = 0; }
00030 
00031 inline int umfpack_symbolic(int n_row,int n_col,
00032                             const int Ap[], const int Ai[], const double Ax[], void **Symbolic,
00033                             const double Control [UMFPACK_CONTROL], double Info [UMFPACK_INFO])
00034 {
00035   return umfpack_di_symbolic(n_row,n_col,Ap,Ai,Ax,Symbolic,Control,Info);
00036 }
00037 
00038 inline int umfpack_symbolic(int n_row,int n_col,
00039                             const int Ap[], const int Ai[], const std::complex<double> Ax[], void **Symbolic,
00040                             const double Control [UMFPACK_CONTROL], double Info [UMFPACK_INFO])
00041 {
00042   return umfpack_zi_symbolic(n_row,n_col,Ap,Ai,&internal::real_ref(Ax[0]),0,Symbolic,Control,Info);
00043 }
00044 
00045 inline int umfpack_numeric( const int Ap[], const int Ai[], const double Ax[],
00046                             void *Symbolic, void **Numeric,
00047                             const double Control[UMFPACK_CONTROL],double Info [UMFPACK_INFO])
00048 {
00049   return umfpack_di_numeric(Ap,Ai,Ax,Symbolic,Numeric,Control,Info);
00050 }
00051 
00052 inline int umfpack_numeric( const int Ap[], const int Ai[], const std::complex<double> Ax[],
00053                             void *Symbolic, void **Numeric,
00054                             const double Control[UMFPACK_CONTROL],double Info [UMFPACK_INFO])
00055 {
00056   return umfpack_zi_numeric(Ap,Ai,&internal::real_ref(Ax[0]),0,Symbolic,Numeric,Control,Info);
00057 }
00058 
00059 inline int umfpack_solve( int sys, const int Ap[], const int Ai[], const double Ax[],
00060                           double X[], const double B[], void *Numeric,
00061                           const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
00062 {
00063   return umfpack_di_solve(sys,Ap,Ai,Ax,X,B,Numeric,Control,Info);
00064 }
00065 
00066 inline int umfpack_solve( int sys, const int Ap[], const int Ai[], const std::complex<double> Ax[],
00067                           std::complex<double> X[], const std::complex<double> B[], void *Numeric,
00068                           const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
00069 {
00070   return umfpack_zi_solve(sys,Ap,Ai,&internal::real_ref(Ax[0]),0,&internal::real_ref(X[0]),0,&internal::real_ref(B[0]),0,Numeric,Control,Info);
00071 }
00072 
00073 inline int umfpack_get_lunz(int *lnz, int *unz, int *n_row, int *n_col, int *nz_udiag, void *Numeric, double)
00074 {
00075   return umfpack_di_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
00076 }
00077 
00078 inline int umfpack_get_lunz(int *lnz, int *unz, int *n_row, int *n_col, int *nz_udiag, void *Numeric, std::complex<double>)
00079 {
00080   return umfpack_zi_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
00081 }
00082 
00083 inline int umfpack_get_numeric(int Lp[], int Lj[], double Lx[], int Up[], int Ui[], double Ux[],
00084                                int P[], int Q[], double Dx[], int *do_recip, double Rs[], void *Numeric)
00085 {
00086   return umfpack_di_get_numeric(Lp,Lj,Lx,Up,Ui,Ux,P,Q,Dx,do_recip,Rs,Numeric);
00087 }
00088 
00089 inline int umfpack_get_numeric(int Lp[], int Lj[], std::complex<double> Lx[], int Up[], int Ui[], std::complex<double> Ux[],
00090                                int P[], int Q[], std::complex<double> Dx[], int *do_recip, double Rs[], void *Numeric)
00091 {
00092   double& lx0_real = internal::real_ref(Lx[0]);
00093   double& ux0_real = internal::real_ref(Ux[0]);
00094   double& dx0_real = internal::real_ref(Dx[0]);
00095   return umfpack_zi_get_numeric(Lp,Lj,Lx?&lx0_real:0,0,Up,Ui,Ux?&ux0_real:0,0,P,Q,
00096                                 Dx?&dx0_real:0,0,do_recip,Rs,Numeric);
00097 }
00098 
00099 inline int umfpack_get_determinant(double *Mx, double *Ex, void *NumericHandle, double User_Info [UMFPACK_INFO])
00100 {
00101   return umfpack_di_get_determinant(Mx,Ex,NumericHandle,User_Info);
00102 }
00103 
00104 inline int umfpack_get_determinant(std::complex<double> *Mx, double *Ex, void *NumericHandle, double User_Info [UMFPACK_INFO])
00105 {
00106   double& mx_real = internal::real_ref(*Mx);
00107   return umfpack_zi_get_determinant(&mx_real,0,Ex,NumericHandle,User_Info);
00108 }
00109 
00123 template<typename _MatrixType>
00124 class UmfPackLU : internal::noncopyable
00125 {
00126   public:
00127     typedef _MatrixType MatrixType;
00128     typedef typename MatrixType::Scalar Scalar;
00129     typedef typename MatrixType::RealScalar RealScalar;
00130     typedef typename MatrixType::Index Index;
00131     typedef Matrix<Scalar,Dynamic,1> Vector;
00132     typedef Matrix<int, 1, MatrixType::ColsAtCompileTime> IntRowVectorType;
00133     typedef Matrix<int, MatrixType::RowsAtCompileTime, 1> IntColVectorType;
00134     typedef SparseMatrix<Scalar> LUMatrixType;
00135     typedef SparseMatrix<Scalar,ColMajor,int> UmfpackMatrixType;
00136 
00137   public:
00138 
00139     UmfPackLU() { init(); }
00140 
00141     UmfPackLU(const MatrixType& matrix)
00142     {
00143       init();
00144       compute(matrix);
00145     }
00146 
00147     ~UmfPackLU()
00148     {
00149       if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
00150       if(m_numeric)  umfpack_free_numeric(&m_numeric,Scalar());
00151     }
00152 
00153     inline Index rows() const { return m_copyMatrix.rows(); }
00154     inline Index cols() const { return m_copyMatrix.cols(); }
00155 
00161     ComputationInfo info() const
00162     {
00163       eigen_assert(m_isInitialized && "Decomposition is not initialized.");
00164       return m_info;
00165     }
00166 
00167     inline const LUMatrixType& matrixL() const
00168     {
00169       if (m_extractedDataAreDirty) extractData();
00170       return m_l;
00171     }
00172 
00173     inline const LUMatrixType& matrixU() const
00174     {
00175       if (m_extractedDataAreDirty) extractData();
00176       return m_u;
00177     }
00178 
00179     inline const IntColVectorType& permutationP() const
00180     {
00181       if (m_extractedDataAreDirty) extractData();
00182       return m_p;
00183     }
00184 
00185     inline const IntRowVectorType& permutationQ() const
00186     {
00187       if (m_extractedDataAreDirty) extractData();
00188       return m_q;
00189     }
00190 
00195     void compute(const MatrixType& matrix)
00196     {
00197       analyzePattern(matrix);
00198       factorize(matrix);
00199     }
00200 
00205     template<typename Rhs>
00206     inline const internal::solve_retval<UmfPackLU, Rhs> solve(const MatrixBase<Rhs>& b) const
00207     {
00208       eigen_assert(m_isInitialized && "UmfPackLU is not initialized.");
00209       eigen_assert(rows()==b.rows()
00210                 && "UmfPackLU::solve(): invalid number of rows of the right hand side matrix b");
00211       return internal::solve_retval<UmfPackLU, Rhs>(*this, b.derived());
00212     }
00213 
00218 //     template<typename Rhs>
00219 //     inline const internal::sparse_solve_retval<UmfPAckLU, Rhs> solve(const SparseMatrixBase<Rhs>& b) const
00220 //     {
00221 //       eigen_assert(m_isInitialized && "UmfPAckLU is not initialized.");
00222 //       eigen_assert(rows()==b.rows()
00223 //                 && "UmfPAckLU::solve(): invalid number of rows of the right hand side matrix b");
00224 //       return internal::sparse_solve_retval<UmfPAckLU, Rhs>(*this, b.derived());
00225 //     }
00226 
00233     void analyzePattern(const MatrixType& matrix)
00234     {
00235       if(m_symbolic)
00236         umfpack_free_symbolic(&m_symbolic,Scalar());
00237       if(m_numeric)
00238         umfpack_free_numeric(&m_numeric,Scalar());
00239       
00240       grapInput(matrix);
00241 
00242       int errorCode = 0;
00243       errorCode = umfpack_symbolic(matrix.rows(), matrix.cols(), m_outerIndexPtr, m_innerIndexPtr, m_valuePtr,
00244                                    &m_symbolic, 0, 0);
00245 
00246       m_isInitialized = true;
00247       m_info = errorCode ? InvalidInput : Success;
00248       m_analysisIsOk = true;
00249       m_factorizationIsOk = false;
00250     }
00251 
00258     void factorize(const MatrixType& matrix)
00259     {
00260       eigen_assert(m_analysisIsOk && "UmfPackLU: you must first call analyzePattern()");
00261       if(m_numeric)
00262         umfpack_free_numeric(&m_numeric,Scalar());
00263 
00264       grapInput(matrix);
00265 
00266       int errorCode;
00267       errorCode = umfpack_numeric(m_outerIndexPtr, m_innerIndexPtr, m_valuePtr,
00268                                   m_symbolic, &m_numeric, 0, 0);
00269 
00270       m_info = errorCode ? NumericalIssue : Success;
00271       m_factorizationIsOk = true;
00272     }
00273 
00274     #ifndef EIGEN_PARSED_BY_DOXYGEN
00275 
00276     template<typename BDerived,typename XDerived>
00277     bool _solve(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const;
00278     #endif
00279 
00280     Scalar determinant() const;
00281 
00282     void extractData() const;
00283 
00284   protected:
00285 
00286 
00287     void init()
00288     {
00289       m_info = InvalidInput;
00290       m_isInitialized = false;
00291       m_numeric = 0;
00292       m_symbolic = 0;
00293       m_outerIndexPtr = 0;
00294       m_innerIndexPtr = 0;
00295       m_valuePtr      = 0;
00296     }
00297     
00298     void grapInput(const MatrixType& mat)
00299     {
00300       m_copyMatrix.resize(mat.rows(), mat.cols());
00301       if( ((MatrixType::Flags&RowMajorBit)==RowMajorBit) || sizeof(typename MatrixType::Index)!=sizeof(int) || !mat.isCompressed() )
00302       {
00303         // non supported input -> copy
00304         m_copyMatrix = mat;
00305         m_outerIndexPtr = m_copyMatrix.outerIndexPtr();
00306         m_innerIndexPtr = m_copyMatrix.innerIndexPtr();
00307         m_valuePtr      = m_copyMatrix.valuePtr();
00308       }
00309       else
00310       {
00311         m_outerIndexPtr = mat.outerIndexPtr();
00312         m_innerIndexPtr = mat.innerIndexPtr();
00313         m_valuePtr      = mat.valuePtr();
00314       }
00315     }
00316 
00317     // cached data to reduce reallocation, etc.
00318     mutable LUMatrixType m_l;
00319     mutable LUMatrixType m_u;
00320     mutable IntColVectorType m_p;
00321     mutable IntRowVectorType m_q;
00322 
00323     UmfpackMatrixType m_copyMatrix;
00324     const Scalar* m_valuePtr;
00325     const int* m_outerIndexPtr;
00326     const int* m_innerIndexPtr;
00327     void* m_numeric;
00328     void* m_symbolic;
00329 
00330     mutable ComputationInfo m_info;
00331     bool m_isInitialized;
00332     int m_factorizationIsOk;
00333     int m_analysisIsOk;
00334     mutable bool m_extractedDataAreDirty;
00335     
00336   private:
00337     UmfPackLU(UmfPackLU& ) { }
00338 };
00339 
00340 
00341 template<typename MatrixType>
00342 void UmfPackLU<MatrixType>::extractData() const
00343 {
00344   if (m_extractedDataAreDirty)
00345   {
00346     // get size of the data
00347     int lnz, unz, rows, cols, nz_udiag;
00348     umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());
00349 
00350     // allocate data
00351     m_l.resize(rows,(std::min)(rows,cols));
00352     m_l.resizeNonZeros(lnz);
00353 
00354     m_u.resize((std::min)(rows,cols),cols);
00355     m_u.resizeNonZeros(unz);
00356 
00357     m_p.resize(rows);
00358     m_q.resize(cols);
00359 
00360     // extract
00361     umfpack_get_numeric(m_l.outerIndexPtr(), m_l.innerIndexPtr(), m_l.valuePtr(),
00362                         m_u.outerIndexPtr(), m_u.innerIndexPtr(), m_u.valuePtr(),
00363                         m_p.data(), m_q.data(), 0, 0, 0, m_numeric);
00364 
00365     m_extractedDataAreDirty = false;
00366   }
00367 }
00368 
00369 template<typename MatrixType>
00370 typename UmfPackLU<MatrixType>::Scalar UmfPackLU<MatrixType>::determinant() const
00371 {
00372   Scalar det;
00373   umfpack_get_determinant(&det, 0, m_numeric, 0);
00374   return det;
00375 }
00376 
00377 template<typename MatrixType>
00378 template<typename BDerived,typename XDerived>
00379 bool UmfPackLU<MatrixType>::_solve(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const
00380 {
00381   const int rhsCols = b.cols();
00382   eigen_assert((BDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major rhs yet");
00383   eigen_assert((XDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major result yet");
00384 
00385   int errorCode;
00386   for (int j=0; j<rhsCols; ++j)
00387   {
00388     errorCode = umfpack_solve(UMFPACK_A,
00389         m_outerIndexPtr, m_innerIndexPtr, m_valuePtr,
00390         &x.col(j).coeffRef(0), &b.const_cast_derived().col(j).coeffRef(0), m_numeric, 0, 0);
00391     if (errorCode!=0)
00392       return false;
00393   }
00394 
00395   return true;
00396 }
00397 
00398 
00399 namespace internal {
00400 
00401 template<typename _MatrixType, typename Rhs>
00402 struct solve_retval<UmfPackLU<_MatrixType>, Rhs>
00403   : solve_retval_base<UmfPackLU<_MatrixType>, Rhs>
00404 {
00405   typedef UmfPackLU<_MatrixType> Dec;
00406   EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
00407 
00408   template<typename Dest> void evalTo(Dest& dst) const
00409   {
00410     dec()._solve(rhs(),dst);
00411   }
00412 };
00413 
00414 template<typename _MatrixType, typename Rhs>
00415 struct sparse_solve_retval<UmfPackLU<_MatrixType>, Rhs>
00416   : sparse_solve_retval_base<UmfPackLU<_MatrixType>, Rhs>
00417 {
00418   typedef UmfPackLU<_MatrixType> Dec;
00419   EIGEN_MAKE_SPARSE_SOLVE_HELPERS(Dec,Rhs)
00420 
00421   template<typename Dest> void evalTo(Dest& dst) const
00422   {
00423     dec()._solve(rhs(),dst);
00424   }
00425 };
00426 
00427 } // end namespace internal
00428 
00429 } // end namespace Eigen
00430 
00431 #endif // EIGEN_UMFPACKSUPPORT_H


win_eigen
Author(s): Daniel Stonier
autogenerated on Wed Sep 16 2015 07:12:36