Eigen  3.2.91
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 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 
124 template<typename _MatrixType>
125 class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
126 {
127  protected:
129  using Base::m_isInitialized;
130  public:
131  using Base::_solve_impl;
132  typedef _MatrixType MatrixType;
133  typedef typename MatrixType::Scalar Scalar;
134  typedef typename MatrixType::RealScalar RealScalar;
135  typedef typename MatrixType::StorageIndex StorageIndex;
142 
143  public:
144 
145  UmfPackLU()
146  : m_dummy(0,0), mp_matrix(m_dummy)
147  {
148  init();
149  }
150 
151  explicit UmfPackLU(const MatrixType& matrix)
152  : mp_matrix(matrix)
153  {
154  init();
155  compute(matrix);
156  }
157 
158  ~UmfPackLU()
159  {
160  if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
161  if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
162  }
163 
164  inline Index rows() const { return mp_matrix.rows(); }
165  inline Index cols() const { return mp_matrix.cols(); }
166 
173  {
174  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
175  return m_info;
176  }
177 
178  inline const LUMatrixType& matrixL() const
179  {
180  if (m_extractedDataAreDirty) extractData();
181  return m_l;
182  }
183 
184  inline const LUMatrixType& matrixU() const
185  {
186  if (m_extractedDataAreDirty) extractData();
187  return m_u;
188  }
189 
190  inline const IntColVectorType& permutationP() const
191  {
192  if (m_extractedDataAreDirty) extractData();
193  return m_p;
194  }
195 
196  inline const IntRowVectorType& permutationQ() const
197  {
198  if (m_extractedDataAreDirty) extractData();
199  return m_q;
200  }
201 
206  template<typename InputMatrixType>
207  void compute(const InputMatrixType& matrix)
208  {
209  if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
210  if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
211  grab(matrix.derived());
212  analyzePattern_impl();
213  factorize_impl();
214  }
215 
222  template<typename InputMatrixType>
223  void analyzePattern(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 
228  grab(matrix.derived());
229 
230  analyzePattern_impl();
231  }
232 
239  template<typename InputMatrixType>
240  void factorize(const InputMatrixType& matrix)
241  {
242  eigen_assert(m_analysisIsOk && "UmfPackLU: you must first call analyzePattern()");
243  if(m_numeric)
244  umfpack_free_numeric(&m_numeric,Scalar());
245 
246  grab(matrix.derived());
247 
248  factorize_impl();
249  }
250 
252  template<typename BDerived,typename XDerived>
253  bool _solve_impl(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const;
254 
255  Scalar determinant() const;
256 
257  void extractData() const;
258 
259  protected:
260 
261  void init()
262  {
263  m_info = InvalidInput;
264  m_isInitialized = false;
265  m_numeric = 0;
266  m_symbolic = 0;
267  m_extractedDataAreDirty = true;
268  }
269 
270  void analyzePattern_impl()
271  {
272  int errorCode = 0;
273  errorCode = umfpack_symbolic(internal::convert_index<int>(mp_matrix.rows()),
274  internal::convert_index<int>(mp_matrix.cols()),
275  mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
276  &m_symbolic, 0, 0);
277 
278  m_isInitialized = true;
279  m_info = errorCode ? InvalidInput : Success;
280  m_analysisIsOk = true;
281  m_factorizationIsOk = false;
282  m_extractedDataAreDirty = true;
283  }
284 
285  void factorize_impl()
286  {
287  int errorCode;
288  errorCode = umfpack_numeric(mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
289  m_symbolic, &m_numeric, 0, 0);
290 
291  m_info = errorCode ? NumericalIssue : Success;
292  m_factorizationIsOk = true;
293  m_extractedDataAreDirty = true;
294  }
295 
296  template<typename MatrixDerived>
297  void grab(const EigenBase<MatrixDerived> &A)
298  {
299  mp_matrix.~UmfpackMatrixRef();
300  ::new (&mp_matrix) UmfpackMatrixRef(A.derived());
301  }
302 
303  void grab(const UmfpackMatrixRef &A)
304  {
305  if(&(A.derived()) != &mp_matrix)
306  {
307  mp_matrix.~UmfpackMatrixRef();
308  ::new (&mp_matrix) UmfpackMatrixRef(A);
309  }
310  }
311 
312  // cached data to reduce reallocation, etc.
313  mutable LUMatrixType m_l;
314  mutable LUMatrixType m_u;
315  mutable IntColVectorType m_p;
316  mutable IntRowVectorType m_q;
317 
318  UmfpackMatrixType m_dummy;
319  UmfpackMatrixRef mp_matrix;
320 
321  void* m_numeric;
322  void* m_symbolic;
323 
324  mutable ComputationInfo m_info;
325  int m_factorizationIsOk;
326  int m_analysisIsOk;
327  mutable bool m_extractedDataAreDirty;
328 
329  private:
330  UmfPackLU(UmfPackLU& ) { }
331 };
332 
333 
334 template<typename MatrixType>
335 void UmfPackLU<MatrixType>::extractData() const
336 {
337  if (m_extractedDataAreDirty)
338  {
339  // get size of the data
340  int lnz, unz, rows, cols, nz_udiag;
341  umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());
342 
343  // allocate data
344  m_l.resize(rows,(std::min)(rows,cols));
345  m_l.resizeNonZeros(lnz);
346 
347  m_u.resize((std::min)(rows,cols),cols);
348  m_u.resizeNonZeros(unz);
349 
350  m_p.resize(rows);
351  m_q.resize(cols);
352 
353  // extract
354  umfpack_get_numeric(m_l.outerIndexPtr(), m_l.innerIndexPtr(), m_l.valuePtr(),
355  m_u.outerIndexPtr(), m_u.innerIndexPtr(), m_u.valuePtr(),
356  m_p.data(), m_q.data(), 0, 0, 0, m_numeric);
357 
358  m_extractedDataAreDirty = false;
359  }
360 }
361 
362 template<typename MatrixType>
363 typename UmfPackLU<MatrixType>::Scalar UmfPackLU<MatrixType>::determinant() const
364 {
365  Scalar det;
366  umfpack_get_determinant(&det, 0, m_numeric, 0);
367  return det;
368 }
369 
370 template<typename MatrixType>
371 template<typename BDerived,typename XDerived>
372 bool UmfPackLU<MatrixType>::_solve_impl(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const
373 {
374  Index rhsCols = b.cols();
375  eigen_assert((BDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major rhs yet");
376  eigen_assert((XDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major result yet");
377  eigen_assert(b.derived().data() != x.derived().data() && " Umfpack does not support inplace solve");
378 
379  int errorCode;
380  Scalar* x_ptr = 0;
381  Matrix<Scalar,Dynamic,1> x_tmp;
382  if(x.innerStride()!=1)
383  {
384  x_tmp.resize(x.rows());
385  x_ptr = x_tmp.data();
386  }
387  for (int j=0; j<rhsCols; ++j)
388  {
389  if(x.innerStride()==1)
390  x_ptr = &x.col(j).coeffRef(0);
391  errorCode = umfpack_solve(UMFPACK_A,
392  mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
393  x_ptr, &b.const_cast_derived().col(j).coeffRef(0), m_numeric, 0, 0);
394  if(x.innerStride()!=1)
395  x.col(j) = x_tmp;
396  if (errorCode!=0)
397  return false;
398  }
399 
400  return true;
401 }
402 
403 } // end namespace Eigen
404 
405 #endif // EIGEN_UMFPACKSUPPORT_H
A sparse LU factorization and solver based on UmfPack.
Definition: UmfPackSupport.h:125
void factorize(const InputMatrixType &matrix)
Definition: UmfPackSupport.h:240
void compute(const InputMatrixType &matrix)
Definition: UmfPackSupport.h:207
void analyzePattern(const InputMatrixType &matrix)
Definition: UmfPackSupport.h:223
A base class for sparse solvers.
Definition: SparseSolverBase.h:53
Definition: LDLT.h:16
const StorageIndex * outerIndexPtr() const
Definition: SparseMatrix.h:160
const unsigned int RowMajorBit
Definition: Constants.h:53
void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:252
Definition: Constants.h:426
void resize(Index rows, Index cols)
Definition: SparseMatrix.h:605
Definition: Constants.h:431
Definition: Constants.h:424
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: UmfPackSupport.h:172
const StorageIndex * innerIndexPtr() const
Definition: SparseMatrix.h:151
const Scalar * data() const
Definition: PlainObjectBase.h:228
const Scalar * valuePtr() const
Definition: SparseMatrix.h:142
ComputationInfo
Definition: Constants.h:422
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48