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-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 //
00006 // Eigen is free software; you can redistribute it and/or
00007 // modify it under the terms of the GNU Lesser General Public
00008 // License as published by the Free Software Foundation; either
00009 // version 3 of the License, or (at your option) any later version.
00010 //
00011 // Alternatively, you can redistribute it and/or
00012 // modify it under the terms of the GNU General Public License as
00013 // published by the Free Software Foundation; either version 2 of
00014 // the License, or (at your option) any later version.
00015 //
00016 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
00017 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00018 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
00019 // GNU General Public License for more details.
00020 //
00021 // You should have received a copy of the GNU Lesser General Public
00022 // License and a copy of the GNU General Public License along with
00023 // Eigen. If not, see <http://www.gnu.org/licenses/>.
00024 
00025 #ifndef EIGEN_UMFPACKSUPPORT_H
00026 #define EIGEN_UMFPACKSUPPORT_H
00027 
00028 /* TODO extract L, extract U, compute det, etc... */
00029 
00030 // generic double/complex<double> wrapper functions:
00031 
00032 inline void umfpack_free_numeric(void **Numeric, double)
00033 { umfpack_di_free_numeric(Numeric); }
00034 
00035 inline void umfpack_free_numeric(void **Numeric, std::complex<double>)
00036 { umfpack_zi_free_numeric(Numeric); }
00037 
00038 inline void umfpack_free_symbolic(void **Symbolic, double)
00039 { umfpack_di_free_symbolic(Symbolic); }
00040 
00041 inline void umfpack_free_symbolic(void **Symbolic, std::complex<double>)
00042 { umfpack_zi_free_symbolic(Symbolic); }
00043 
00044 inline int umfpack_symbolic(int n_row,int n_col,
00045                             const int Ap[], const int Ai[], const double Ax[], void **Symbolic,
00046                             const double Control [UMFPACK_CONTROL], double Info [UMFPACK_INFO])
00047 {
00048   return umfpack_di_symbolic(n_row,n_col,Ap,Ai,Ax,Symbolic,Control,Info);
00049 }
00050 
00051 inline int umfpack_symbolic(int n_row,int n_col,
00052                             const int Ap[], const int Ai[], const std::complex<double> Ax[], void **Symbolic,
00053                             const double Control [UMFPACK_CONTROL], double Info [UMFPACK_INFO])
00054 {
00055   return umfpack_zi_symbolic(n_row,n_col,Ap,Ai,&internal::real_ref(Ax[0]),0,Symbolic,Control,Info);
00056 }
00057 
00058 inline int umfpack_numeric( const int Ap[], const int Ai[], const double Ax[],
00059                             void *Symbolic, void **Numeric,
00060                             const double Control[UMFPACK_CONTROL],double Info [UMFPACK_INFO])
00061 {
00062   return umfpack_di_numeric(Ap,Ai,Ax,Symbolic,Numeric,Control,Info);
00063 }
00064 
00065 inline int umfpack_numeric( const int Ap[], const int Ai[], const std::complex<double> Ax[],
00066                             void *Symbolic, void **Numeric,
00067                             const double Control[UMFPACK_CONTROL],double Info [UMFPACK_INFO])
00068 {
00069   return umfpack_zi_numeric(Ap,Ai,&internal::real_ref(Ax[0]),0,Symbolic,Numeric,Control,Info);
00070 }
00071 
00072 inline int umfpack_solve( int sys, const int Ap[], const int Ai[], const double Ax[],
00073                           double X[], const double B[], void *Numeric,
00074                           const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
00075 {
00076   return umfpack_di_solve(sys,Ap,Ai,Ax,X,B,Numeric,Control,Info);
00077 }
00078 
00079 inline int umfpack_solve( int sys, const int Ap[], const int Ai[], const std::complex<double> Ax[],
00080                           std::complex<double> X[], const std::complex<double> B[], void *Numeric,
00081                           const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
00082 {
00083   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);
00084 }
00085 
00086 inline int umfpack_get_lunz(int *lnz, int *unz, int *n_row, int *n_col, int *nz_udiag, void *Numeric, double)
00087 {
00088   return umfpack_di_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
00089 }
00090 
00091 inline int umfpack_get_lunz(int *lnz, int *unz, int *n_row, int *n_col, int *nz_udiag, void *Numeric, std::complex<double>)
00092 {
00093   return umfpack_zi_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
00094 }
00095 
00096 inline int umfpack_get_numeric(int Lp[], int Lj[], double Lx[], int Up[], int Ui[], double Ux[],
00097                                int P[], int Q[], double Dx[], int *do_recip, double Rs[], void *Numeric)
00098 {
00099   return umfpack_di_get_numeric(Lp,Lj,Lx,Up,Ui,Ux,P,Q,Dx,do_recip,Rs,Numeric);
00100 }
00101 
00102 inline int umfpack_get_numeric(int Lp[], int Lj[], std::complex<double> Lx[], int Up[], int Ui[], std::complex<double> Ux[],
00103                                int P[], int Q[], std::complex<double> Dx[], int *do_recip, double Rs[], void *Numeric)
00104 {
00105   double& lx0_real = internal::real_ref(Lx[0]);
00106   double& ux0_real = internal::real_ref(Ux[0]);
00107   double& dx0_real = internal::real_ref(Dx[0]);
00108   return umfpack_zi_get_numeric(Lp,Lj,Lx?&lx0_real:0,0,Up,Ui,Ux?&ux0_real:0,0,P,Q,
00109                                 Dx?&dx0_real:0,0,do_recip,Rs,Numeric);
00110 }
00111 
00112 inline int umfpack_get_determinant(double *Mx, double *Ex, void *NumericHandle, double User_Info [UMFPACK_INFO])
00113 {
00114   return umfpack_di_get_determinant(Mx,Ex,NumericHandle,User_Info);
00115 }
00116 
00117 inline int umfpack_get_determinant(std::complex<double> *Mx, double *Ex, void *NumericHandle, double User_Info [UMFPACK_INFO])
00118 {
00119   double& mx_real = internal::real_ref(*Mx);
00120   return umfpack_zi_get_determinant(&mx_real,0,Ex,NumericHandle,User_Info);
00121 }
00122 
00123 
00124 template<typename _MatrixType>
00125 class SparseLU<_MatrixType,UmfPack> : public SparseLU<_MatrixType>
00126 {
00127   protected:
00128     typedef SparseLU<_MatrixType> Base;
00129     typedef typename Base::Scalar Scalar;
00130     typedef typename Base::RealScalar RealScalar;
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,Lower|UnitDiag> LMatrixType;
00135     typedef SparseMatrix<Scalar,Upper> UMatrixType;
00136     using Base::m_flags;
00137     using Base::m_status;
00138 
00139   public:
00140     typedef _MatrixType MatrixType;
00141     typedef typename MatrixType::Index Index;
00142 
00143     SparseLU(int flags = NaturalOrdering)
00144       : Base(flags), m_numeric(0)
00145     {
00146     }
00147 
00148     SparseLU(const MatrixType& matrix, int flags = NaturalOrdering)
00149       : Base(flags), m_numeric(0)
00150     {
00151       compute(matrix);
00152     }
00153 
00154     ~SparseLU()
00155     {
00156       if (m_numeric)
00157         umfpack_free_numeric(&m_numeric,Scalar());
00158     }
00159 
00160     inline const LMatrixType& matrixL() const
00161     {
00162       if (m_extractedDataAreDirty) extractData();
00163       return m_l;
00164     }
00165 
00166     inline const UMatrixType& matrixU() const
00167     {
00168       if (m_extractedDataAreDirty) extractData();
00169       return m_u;
00170     }
00171 
00172     inline const IntColVectorType& permutationP() const
00173     {
00174       if (m_extractedDataAreDirty) extractData();
00175       return m_p;
00176     }
00177 
00178     inline const IntRowVectorType& permutationQ() const
00179     {
00180       if (m_extractedDataAreDirty) extractData();
00181       return m_q;
00182     }
00183 
00184     Scalar determinant() const;
00185 
00186     template<typename BDerived, typename XDerived>
00187     bool solve(const MatrixBase<BDerived> &b, MatrixBase<XDerived>* x) const;
00188 
00189     template<typename Rhs>
00190       inline const internal::solve_retval<SparseLU<MatrixType, UmfPack>, Rhs>
00191     solve(const MatrixBase<Rhs>& b) const
00192     {
00193       eigen_assert(true && "SparseLU is not initialized.");
00194       return internal::solve_retval<SparseLU<MatrixType, UmfPack>, Rhs>(*this, b.derived());
00195     }
00196 
00197     void compute(const MatrixType& matrix);
00198 
00199     inline Index cols() const { return m_matrixRef->cols(); }
00200     inline Index rows() const { return m_matrixRef->rows(); }
00201 
00202     inline const MatrixType& matrixLU() const
00203     {
00204       //eigen_assert(m_isInitialized && "LU is not initialized.");
00205       return *m_matrixRef;
00206     }
00207 
00208     const void* numeric() const
00209     {
00210       return m_numeric;
00211     }
00212 
00213   protected:
00214 
00215     void extractData() const;
00216   
00217   protected:
00218     // cached data:
00219     void* m_numeric;
00220     const MatrixType* m_matrixRef;
00221     mutable LMatrixType m_l;
00222     mutable UMatrixType m_u;
00223     mutable IntColVectorType m_p;
00224     mutable IntRowVectorType m_q;
00225     mutable bool m_extractedDataAreDirty;
00226 };
00227 
00228 namespace internal {
00229 
00230 template<typename _MatrixType, typename Rhs>
00231   struct solve_retval<SparseLU<_MatrixType, UmfPack>, Rhs>
00232   : solve_retval_base<SparseLU<_MatrixType, UmfPack>, Rhs>
00233 {
00234   typedef SparseLU<_MatrixType, UmfPack> SpLUDecType;
00235   EIGEN_MAKE_SOLVE_HELPERS(SpLUDecType,Rhs)
00236 
00237   template<typename Dest> void evalTo(Dest& dst) const
00238   {
00239     const int rhsCols = rhs().cols();
00240 
00241     eigen_assert((Rhs::Flags&RowMajorBit)==0 && "UmfPack backend does not support non col-major rhs yet");
00242     eigen_assert((Dest::Flags&RowMajorBit)==0 && "UmfPack backend does not support non col-major result yet");
00243 
00244     void* numeric = const_cast<void*>(dec().numeric());
00245 
00246     EIGEN_UNUSED int errorCode = 0;
00247     for (int j=0; j<rhsCols; ++j)
00248     {
00249       errorCode = umfpack_solve(UMFPACK_A,
00250                                 dec().matrixLU()._outerIndexPtr(), dec().matrixLU()._innerIndexPtr(), dec().matrixLU()._valuePtr(),
00251                                 &dst.col(j).coeffRef(0), &rhs().const_cast_derived().col(j).coeffRef(0), numeric, 0, 0);
00252       eigen_assert(!errorCode && "UmfPack could not solve the system.");
00253     }
00254   }
00255     
00256 };
00257 
00258 } // end namespace internal
00259 
00260 template<typename MatrixType>
00261 void SparseLU<MatrixType,UmfPack>::compute(const MatrixType& a)
00262 {
00263   typedef typename MatrixType::Index Index;
00264   const Index rows = a.rows();
00265   const Index cols = a.cols();
00266   eigen_assert((MatrixType::Flags&RowMajorBit)==0 && "Row major matrices are not supported yet");
00267 
00268   m_matrixRef = &a;
00269 
00270   if (m_numeric)
00271     umfpack_free_numeric(&m_numeric,Scalar());
00272 
00273   void* symbolic;
00274   int errorCode = 0;
00275   errorCode = umfpack_symbolic(rows, cols, a._outerIndexPtr(), a._innerIndexPtr(), a._valuePtr(),
00276                                   &symbolic, 0, 0);
00277   if (errorCode==0)
00278     errorCode = umfpack_numeric(a._outerIndexPtr(), a._innerIndexPtr(), a._valuePtr(),
00279                                    symbolic, &m_numeric, 0, 0);
00280 
00281   umfpack_free_symbolic(&symbolic,Scalar());
00282 
00283   m_extractedDataAreDirty = true;
00284 
00285   Base::m_succeeded = (errorCode==0);
00286 }
00287 
00288 template<typename MatrixType>
00289 void SparseLU<MatrixType,UmfPack>::extractData() const
00290 {
00291   if (m_extractedDataAreDirty)
00292   {
00293     // get size of the data
00294     int lnz, unz, rows, cols, nz_udiag;
00295     umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());
00296 
00297     // allocate data
00298     m_l.resize(rows,std::min(rows,cols));
00299     m_l.resizeNonZeros(lnz);
00300     
00301     m_u.resize(std::min(rows,cols),cols);
00302     m_u.resizeNonZeros(unz);
00303 
00304     m_p.resize(rows);
00305     m_q.resize(cols);
00306 
00307     // extract
00308     umfpack_get_numeric(m_l._outerIndexPtr(), m_l._innerIndexPtr(), m_l._valuePtr(),
00309                         m_u._outerIndexPtr(), m_u._innerIndexPtr(), m_u._valuePtr(),
00310                         m_p.data(), m_q.data(), 0, 0, 0, m_numeric);
00311     
00312     m_extractedDataAreDirty = false;
00313   }
00314 }
00315 
00316 template<typename MatrixType>
00317 typename SparseLU<MatrixType,UmfPack>::Scalar SparseLU<MatrixType,UmfPack>::determinant() const
00318 {
00319   Scalar det;
00320   umfpack_get_determinant(&det, 0, m_numeric, 0);
00321   return det;
00322 }
00323 
00324 template<typename MatrixType>
00325 template<typename BDerived,typename XDerived>
00326 bool SparseLU<MatrixType,UmfPack>::solve(const MatrixBase<BDerived> &b, MatrixBase<XDerived> *x) const
00327 {
00328   //const int size = m_matrix.rows();
00329   const int rhsCols = b.cols();
00330 //   eigen_assert(size==b.rows());
00331   eigen_assert((BDerived::Flags&RowMajorBit)==0 && "UmfPack backend does not support non col-major rhs yet");
00332   eigen_assert((XDerived::Flags&RowMajorBit)==0 && "UmfPack backend does not support non col-major result yet");
00333 
00334   int errorCode;
00335   for (int j=0; j<rhsCols; ++j)
00336   {
00337     errorCode = umfpack_solve(UMFPACK_A,
00338         m_matrixRef->_outerIndexPtr(), m_matrixRef->_innerIndexPtr(), m_matrixRef->_valuePtr(),
00339         &x->col(j).coeffRef(0), &b.const_cast_derived().col(j).coeffRef(0), m_numeric, 0, 0);
00340     if (errorCode!=0)
00341       return false;
00342   }
00343 //   errorCode = umfpack_di_solve(UMFPACK_A,
00344 //       m_matrixRef._outerIndexPtr(), m_matrixRef._innerIndexPtr(), m_matrixRef._valuePtr(),
00345 //       x->derived().data(), b.derived().data(), m_numeric, 0, 0);
00346 
00347   return true;
00348 }
00349 
00350 #endif // EIGEN_UMFPACKSUPPORT_H


re_vision
Author(s): Dorian Galvez-Lopez
autogenerated on Sun Jan 5 2014 11:33:35