UmfPackSupport.h
Go to the documentation of this file.
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 // for compatibility with super old version of umfpack,
14 // not sure this is really needed, but this is harmless.
15 #ifndef SuiteSparse_long
16 #ifdef UF_long
17 #define SuiteSparse_long UF_long
18 #else
19 #error neither SuiteSparse_long nor UF_long are defined
20 #endif
21 #endif
22 
23 #include "./InternalHeaderCheck.h"
24 
25 namespace Eigen {
26 
27 /* TODO extract L, extract U, compute det, etc... */
28 
29 // generic double/complex<double> wrapper functions:
30 
31 
32  // Defaults
33 inline void umfpack_defaults(double control[UMFPACK_CONTROL], double, int)
34 { umfpack_di_defaults(control); }
35 
36 inline void umfpack_defaults(double control[UMFPACK_CONTROL], std::complex<double>, int)
37 { umfpack_zi_defaults(control); }
38 
39 inline void umfpack_defaults(double control[UMFPACK_CONTROL], double, SuiteSparse_long)
40 { umfpack_dl_defaults(control); }
41 
42 inline void umfpack_defaults(double control[UMFPACK_CONTROL], std::complex<double>, SuiteSparse_long)
43 { umfpack_zl_defaults(control); }
44 
45 // Report info
46 inline void umfpack_report_info(double control[UMFPACK_CONTROL], double info[UMFPACK_INFO], double, int)
47 { umfpack_di_report_info(control, info);}
48 
49 inline void umfpack_report_info(double control[UMFPACK_CONTROL], double info[UMFPACK_INFO], std::complex<double>, int)
50 { umfpack_zi_report_info(control, info);}
51 
52 inline void umfpack_report_info(double control[UMFPACK_CONTROL], double info[UMFPACK_INFO], double, SuiteSparse_long)
53 { umfpack_dl_report_info(control, info);}
54 
55 inline void umfpack_report_info(double control[UMFPACK_CONTROL], double info[UMFPACK_INFO], std::complex<double>, SuiteSparse_long)
56 { umfpack_zl_report_info(control, info);}
57 
58 // Report status
59 inline void umfpack_report_status(double control[UMFPACK_CONTROL], int status, double, int)
60 { umfpack_di_report_status(control, status);}
61 
62 inline void umfpack_report_status(double control[UMFPACK_CONTROL], int status, std::complex<double>, int)
63 { umfpack_zi_report_status(control, status);}
64 
65 inline void umfpack_report_status(double control[UMFPACK_CONTROL], int status, double, SuiteSparse_long)
66 { umfpack_dl_report_status(control, status);}
67 
68 inline void umfpack_report_status(double control[UMFPACK_CONTROL], int status, std::complex<double>, SuiteSparse_long)
69 { umfpack_zl_report_status(control, status);}
70 
71 // report control
72 inline void umfpack_report_control(double control[UMFPACK_CONTROL], double, int)
73 { umfpack_di_report_control(control);}
74 
75 inline void umfpack_report_control(double control[UMFPACK_CONTROL], std::complex<double>, int)
76 { umfpack_zi_report_control(control);}
77 
78 inline void umfpack_report_control(double control[UMFPACK_CONTROL], double, SuiteSparse_long)
79 { umfpack_dl_report_control(control);}
80 
81 inline void umfpack_report_control(double control[UMFPACK_CONTROL], std::complex<double>, SuiteSparse_long)
82 { umfpack_zl_report_control(control);}
83 
84 // Free numeric
85 inline void umfpack_free_numeric(void **Numeric, double, int)
86 { umfpack_di_free_numeric(Numeric); *Numeric = 0; }
87 
88 inline void umfpack_free_numeric(void **Numeric, std::complex<double>, int)
89 { umfpack_zi_free_numeric(Numeric); *Numeric = 0; }
90 
91 inline void umfpack_free_numeric(void **Numeric, double, SuiteSparse_long)
92 { umfpack_dl_free_numeric(Numeric); *Numeric = 0; }
93 
94 inline void umfpack_free_numeric(void **Numeric, std::complex<double>, SuiteSparse_long)
95 { umfpack_zl_free_numeric(Numeric); *Numeric = 0; }
96 
97 // Free symbolic
98 inline void umfpack_free_symbolic(void **Symbolic, double, int)
99 { umfpack_di_free_symbolic(Symbolic); *Symbolic = 0; }
100 
101 inline void umfpack_free_symbolic(void **Symbolic, std::complex<double>, int)
102 { umfpack_zi_free_symbolic(Symbolic); *Symbolic = 0; }
103 
104 inline void umfpack_free_symbolic(void **Symbolic, double, SuiteSparse_long)
105 { umfpack_dl_free_symbolic(Symbolic); *Symbolic = 0; }
106 
107 inline void umfpack_free_symbolic(void **Symbolic, std::complex<double>, SuiteSparse_long)
108 { umfpack_zl_free_symbolic(Symbolic); *Symbolic = 0; }
109 
110 // Symbolic
111 inline int umfpack_symbolic(int n_row,int n_col,
112  const int Ap[], const int Ai[], const double Ax[], void **Symbolic,
113  const double Control [UMFPACK_CONTROL], double Info [UMFPACK_INFO])
114 {
115  return umfpack_di_symbolic(n_row,n_col,Ap,Ai,Ax,Symbolic,Control,Info);
116 }
117 
118 inline int umfpack_symbolic(int n_row,int n_col,
119  const int Ap[], const int Ai[], const std::complex<double> Ax[], void **Symbolic,
120  const double Control [UMFPACK_CONTROL], double Info [UMFPACK_INFO])
121 {
122  return umfpack_zi_symbolic(n_row,n_col,Ap,Ai,&numext::real_ref(Ax[0]),0,Symbolic,Control,Info);
123 }
124 inline SuiteSparse_long umfpack_symbolic( SuiteSparse_long n_row,SuiteSparse_long n_col,
125  const SuiteSparse_long Ap[], const SuiteSparse_long Ai[], const double Ax[], void **Symbolic,
126  const double Control [UMFPACK_CONTROL], double Info [UMFPACK_INFO])
127 {
128  return umfpack_dl_symbolic(n_row,n_col,Ap,Ai,Ax,Symbolic,Control,Info);
129 }
130 
131 inline SuiteSparse_long umfpack_symbolic( SuiteSparse_long n_row,SuiteSparse_long n_col,
132  const SuiteSparse_long Ap[], const SuiteSparse_long Ai[], const std::complex<double> Ax[], void **Symbolic,
133  const double Control [UMFPACK_CONTROL], double Info [UMFPACK_INFO])
134 {
135  return umfpack_zl_symbolic(n_row,n_col,Ap,Ai,&numext::real_ref(Ax[0]),0,Symbolic,Control,Info);
136 }
137 
138 // Numeric
139 inline int umfpack_numeric( const int Ap[], const int Ai[], const double Ax[],
140  void *Symbolic, void **Numeric,
141  const double Control[UMFPACK_CONTROL],double Info [UMFPACK_INFO])
142 {
143  return umfpack_di_numeric(Ap,Ai,Ax,Symbolic,Numeric,Control,Info);
144 }
145 
146 inline int umfpack_numeric( const int Ap[], const int Ai[], const std::complex<double> Ax[],
147  void *Symbolic, void **Numeric,
148  const double Control[UMFPACK_CONTROL],double Info [UMFPACK_INFO])
149 {
150  return umfpack_zi_numeric(Ap,Ai,&numext::real_ref(Ax[0]),0,Symbolic,Numeric,Control,Info);
151 }
152 inline SuiteSparse_long umfpack_numeric(const SuiteSparse_long Ap[], const SuiteSparse_long Ai[], const double Ax[],
153  void *Symbolic, void **Numeric,
154  const double Control[UMFPACK_CONTROL],double Info [UMFPACK_INFO])
155 {
156  return umfpack_dl_numeric(Ap,Ai,Ax,Symbolic,Numeric,Control,Info);
157 }
158 
159 inline SuiteSparse_long umfpack_numeric(const SuiteSparse_long Ap[], const SuiteSparse_long Ai[], const std::complex<double> Ax[],
160  void *Symbolic, void **Numeric,
161  const double Control[UMFPACK_CONTROL],double Info [UMFPACK_INFO])
162 {
163  return umfpack_zl_numeric(Ap,Ai,&numext::real_ref(Ax[0]),0,Symbolic,Numeric,Control,Info);
164 }
165 
166 // solve
167 inline int umfpack_solve( int sys, const int Ap[], const int Ai[], const double Ax[],
168  double X[], const double B[], void *Numeric,
169  const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
170 {
171  return umfpack_di_solve(sys,Ap,Ai,Ax,X,B,Numeric,Control,Info);
172 }
173 
174 inline int umfpack_solve( int sys, const int Ap[], const int Ai[], const std::complex<double> Ax[],
175  std::complex<double> X[], const std::complex<double> B[], void *Numeric,
176  const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
177 {
178  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);
179 }
180 
181 inline SuiteSparse_long umfpack_solve(int sys, const SuiteSparse_long Ap[], const SuiteSparse_long Ai[], const double Ax[],
182  double X[], const double B[], void *Numeric,
183  const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
184 {
185  return umfpack_dl_solve(sys,Ap,Ai,Ax,X,B,Numeric,Control,Info);
186 }
187 
188 inline SuiteSparse_long umfpack_solve(int sys, const SuiteSparse_long Ap[], const SuiteSparse_long Ai[], const std::complex<double> Ax[],
189  std::complex<double> X[], const std::complex<double> B[], void *Numeric,
190  const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
191 {
192  return umfpack_zl_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);
193 }
194 
195 // Get Lunz
196 inline int umfpack_get_lunz(int *lnz, int *unz, int *n_row, int *n_col, int *nz_udiag, void *Numeric, double)
197 {
198  return umfpack_di_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
199 }
200 
201 inline int umfpack_get_lunz(int *lnz, int *unz, int *n_row, int *n_col, int *nz_udiag, void *Numeric, std::complex<double>)
202 {
203  return umfpack_zi_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
204 }
205 
206 inline SuiteSparse_long umfpack_get_lunz( SuiteSparse_long *lnz, SuiteSparse_long *unz, SuiteSparse_long *n_row, SuiteSparse_long *n_col,
207  SuiteSparse_long *nz_udiag, void *Numeric, double)
208 {
209  return umfpack_dl_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
210 }
211 
212 inline SuiteSparse_long umfpack_get_lunz( SuiteSparse_long *lnz, SuiteSparse_long *unz, SuiteSparse_long *n_row, SuiteSparse_long *n_col,
213  SuiteSparse_long *nz_udiag, void *Numeric, std::complex<double>)
214 {
215  return umfpack_zl_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
216 }
217 
218 // Get Numeric
219 inline int umfpack_get_numeric(int Lp[], int Lj[], double Lx[], int Up[], int Ui[], double Ux[],
220  int P[], int Q[], double Dx[], int *do_recip, double Rs[], void *Numeric)
221 {
222  return umfpack_di_get_numeric(Lp,Lj,Lx,Up,Ui,Ux,P,Q,Dx,do_recip,Rs,Numeric);
223 }
224 
225 inline int umfpack_get_numeric(int Lp[], int Lj[], std::complex<double> Lx[], int Up[], int Ui[], std::complex<double> Ux[],
226  int P[], int Q[], std::complex<double> Dx[], int *do_recip, double Rs[], void *Numeric)
227 {
228  double& lx0_real = numext::real_ref(Lx[0]);
229  double& ux0_real = numext::real_ref(Ux[0]);
230  double& dx0_real = numext::real_ref(Dx[0]);
231  return umfpack_zi_get_numeric(Lp,Lj,Lx?&lx0_real:0,0,Up,Ui,Ux?&ux0_real:0,0,P,Q,
232  Dx?&dx0_real:0,0,do_recip,Rs,Numeric);
233 }
234 inline SuiteSparse_long umfpack_get_numeric(SuiteSparse_long Lp[], SuiteSparse_long Lj[], double Lx[], SuiteSparse_long Up[], SuiteSparse_long Ui[], double Ux[],
235  SuiteSparse_long P[], SuiteSparse_long Q[], double Dx[], SuiteSparse_long *do_recip, double Rs[], void *Numeric)
236 {
237  return umfpack_dl_get_numeric(Lp,Lj,Lx,Up,Ui,Ux,P,Q,Dx,do_recip,Rs,Numeric);
238 }
239 
240 inline SuiteSparse_long umfpack_get_numeric(SuiteSparse_long Lp[], SuiteSparse_long Lj[], std::complex<double> Lx[], SuiteSparse_long Up[], SuiteSparse_long Ui[], std::complex<double> Ux[],
241  SuiteSparse_long P[], SuiteSparse_long Q[], std::complex<double> Dx[], SuiteSparse_long *do_recip, double Rs[], void *Numeric)
242 {
243  double& lx0_real = numext::real_ref(Lx[0]);
244  double& ux0_real = numext::real_ref(Ux[0]);
245  double& dx0_real = numext::real_ref(Dx[0]);
246  return umfpack_zl_get_numeric(Lp,Lj,Lx?&lx0_real:0,0,Up,Ui,Ux?&ux0_real:0,0,P,Q,
247  Dx?&dx0_real:0,0,do_recip,Rs,Numeric);
248 }
249 
250 // Get Determinant
251 inline int umfpack_get_determinant(double *Mx, double *Ex, void *NumericHandle, double User_Info [UMFPACK_INFO], int)
252 {
253  return umfpack_di_get_determinant(Mx,Ex,NumericHandle,User_Info);
254 }
255 
256 inline int umfpack_get_determinant(std::complex<double> *Mx, double *Ex, void *NumericHandle, double User_Info [UMFPACK_INFO], int)
257 {
258  double& mx_real = numext::real_ref(*Mx);
259  return umfpack_zi_get_determinant(&mx_real,0,Ex,NumericHandle,User_Info);
260 }
261 
262 inline SuiteSparse_long umfpack_get_determinant(double *Mx, double *Ex, void *NumericHandle, double User_Info [UMFPACK_INFO], SuiteSparse_long)
263 {
264  return umfpack_dl_get_determinant(Mx,Ex,NumericHandle,User_Info);
265 }
266 
267 inline SuiteSparse_long umfpack_get_determinant(std::complex<double> *Mx, double *Ex, void *NumericHandle, double User_Info [UMFPACK_INFO], SuiteSparse_long)
268 {
269  double& mx_real = numext::real_ref(*Mx);
270  return umfpack_zl_get_determinant(&mx_real,0,Ex,NumericHandle,User_Info);
271 }
272 
273 
289 template<typename MatrixType_>
290 class UmfPackLU : public SparseSolverBase<UmfPackLU<MatrixType_> >
291 {
292  protected:
294  using Base::m_isInitialized;
295  public:
296  using Base::_solve_impl;
297  typedef MatrixType_ MatrixType;
298  typedef typename MatrixType::Scalar Scalar;
307  enum {
308  ColsAtCompileTime = MatrixType::ColsAtCompileTime,
309  MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
310  };
311 
312  public:
313 
316 
318  : m_dummy(0,0), mp_matrix(m_dummy)
319  {
320  init();
321  }
322 
323  template<typename InputMatrixType>
324  explicit UmfPackLU(const InputMatrixType& matrix)
325  : mp_matrix(matrix)
326  {
327  init();
328  compute(matrix);
329  }
330 
332  {
335  }
336 
337  inline Index rows() const { return mp_matrix.rows(); }
338  inline Index cols() const { return mp_matrix.cols(); }
339 
346  {
347  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
348  return m_info;
349  }
350 
351  inline const LUMatrixType& matrixL() const
352  {
354  return m_l;
355  }
356 
357  inline const LUMatrixType& matrixU() const
358  {
360  return m_u;
361  }
362 
363  inline const IntColVectorType& permutationP() const
364  {
366  return m_p;
367  }
368 
369  inline const IntRowVectorType& permutationQ() const
370  {
372  return m_q;
373  }
374 
379  template<typename InputMatrixType>
380  void compute(const InputMatrixType& matrix)
381  {
384  grab(matrix.derived());
386  factorize_impl();
387  }
388 
395  template<typename InputMatrixType>
396  void analyzePattern(const InputMatrixType& matrix)
397  {
400 
401  grab(matrix.derived());
402 
404  }
405 
411  inline int umfpackFactorizeReturncode() const
412  {
413  eigen_assert(m_numeric && "UmfPackLU: you must first call factorize()");
414  return m_fact_errorCode;
415  }
416 
423  inline const UmfpackControl& umfpackControl() const
424  {
425  return m_control;
426  }
427 
435  {
436  return m_control;
437  }
438 
445  template<typename InputMatrixType>
446  void factorize(const InputMatrixType& matrix)
447  {
448  eigen_assert(m_analysisIsOk && "UmfPackLU: you must first call analyzePattern()");
449  if(m_numeric)
451 
452  grab(matrix.derived());
453 
454  factorize_impl();
455  }
456 
462  {
464  }
465 
471  {
472  eigen_assert(m_analysisIsOk && "UmfPackLU: you must first call analyzePattern()");
474  }
475 
481  eigen_assert(m_analysisIsOk && "UmfPackLU: you must first call analyzePattern()");
483  }
484 
486  template<typename BDerived,typename XDerived>
488 
489  Scalar determinant() const;
490 
491  void extractData() const;
492 
493  protected:
494 
495  void init()
496  {
498  m_isInitialized = false;
499  m_numeric = 0;
500  m_symbolic = 0;
502 
504  }
505 
507  {
508  m_fact_errorCode = umfpack_symbolic(internal::convert_index<StorageIndex>(mp_matrix.rows()),
509  internal::convert_index<StorageIndex>(mp_matrix.cols()),
510  mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
512 
513  m_isInitialized = true;
515  m_analysisIsOk = true;
516  m_factorizationIsOk = false;
518  }
519 
521  {
522 
523  m_fact_errorCode = umfpack_numeric(mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
525 
526  m_info = m_fact_errorCode == UMFPACK_OK ? Success : NumericalIssue;
527  m_factorizationIsOk = true;
529  }
530 
531  template<typename MatrixDerived>
533  {
535  internal::construct_at(&mp_matrix, A.derived());
536  }
537 
538  void grab(const UmfpackMatrixRef &A)
539  {
540  if(&(A.derived()) != &mp_matrix)
541  {
544  }
545  }
546 
547  // cached data to reduce reallocation, etc.
548  mutable LUMatrixType m_l;
552 
553  mutable LUMatrixType m_u;
556 
559 
560  void* m_numeric;
561  void* m_symbolic;
562 
567 
568  private:
569  UmfPackLU(const UmfPackLU& ) { }
570 };
571 
572 
573 template<typename MatrixType>
575 {
576  if (m_extractedDataAreDirty)
577  {
578  // get size of the data
579  StorageIndex lnz, unz, rows, cols, nz_udiag;
580  umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());
581 
582  // allocate data
583  m_l.resize(rows,(std::min)(rows,cols));
584  m_l.resizeNonZeros(lnz);
585 
586  m_u.resize((std::min)(rows,cols),cols);
587  m_u.resizeNonZeros(unz);
588 
589  m_p.resize(rows);
590  m_q.resize(cols);
591 
592  // extract
593  umfpack_get_numeric(m_l.outerIndexPtr(), m_l.innerIndexPtr(), m_l.valuePtr(),
594  m_u.outerIndexPtr(), m_u.innerIndexPtr(), m_u.valuePtr(),
595  m_p.data(), m_q.data(), 0, 0, 0, m_numeric);
596 
597  m_extractedDataAreDirty = false;
598  }
599 }
600 
601 template<typename MatrixType>
603 {
604  Scalar det;
605  umfpack_get_determinant(&det, 0, m_numeric, 0, StorageIndex());
606  return det;
607 }
608 
609 template<typename MatrixType>
610 template<typename BDerived,typename XDerived>
612 {
613  Index rhsCols = b.cols();
614  eigen_assert((BDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major rhs yet");
615  eigen_assert((XDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major result yet");
616  eigen_assert(b.derived().data() != x.derived().data() && " Umfpack does not support inplace solve");
617 
618  Scalar* x_ptr = 0;
620  if(x.innerStride()!=1)
621  {
622  x_tmp.resize(x.rows());
623  x_ptr = x_tmp.data();
624  }
625  for (int j=0; j<rhsCols; ++j)
626  {
627  if(x.innerStride()==1)
628  x_ptr = &x.col(j).coeffRef(0);
629  StorageIndex errorCode = umfpack_solve(UMFPACK_A,
630  mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
631  x_ptr, &b.const_cast_derived().col(j).coeffRef(0),
632  m_numeric, m_control.data(), m_umfpackInfo.data());
633  if(x.innerStride()!=1)
634  x.col(j) = x_tmp;
635  if (errorCode!=0)
636  return false;
637  }
638 
639  return true;
640 }
641 
642 } // end namespace Eigen
643 
644 #endif // EIGEN_UMFPACKSUPPORT_H
Array< int, 3, 1 > b
MatrixXcf A
Projective3d P(Matrix4d::Random())
MatrixXf B
MatrixXf Q
#define eigen_assert(x)
Definition: Macros.h:902
NumTraits< Scalar >::Real RealScalar
Definition: DenseBase.h:68
internal::traits< Derived >::StorageIndex StorageIndex
The type used to store indices.
Definition: DenseBase.h:58
internal::traits< Derived >::Scalar Scalar
Definition: DenseBase.h:61
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
const Scalar * data() const
constexpr void resize(Index rows, Index cols)
A base class for sparse solvers.
A sparse LU factorization and solver based on UmfPack.
Index cols() const
void analyzePattern_impl()
UmfPackLU(const UmfPackLU &)
const LUMatrixType & matrixL() const
Index rows() const
Array< double, UMFPACK_CONTROL, 1 > UmfpackControl
UmfpackMatrixRef mp_matrix
Matrix< int, MatrixType::RowsAtCompileTime, 1 > IntColVectorType
void printUmfpackControl()
const IntColVectorType & permutationP() const
MatrixType::RealScalar RealScalar
void compute(const InputMatrixType &matrix)
UmfPackLU(const InputMatrixType &matrix)
UmfpackMatrixType m_dummy
UmfpackInfo m_umfpackInfo
SparseMatrix< Scalar, ColMajor, StorageIndex > UmfpackMatrixType
Scalar determinant() const
const IntRowVectorType & permutationQ() const
void factorize(const InputMatrixType &matrix)
UmfpackControl & umfpackControl()
ComputationInfo m_info
MatrixType::Scalar Scalar
void grab(const UmfpackMatrixRef &A)
SparseSolverBase< UmfPackLU< MatrixType_ > > Base
StorageIndex m_fact_errorCode
IntColVectorType m_p
void extractData() const
SparseMatrix< Scalar > LUMatrixType
int umfpackFactorizeReturncode() const
bool _solve_impl(const MatrixBase< BDerived > &b, MatrixBase< XDerived > &x) const
Matrix< int, 1, MatrixType::ColsAtCompileTime > IntRowVectorType
MatrixType_ MatrixType
ComputationInfo info() const
Reports whether previous computation was successful.
Array< double, UMFPACK_INFO, 1 > UmfpackInfo
UmfpackControl m_control
void grab(const EigenBase< MatrixDerived > &A)
const LUMatrixType & matrixU() const
Ref< const UmfpackMatrixType, StandardCompressedFormat > UmfpackMatrixRef
LUMatrixType m_l
void analyzePattern(const InputMatrixType &matrix)
const UmfpackControl & umfpackControl() const
Matrix< Scalar, Dynamic, 1 > Vector
IntRowVectorType m_q
LUMatrixType m_u
MatrixType::StorageIndex StorageIndex
ComputationInfo
Definition: Constants.h:444
@ NumericalIssue
Definition: Constants.h:448
@ InvalidInput
Definition: Constants.h:453
@ Success
Definition: Constants.h:446
const unsigned int RowMajorBit
Definition: Constants.h:68
bfloat16() min(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:684
void destroy_at(T *p)
Definition: Memory.h:1264
T * construct_at(T *p, Args &&... args)
Definition: Memory.h:1248
internal::add_const_on_value_type_t< EIGEN_MATHFUNC_RETVAL(real_ref, Scalar) > real_ref(const Scalar &x)
: InteropHeaders
Definition: Core:139
int umfpack_get_numeric(int Lp[], int Lj[], double Lx[], int Up[], int Ui[], double Ux[], int P[], int Q[], double Dx[], int *do_recip, double Rs[], void *Numeric)
int umfpack_get_determinant(double *Mx, double *Ex, void *NumericHandle, double User_Info[UMFPACK_INFO], int)
int umfpack_solve(int sys, const int Ap[], const int Ai[], const double Ax[], double X[], const double B[], void *Numeric, const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
void umfpack_report_control(double control[UMFPACK_CONTROL], double, int)
void umfpack_free_numeric(void **Numeric, double, int)
void umfpack_defaults(double control[UMFPACK_CONTROL], double, int)
int umfpack_symbolic(int n_row, int n_col, const int Ap[], const int Ai[], const double Ax[], void **Symbolic, const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
void umfpack_report_status(double control[UMFPACK_CONTROL], int status, double, int)
void umfpack_report_info(double control[UMFPACK_CONTROL], double info[UMFPACK_INFO], double, int)
void umfpack_free_symbolic(void **Symbolic, double, int)
int umfpack_numeric(const int Ap[], const int Ai[], const double Ax[], void *Symbolic, void **Numeric, const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
int umfpack_get_lunz(int *lnz, int *unz, int *n_row, int *n_col, int *nz_udiag, void *Numeric, double)
std::ptrdiff_t j