1 // This file is part of Eigen, a lightweight C++ template library 2 // for linear algebra. 3 // 4 // Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud (at) inria.fr> 5 // 6 // This Source Code Form is subject to the terms of the Mozilla 7 // Public License v. 2.0. If a copy of the MPL was not distributed 8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 9 10 #ifndef EIGEN_UMFPACKSUPPORT_H 11 #define EIGEN_UMFPACKSUPPORT_H 12 13 namespace Eigen { 14 15 /* TODO extract L, extract U, compute det, etc... */ 16 17 // generic double/complex<double> wrapper functions: 18 19 inline void umfpack_free_numeric(void **Numeric, double) 20 { umfpack_di_free_numeric(Numeric); *Numeric = 0; } 21 22 inline void umfpack_free_numeric(void **Numeric, std::complex<double>) 23 { umfpack_zi_free_numeric(Numeric); *Numeric = 0; } 24 25 inline void umfpack_free_symbolic(void **Symbolic, double) 26 { umfpack_di_free_symbolic(Symbolic); *Symbolic = 0; } 27 28 inline void umfpack_free_symbolic(void **Symbolic, std::complex<double>) 29 { umfpack_zi_free_symbolic(Symbolic); *Symbolic = 0; } 30 31 inline int umfpack_symbolic(int n_row,int n_col, 32 const int Ap[], const int Ai[], const double Ax[], void **Symbolic, 33 const double Control [UMFPACK_CONTROL], double Info [UMFPACK_INFO]) 34 { 35 return umfpack_di_symbolic(n_row,n_col,Ap,Ai,Ax,Symbolic,Control,Info); 36 } 37 38 inline int umfpack_symbolic(int n_row,int n_col, 39 const int Ap[], const int Ai[], const std::complex<double> Ax[], void **Symbolic, 40 const double Control [UMFPACK_CONTROL], double Info [UMFPACK_INFO]) 41 { 42 return umfpack_zi_symbolic(n_row,n_col,Ap,Ai,&numext::real_ref(Ax[0]),0,Symbolic,Control,Info); 43 } 44 45 inline int umfpack_numeric( const int Ap[], const int Ai[], const double Ax[], 46 void *Symbolic, void **Numeric, 47 const double Control[UMFPACK_CONTROL],double Info [UMFPACK_INFO]) 48 { 49 return umfpack_di_numeric(Ap,Ai,Ax,Symbolic,Numeric,Control,Info); 50 } 51 52 inline int umfpack_numeric( const int Ap[], const int Ai[], const std::complex<double> Ax[], 53 void *Symbolic, void **Numeric, 54 const double Control[UMFPACK_CONTROL],double Info [UMFPACK_INFO]) 55 { 56 return umfpack_zi_numeric(Ap,Ai,&numext::real_ref(Ax[0]),0,Symbolic,Numeric,Control,Info); 57 } 58 59 inline int umfpack_solve( int sys, const int Ap[], const int Ai[], const double Ax[], 60 double X[], const double B[], void *Numeric, 61 const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO]) 62 { 63 return umfpack_di_solve(sys,Ap,Ai,Ax,X,B,Numeric,Control,Info); 64 } 65 66 inline int umfpack_solve( int sys, const int Ap[], const int Ai[], const std::complex<double> Ax[], 67 std::complex<double> X[], const std::complex<double> B[], void *Numeric, 68 const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO]) 69 { 70 return umfpack_zi_solve(sys,Ap,Ai,&numext::real_ref(Ax[0]),0,&numext::real_ref(X[0]),0,&numext::real_ref(B[0]),0,Numeric,Control,Info); 71 } 72 73 inline int umfpack_get_lunz(int *lnz, int *unz, int *n_row, int *n_col, int *nz_udiag, void *Numeric, double) 74 { 75 return umfpack_di_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric); 76 } 77 78 inline int umfpack_get_lunz(int *lnz, int *unz, int *n_row, int *n_col, int *nz_udiag, void *Numeric, std::complex<double>) 79 { 80 return umfpack_zi_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric); 81 } 82 83 inline int umfpack_get_numeric(int Lp[], int Lj[], double Lx[], int Up[], int Ui[], double Ux[], 84 int P[], int Q[], double Dx[], int *do_recip, double Rs[], void *Numeric) 85 { 86 return umfpack_di_get_numeric(Lp,Lj,Lx,Up,Ui,Ux,P,Q,Dx,do_recip,Rs,Numeric); 87 } 88 89 inline int umfpack_get_numeric(int Lp[], int Lj[], std::complex<double> Lx[], int Up[], int Ui[], std::complex<double> Ux[], 90 int P[], int Q[], std::complex<double> Dx[], int *do_recip, double Rs[], void *Numeric) 91 { 92 double& lx0_real = numext::real_ref(Lx[0]); 93 double& ux0_real = numext::real_ref(Ux[0]); 94 double& dx0_real = numext::real_ref(Dx[0]); 95 return umfpack_zi_get_numeric(Lp,Lj,Lx?&lx0_real:0,0,Up,Ui,Ux?&ux0_real:0,0,P,Q, 96 Dx?&dx0_real:0,0,do_recip,Rs,Numeric); 97 } 98 99 inline int umfpack_get_determinant(double *Mx, double *Ex, void *NumericHandle, double User_Info [UMFPACK_INFO]) 100 { 101 return umfpack_di_get_determinant(Mx,Ex,NumericHandle,User_Info); 102 } 103 104 inline int umfpack_get_determinant(std::complex<double> *Mx, double *Ex, void *NumericHandle, double User_Info [UMFPACK_INFO]) 105 { 106 double& mx_real = numext::real_ref(*Mx); 107 return umfpack_zi_get_determinant(&mx_real,0,Ex,NumericHandle,User_Info); 108 } 109 110 namespace internal { 111 template<typename T> struct umfpack_helper_is_sparse_plain : false_type {}; 112 template<typename Scalar, int Options, typename StorageIndex> 113 struct umfpack_helper_is_sparse_plain<SparseMatrix<Scalar,Options,StorageIndex> > 114 : true_type {}; 115 template<typename Scalar, int Options, typename StorageIndex> 116 struct umfpack_helper_is_sparse_plain<MappedSparseMatrix<Scalar,Options,StorageIndex> > 117 : true_type {}; 118 } 119 120 /** \ingroup UmfPackSupport_Module 121 * \brief A sparse LU factorization and solver based on UmfPack 122 * 123 * This class allows to solve for A.X = B sparse linear problems via a LU factorization 124 * using the UmfPack library. The sparse matrix A must be squared and full rank. 125 * The vectors or matrices X and B can be either dense or sparse. 126 * 127 * \warning The input matrix A should be in a \b compressed and \b column-major form. 128 * Otherwise an expensive copy will be made. You can call the inexpensive makeCompressed() to get a compressed matrix. 129 * \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<> 130 * 131 * \sa \ref TutorialSparseDirectSolvers 132 */ 133 template<typename _MatrixType> 134 class UmfPackLU : internal::noncopyable 135 { 136 public: 137 typedef _MatrixType MatrixType; 138 typedef typename MatrixType::Scalar Scalar; 139 typedef typename MatrixType::RealScalar RealScalar; 140 typedef typename MatrixType::Index Index; 141 typedef Matrix<Scalar,Dynamic,1> Vector; 142 typedef Matrix<int, 1, MatrixType::ColsAtCompileTime> IntRowVectorType; 143 typedef Matrix<int, MatrixType::RowsAtCompileTime, 1> IntColVectorType; 144 typedef SparseMatrix<Scalar> LUMatrixType; 145 typedef SparseMatrix<Scalar,ColMajor,int> UmfpackMatrixType; 146 147 public: 148 149 UmfPackLU() { init(); } 150 151 UmfPackLU(const MatrixType& matrix) 152 { 153 init(); 154 compute(matrix); 155 } 156 157 ~UmfPackLU() 158 { 159 if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar()); 160 if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar()); 161 } 162 163 inline Index rows() const { return m_copyMatrix.rows(); } 164 inline Index cols() const { return m_copyMatrix.cols(); } 165 166 /** \brief Reports whether previous computation was successful. 167 * 168 * \returns \c Success if computation was succesful, 169 * \c NumericalIssue if the matrix.appears to be negative. 170 */ 171 ComputationInfo info() const 172 { 173 eigen_assert(m_isInitialized && "Decomposition is not initialized."); 174 return m_info; 175 } 176 177 inline const LUMatrixType& matrixL() const 178 { 179 if (m_extractedDataAreDirty) extractData(); 180 return m_l; 181 } 182 183 inline const LUMatrixType& matrixU() const 184 { 185 if (m_extractedDataAreDirty) extractData(); 186 return m_u; 187 } 188 189 inline const IntColVectorType& permutationP() const 190 { 191 if (m_extractedDataAreDirty) extractData(); 192 return m_p; 193 } 194 195 inline const IntRowVectorType& permutationQ() const 196 { 197 if (m_extractedDataAreDirty) extractData(); 198 return m_q; 199 } 200 201 /** Computes the sparse Cholesky decomposition of \a matrix 202 * Note that the matrix should be column-major, and in compressed format for best performance. 203 * \sa SparseMatrix::makeCompressed(). 204 */ 205 template<typename InputMatrixType> 206 void compute(const InputMatrixType& matrix) 207 { 208 if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar()); 209 if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar()); 210 grapInput(matrix.derived()); 211 analyzePattern_impl(); 212 factorize_impl(); 213 } 214 215 /** \returns the solution x of \f$ A x = b \f$ using the current decomposition of A. 216 * 217 * \sa compute() 218 */ 219 template<typename Rhs> 220 inline const internal::solve_retval<UmfPackLU, Rhs> solve(const MatrixBase<Rhs>& b) const 221 { 222 eigen_assert(m_isInitialized && "UmfPackLU is not initialized."); 223 eigen_assert(rows()==b.rows() 224 && "UmfPackLU::solve(): invalid number of rows of the right hand side matrix b"); 225 return internal::solve_retval<UmfPackLU, Rhs>(*this, b.derived()); 226 } 227 228 /** \returns the solution x of \f$ A x = b \f$ using the current decomposition of A. 229 * 230 * \sa compute() 231 */ 232 template<typename Rhs> 233 inline const internal::sparse_solve_retval<UmfPackLU, Rhs> solve(const SparseMatrixBase<Rhs>& b) const 234 { 235 eigen_assert(m_isInitialized && "UmfPackLU is not initialized."); 236 eigen_assert(rows()==b.rows() 237 && "UmfPackLU::solve(): invalid number of rows of the right hand side matrix b"); 238 return internal::sparse_solve_retval<UmfPackLU, Rhs>(*this, b.derived()); 239 } 240 241 /** Performs a symbolic decomposition on the sparcity of \a matrix. 242 * 243 * This function is particularly useful when solving for several problems having the same structure. 244 * 245 * \sa factorize(), compute() 246 */ 247 template<typename InputMatrixType> 248 void analyzePattern(const InputMatrixType& matrix) 249 { 250 if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar()); 251 if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar()); 252 253 grapInput(matrix.derived()); 254 255 analyzePattern_impl(); 256 } 257 258 /** Performs a numeric decomposition of \a matrix 259 * 260 * The given matrix must has the same sparcity than the matrix on which the pattern anylysis has been performed. 261 * 262 * \sa analyzePattern(), compute() 263 */ 264 template<typename InputMatrixType> 265 void factorize(const InputMatrixType& matrix) 266 { 267 eigen_assert(m_analysisIsOk && "UmfPackLU: you must first call analyzePattern()"); 268 if(m_numeric) 269 umfpack_free_numeric(&m_numeric,Scalar()); 270 271 grapInput(matrix.derived()); 272 273 factorize_impl(); 274 } 275 276 #ifndef EIGEN_PARSED_BY_DOXYGEN 277 /** \internal */ 278 template<typename BDerived,typename XDerived> 279 bool _solve(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const; 280 #endif 281 282 Scalar determinant() const; 283 284 void extractData() const; 285 286 protected: 287 288 void init() 289 { 290 m_info = InvalidInput; 291 m_isInitialized = false; 292 m_numeric = 0; 293 m_symbolic = 0; 294 m_outerIndexPtr = 0; 295 m_innerIndexPtr = 0; 296 m_valuePtr = 0; 297 m_extractedDataAreDirty = true; 298 } 299 300 template<typename InputMatrixType> 301 void grapInput_impl(const InputMatrixType& mat, internal::true_type) 302 { 303 m_copyMatrix.resize(mat.rows(), mat.cols()); 304 if( ((MatrixType::Flags&RowMajorBit)==RowMajorBit) || sizeof(typename MatrixType::Index)!=sizeof(int) || !mat.isCompressed() ) 305 { 306 // non supported input -> copy 307 m_copyMatrix = mat; 308 m_outerIndexPtr = m_copyMatrix.outerIndexPtr(); 309 m_innerIndexPtr = m_copyMatrix.innerIndexPtr(); 310 m_valuePtr = m_copyMatrix.valuePtr(); 311 } 312 else 313 { 314 m_outerIndexPtr = mat.outerIndexPtr(); 315 m_innerIndexPtr = mat.innerIndexPtr(); 316 m_valuePtr = mat.valuePtr(); 317 } 318 } 319 320 template<typename InputMatrixType> 321 void grapInput_impl(const InputMatrixType& mat, internal::false_type) 322 { 323 m_copyMatrix = mat; 324 m_outerIndexPtr = m_copyMatrix.outerIndexPtr(); 325 m_innerIndexPtr = m_copyMatrix.innerIndexPtr(); 326 m_valuePtr = m_copyMatrix.valuePtr(); 327 } 328 329 template<typename InputMatrixType> 330 void grapInput(const InputMatrixType& mat) 331 { 332 grapInput_impl(mat, internal::umfpack_helper_is_sparse_plain<InputMatrixType>()); 333 } 334 335 void analyzePattern_impl() 336 { 337 int errorCode = 0; 338 errorCode = umfpack_symbolic(m_copyMatrix.rows(), m_copyMatrix.cols(), m_outerIndexPtr, m_innerIndexPtr, m_valuePtr, 339 &m_symbolic, 0, 0); 340 341 m_isInitialized = true; 342 m_info = errorCode ? InvalidInput : Success; 343 m_analysisIsOk = true; 344 m_factorizationIsOk = false; 345 m_extractedDataAreDirty = true; 346 } 347 348 void factorize_impl() 349 { 350 int errorCode; 351 errorCode = umfpack_numeric(m_outerIndexPtr, m_innerIndexPtr, m_valuePtr, 352 m_symbolic, &m_numeric, 0, 0); 353 354 m_info = errorCode ? NumericalIssue : Success; 355 m_factorizationIsOk = true; 356 m_extractedDataAreDirty = true; 357 } 358 359 // cached data to reduce reallocation, etc. 360 mutable LUMatrixType m_l; 361 mutable LUMatrixType m_u; 362 mutable IntColVectorType m_p; 363 mutable IntRowVectorType m_q; 364 365 UmfpackMatrixType m_copyMatrix; 366 const Scalar* m_valuePtr; 367 const int* m_outerIndexPtr; 368 const int* m_innerIndexPtr; 369 void* m_numeric; 370 void* m_symbolic; 371 372 mutable ComputationInfo m_info; 373 bool m_isInitialized; 374 int m_factorizationIsOk; 375 int m_analysisIsOk; 376 mutable bool m_extractedDataAreDirty; 377 378 private: 379 UmfPackLU(UmfPackLU& ) { } 380 }; 381 382 383 template<typename MatrixType> 384 void UmfPackLU<MatrixType>::extractData() const 385 { 386 if (m_extractedDataAreDirty) 387 { 388 // get size of the data 389 int lnz, unz, rows, cols, nz_udiag; 390 umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar()); 391 392 // allocate data 393 m_l.resize(rows,(std::min)(rows,cols)); 394 m_l.resizeNonZeros(lnz); 395 396 m_u.resize((std::min)(rows,cols),cols); 397 m_u.resizeNonZeros(unz); 398 399 m_p.resize(rows); 400 m_q.resize(cols); 401 402 // extract 403 umfpack_get_numeric(m_l.outerIndexPtr(), m_l.innerIndexPtr(), m_l.valuePtr(), 404 m_u.outerIndexPtr(), m_u.innerIndexPtr(), m_u.valuePtr(), 405 m_p.data(), m_q.data(), 0, 0, 0, m_numeric); 406 407 m_extractedDataAreDirty = false; 408 } 409 } 410 411 template<typename MatrixType> 412 typename UmfPackLU<MatrixType>::Scalar UmfPackLU<MatrixType>::determinant() const 413 { 414 Scalar det; 415 umfpack_get_determinant(&det, 0, m_numeric, 0); 416 return det; 417 } 418 419 template<typename MatrixType> 420 template<typename BDerived,typename XDerived> 421 bool UmfPackLU<MatrixType>::_solve(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const 422 { 423 const int rhsCols = b.cols(); 424 eigen_assert((BDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major rhs yet"); 425 eigen_assert((XDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major result yet"); 426 eigen_assert(b.derived().data() != x.derived().data() && " Umfpack does not support inplace solve"); 427 428 int errorCode; 429 for (int j=0; j<rhsCols; ++j) 430 { 431 errorCode = umfpack_solve(UMFPACK_A, 432 m_outerIndexPtr, m_innerIndexPtr, m_valuePtr, 433 &x.col(j).coeffRef(0), &b.const_cast_derived().col(j).coeffRef(0), m_numeric, 0, 0); 434 if (errorCode!=0) 435 return false; 436 } 437 438 return true; 439 } 440 441 442 namespace internal { 443 444 template<typename _MatrixType, typename Rhs> 445 struct solve_retval<UmfPackLU<_MatrixType>, Rhs> 446 : solve_retval_base<UmfPackLU<_MatrixType>, Rhs> 447 { 448 typedef UmfPackLU<_MatrixType> Dec; 449 EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs) 450 451 template<typename Dest> void evalTo(Dest& dst) const 452 { 453 dec()._solve(rhs(),dst); 454 } 455 }; 456 457 template<typename _MatrixType, typename Rhs> 458 struct sparse_solve_retval<UmfPackLU<_MatrixType>, Rhs> 459 : sparse_solve_retval_base<UmfPackLU<_MatrixType>, Rhs> 460 { 461 typedef UmfPackLU<_MatrixType> Dec; 462 EIGEN_MAKE_SPARSE_SOLVE_HELPERS(Dec,Rhs) 463 464 template<typename Dest> void evalTo(Dest& dst) const 465 { 466 this->defaultEvalTo(dst); 467 } 468 }; 469 470 } // end namespace internal 471 472 } // end namespace Eigen 473 474 #endif // EIGEN_UMFPACKSUPPORT_H 475