compbio
UmfPackSupport.h
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@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 
20 inline void umfpack_defaults(double control[UMFPACK_CONTROL], double)
21 { umfpack_di_defaults(control); }
22 
23 inline void umfpack_defaults(double control[UMFPACK_CONTROL], std::complex<double>)
24 { umfpack_zi_defaults(control); }
25 
26 inline void umfpack_free_numeric(void **Numeric, double)
27 { umfpack_di_free_numeric(Numeric); *Numeric = 0; }
28 
29 inline void umfpack_free_numeric(void **Numeric, std::complex<double>)
30 { umfpack_zi_free_numeric(Numeric); *Numeric = 0; }
31 
32 inline void umfpack_free_symbolic(void **Symbolic, double)
33 { umfpack_di_free_symbolic(Symbolic); *Symbolic = 0; }
34 
35 inline void umfpack_free_symbolic(void **Symbolic, std::complex<double>)
36 { umfpack_zi_free_symbolic(Symbolic); *Symbolic = 0; }
37 
38 inline int umfpack_symbolic(int n_row,int n_col,
39  const int Ap[], const int Ai[], const double Ax[], void **Symbolic,
40  const double Control [UMFPACK_CONTROL], double Info [UMFPACK_INFO])
41 {
42  return umfpack_di_symbolic(n_row,n_col,Ap,Ai,Ax,Symbolic,Control,Info);
43 }
44 
45 inline int umfpack_symbolic(int n_row,int n_col,
46  const int Ap[], const int Ai[], const std::complex<double> Ax[], void **Symbolic,
47  const double Control [UMFPACK_CONTROL], double Info [UMFPACK_INFO])
48 {
49  return umfpack_zi_symbolic(n_row,n_col,Ap,Ai,&numext::real_ref(Ax[0]),0,Symbolic,Control,Info);
50 }
51 
52 inline int umfpack_numeric( const int Ap[], const int Ai[], const double Ax[],
53  void *Symbolic, void **Numeric,
54  const double Control[UMFPACK_CONTROL],double Info [UMFPACK_INFO])
55 {
56  return umfpack_di_numeric(Ap,Ai,Ax,Symbolic,Numeric,Control,Info);
57 }
58 
59 inline int umfpack_numeric( const int Ap[], const int Ai[], const std::complex<double> Ax[],
60  void *Symbolic, void **Numeric,
61  const double Control[UMFPACK_CONTROL],double Info [UMFPACK_INFO])
62 {
63  return umfpack_zi_numeric(Ap,Ai,&numext::real_ref(Ax[0]),0,Symbolic,Numeric,Control,Info);
64 }
65 
66 inline int umfpack_solve( int sys, const int Ap[], const int Ai[], const double Ax[],
67  double X[], const double B[], void *Numeric,
68  const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
69 {
70  return umfpack_di_solve(sys,Ap,Ai,Ax,X,B,Numeric,Control,Info);
71 }
72 
73 inline int umfpack_solve( int sys, const int Ap[], const int Ai[], const std::complex<double> Ax[],
74  std::complex<double> X[], const std::complex<double> B[], void *Numeric,
75  const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
76 {
77  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);
78 }
79 
80 inline int umfpack_get_lunz(int *lnz, int *unz, int *n_row, int *n_col, int *nz_udiag, void *Numeric, double)
81 {
82  return umfpack_di_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
83 }
84 
85 inline int umfpack_get_lunz(int *lnz, int *unz, int *n_row, int *n_col, int *nz_udiag, void *Numeric, std::complex<double>)
86 {
87  return umfpack_zi_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
88 }
89 
90 inline int umfpack_get_numeric(int Lp[], int Lj[], double Lx[], int Up[], int Ui[], double Ux[],
91  int P[], int Q[], double Dx[], int *do_recip, double Rs[], void *Numeric)
92 {
93  return umfpack_di_get_numeric(Lp,Lj,Lx,Up,Ui,Ux,P,Q,Dx,do_recip,Rs,Numeric);
94 }
95 
96 inline int umfpack_get_numeric(int Lp[], int Lj[], std::complex<double> Lx[], int Up[], int Ui[], std::complex<double> Ux[],
97  int P[], int Q[], std::complex<double> Dx[], int *do_recip, double Rs[], void *Numeric)
98 {
99  double& lx0_real = numext::real_ref(Lx[0]);
100  double& ux0_real = numext::real_ref(Ux[0]);
101  double& dx0_real = numext::real_ref(Dx[0]);
102  return umfpack_zi_get_numeric(Lp,Lj,Lx?&lx0_real:0,0,Up,Ui,Ux?&ux0_real:0,0,P,Q,
103  Dx?&dx0_real:0,0,do_recip,Rs,Numeric);
104 }
105 
106 inline int umfpack_get_determinant(double *Mx, double *Ex, void *NumericHandle, double User_Info [UMFPACK_INFO])
107 {
108  return umfpack_di_get_determinant(Mx,Ex,NumericHandle,User_Info);
109 }
110 
111 inline int umfpack_get_determinant(std::complex<double> *Mx, double *Ex, void *NumericHandle, double User_Info [UMFPACK_INFO])
112 {
113  double& mx_real = numext::real_ref(*Mx);
114  return umfpack_zi_get_determinant(&mx_real,0,Ex,NumericHandle,User_Info);
115 }
116 
117 
133 template<typename _MatrixType>
134 class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
135 {
136  protected:
138  using Base::m_isInitialized;
139  public:
140  using Base::_solve_impl;
141  typedef _MatrixType MatrixType;
142  typedef typename MatrixType::Scalar Scalar;
143  typedef typename MatrixType::RealScalar RealScalar;
144  typedef typename MatrixType::StorageIndex StorageIndex;
151  enum {
152  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
153  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
154  };
155 
156  public:
157 
159 
160  UmfPackLU()
161  : m_dummy(0,0), mp_matrix(m_dummy)
162  {
163  init();
164  }
165 
166  template<typename InputMatrixType>
167  explicit UmfPackLU(const InputMatrixType& matrix)
168  : mp_matrix(matrix)
169  {
170  init();
171  compute(matrix);
172  }
173 
174  ~UmfPackLU()
175  {
176  if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
177  if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
178  }
179 
180  inline Index rows() const { return mp_matrix.rows(); }
181  inline Index cols() const { return mp_matrix.cols(); }
182 
189  {
190  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
191  return m_info;
192  }
193 
194  inline const LUMatrixType& matrixL() const
195  {
196  if (m_extractedDataAreDirty) extractData();
197  return m_l;
198  }
199 
200  inline const LUMatrixType& matrixU() const
201  {
202  if (m_extractedDataAreDirty) extractData();
203  return m_u;
204  }
205 
206  inline const IntColVectorType& permutationP() const
207  {
208  if (m_extractedDataAreDirty) extractData();
209  return m_p;
210  }
211 
212  inline const IntRowVectorType& permutationQ() const
213  {
214  if (m_extractedDataAreDirty) extractData();
215  return m_q;
216  }
217 
222  template<typename InputMatrixType>
223  void compute(const InputMatrixType& matrix)
224  {
225  if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
226  if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
227  grab(matrix.derived());
228  analyzePattern_impl();
229  factorize_impl();
230  }
231 
238  template<typename InputMatrixType>
239  void analyzePattern(const InputMatrixType& matrix)
240  {
241  if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
242  if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
243 
244  grab(matrix.derived());
245 
246  analyzePattern_impl();
247  }
248 
254  inline int umfpackFactorizeReturncode() const
255  {
256  eigen_assert(m_numeric && "UmfPackLU: you must first call factorize()");
257  return m_fact_errorCode;
258  }
259 
266  inline const UmfpackControl& umfpackControl() const
267  {
268  return m_control;
269  }
270 
277  inline UmfpackControl& umfpackControl()
278  {
279  return m_control;
280  }
281 
288  template<typename InputMatrixType>
289  void factorize(const InputMatrixType& matrix)
290  {
291  eigen_assert(m_analysisIsOk && "UmfPackLU: you must first call analyzePattern()");
292  if(m_numeric)
293  umfpack_free_numeric(&m_numeric,Scalar());
294 
295  grab(matrix.derived());
296 
297  factorize_impl();
298  }
299 
301  template<typename BDerived,typename XDerived>
302  bool _solve_impl(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const;
303 
304  Scalar determinant() const;
305 
306  void extractData() const;
307 
308  protected:
309 
310  void init()
311  {
312  m_info = InvalidInput;
313  m_isInitialized = false;
314  m_numeric = 0;
315  m_symbolic = 0;
316  m_extractedDataAreDirty = true;
317  }
318 
319  void analyzePattern_impl()
320  {
321  umfpack_defaults(m_control.data(), Scalar());
322  int errorCode = 0;
323  errorCode = umfpack_symbolic(internal::convert_index<int>(mp_matrix.rows()),
324  internal::convert_index<int>(mp_matrix.cols()),
325  mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
326  &m_symbolic, m_control.data(), 0);
327 
328  m_isInitialized = true;
329  m_info = errorCode ? InvalidInput : Success;
330  m_analysisIsOk = true;
331  m_factorizationIsOk = false;
332  m_extractedDataAreDirty = true;
333  }
334 
335  void factorize_impl()
336  {
337  m_fact_errorCode = umfpack_numeric(mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
338  m_symbolic, &m_numeric, m_control.data(), 0);
339 
340  m_info = m_fact_errorCode == UMFPACK_OK ? Success : NumericalIssue;
341  m_factorizationIsOk = true;
342  m_extractedDataAreDirty = true;
343  }
344 
345  template<typename MatrixDerived>
346  void grab(const EigenBase<MatrixDerived> &A)
347  {
348  mp_matrix.~UmfpackMatrixRef();
349  ::new (&mp_matrix) UmfpackMatrixRef(A.derived());
350  }
351 
352  void grab(const UmfpackMatrixRef &A)
353  {
354  if(&(A.derived()) != &mp_matrix)
355  {
356  mp_matrix.~UmfpackMatrixRef();
357  ::new (&mp_matrix) UmfpackMatrixRef(A);
358  }
359  }
360 
361  // cached data to reduce reallocation, etc.
362  mutable LUMatrixType m_l;
363  int m_fact_errorCode;
364  UmfpackControl m_control;
365 
366  mutable LUMatrixType m_u;
367  mutable IntColVectorType m_p;
368  mutable IntRowVectorType m_q;
369 
370  UmfpackMatrixType m_dummy;
371  UmfpackMatrixRef mp_matrix;
372 
373  void* m_numeric;
374  void* m_symbolic;
375 
376  mutable ComputationInfo m_info;
377  int m_factorizationIsOk;
378  int m_analysisIsOk;
379  mutable bool m_extractedDataAreDirty;
380 
381  private:
382  UmfPackLU(const UmfPackLU& ) { }
383 };
384 
385 
386 template<typename MatrixType>
388 {
389  if (m_extractedDataAreDirty)
390  {
391  // get size of the data
392  int lnz, unz, rows, cols, nz_udiag;
393  umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());
394 
395  // allocate data
396  m_l.resize(rows,(std::min)(rows,cols));
397  m_l.resizeNonZeros(lnz);
398 
399  m_u.resize((std::min)(rows,cols),cols);
400  m_u.resizeNonZeros(unz);
401 
402  m_p.resize(rows);
403  m_q.resize(cols);
404 
405  // extract
406  umfpack_get_numeric(m_l.outerIndexPtr(), m_l.innerIndexPtr(), m_l.valuePtr(),
407  m_u.outerIndexPtr(), m_u.innerIndexPtr(), m_u.valuePtr(),
408  m_p.data(), m_q.data(), 0, 0, 0, m_numeric);
409 
410  m_extractedDataAreDirty = false;
411  }
412 }
413 
414 template<typename MatrixType>
415 typename UmfPackLU<MatrixType>::Scalar UmfPackLU<MatrixType>::determinant() const
416 {
417  Scalar det;
418  umfpack_get_determinant(&det, 0, m_numeric, 0);
419  return det;
420 }
421 
422 template<typename MatrixType>
423 template<typename BDerived,typename XDerived>
425 {
426  Index rhsCols = b.cols();
427  eigen_assert((BDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major rhs yet");
428  eigen_assert((XDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major result yet");
429  eigen_assert(b.derived().data() != x.derived().data() && " Umfpack does not support inplace solve");
430 
431  int errorCode;
432  Scalar* x_ptr = 0;
434  if(x.innerStride()!=1)
435  {
436  x_tmp.resize(x.rows());
437  x_ptr = x_tmp.data();
438  }
439  for (int j=0; j<rhsCols; ++j)
440  {
441  if(x.innerStride()==1)
442  x_ptr = &x.col(j).coeffRef(0);
443  errorCode = umfpack_solve(UMFPACK_A,
444  mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
445  x_ptr, &b.const_cast_derived().col(j).coeffRef(0), m_numeric, m_control.data(), 0);
446  if(x.innerStride()!=1)
447  x.col(j) = x_tmp;
448  if (errorCode!=0)
449  return false;
450  }
451 
452  return true;
453 }
454 
455 } // end namespace Eigen
456 
457 #endif // EIGEN_UMFPACKSUPPORT_H
EIGEN_DEVICE_FUNC ColXpr col(Index i)
Definition: DenseBase.h:839
A sparse LU factorization and solver based on UmfPack.
Definition: UmfPackSupport.h:134
const Scalar * valuePtr() const
Definition: SparseMatrix.h:144
void factorize(const InputMatrixType &matrix)
Performs a numeric decomposition of matrix.
Definition: UmfPackSupport.h:289
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar * data() const
Definition: PlainObjectBase.h:249
void compute(const InputMatrixType &matrix)
Computes the sparse Cholesky decomposition of matrix Note that the matrix should be column-major...
Definition: UmfPackSupport.h:223
void analyzePattern(const InputMatrixType &matrix)
Performs a symbolic decomposition on the sparcity of matrix.
Definition: UmfPackSupport.h:239
Definition: mainwindow.h:6
A base class for sparse solvers.
Definition: SparseSolverBase.h:67
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:85
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: UmfPackSupport.h:188
const unsigned int RowMajorBit
for a matrix, this means that the storage order is row-major.
Definition: Constants.h:61
Common base class for all classes T such that MatrixBase has an operator=(T) and a constructor Matrix...
Definition: EigenBase.h:28
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Resizes *this to a rows x cols matrix.
Definition: PlainObjectBase.h:273
The provided data did not satisfy the prerequisites.
Definition: Constants.h:434
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:33
void resize(Index rows, Index cols)
Resizes the matrix to a rows x cols matrix and initializes it to zero.
Definition: SparseMatrix.h:617
const UmfpackControl & umfpackControl() const
Provides access to the control settings array used by UmfPack.
Definition: UmfPackSupport.h:266
The inputs are invalid, or the algorithm has been improperly called.
Definition: Constants.h:439
Computation was successful.
Definition: Constants.h:432
const StorageIndex * innerIndexPtr() const
Definition: SparseMatrix.h:153
Definition: TutorialInplaceLU.cpp:2
int umfpackFactorizeReturncode() const
Provides the return status code returned by UmfPack during the numeric factorization.
Definition: UmfPackSupport.h:254
const StorageIndex * outerIndexPtr() const
Definition: SparseMatrix.h:162
UmfpackControl & umfpackControl()
Provides access to the control settings array used by UmfPack.
Definition: UmfPackSupport.h:277
ComputationInfo
Enum for reporting the status of a computation.
Definition: Constants.h:430
EIGEN_DEVICE_FUNC Derived & derived()
Definition: EigenBase.h:44
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48