10 #ifndef EIGEN_UMFPACKSUPPORT_H 11 #define EIGEN_UMFPACKSUPPORT_H 20 inline void umfpack_defaults(
double control[UMFPACK_CONTROL],
double)
21 { umfpack_di_defaults(control); }
23 inline void umfpack_defaults(
double control[UMFPACK_CONTROL], std::complex<double>)
24 { umfpack_zi_defaults(control); }
26 inline void umfpack_free_numeric(
void **Numeric,
double)
27 { umfpack_di_free_numeric(Numeric); *Numeric = 0; }
29 inline void umfpack_free_numeric(
void **Numeric, std::complex<double>)
30 { umfpack_zi_free_numeric(Numeric); *Numeric = 0; }
32 inline void umfpack_free_symbolic(
void **Symbolic,
double)
33 { umfpack_di_free_symbolic(Symbolic); *Symbolic = 0; }
35 inline void umfpack_free_symbolic(
void **Symbolic, std::complex<double>)
36 { umfpack_zi_free_symbolic(Symbolic); *Symbolic = 0; }
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])
42 return umfpack_di_symbolic(n_row,n_col,Ap,Ai,Ax,Symbolic,Control,Info);
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])
49 return umfpack_zi_symbolic(n_row,n_col,Ap,Ai,&numext::real_ref(Ax[0]),0,Symbolic,Control,Info);
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])
56 return umfpack_di_numeric(Ap,Ai,Ax,Symbolic,Numeric,Control,Info);
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])
63 return umfpack_zi_numeric(Ap,Ai,&numext::real_ref(Ax[0]),0,Symbolic,Numeric,Control,Info);
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])
70 return umfpack_di_solve(sys,Ap,Ai,Ax,X,B,Numeric,Control,Info);
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])
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);
80 inline int umfpack_get_lunz(
int *lnz,
int *unz,
int *n_row,
int *n_col,
int *nz_udiag,
void *Numeric,
double)
82 return umfpack_di_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
85 inline int umfpack_get_lunz(
int *lnz,
int *unz,
int *n_row,
int *n_col,
int *nz_udiag,
void *Numeric, std::complex<double>)
87 return umfpack_zi_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
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)
93 return umfpack_di_get_numeric(Lp,Lj,Lx,Up,Ui,Ux,P,Q,Dx,do_recip,Rs,Numeric);
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)
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);
106 inline int umfpack_get_determinant(
double *Mx,
double *Ex,
void *NumericHandle,
double User_Info [UMFPACK_INFO])
108 return umfpack_di_get_determinant(Mx,Ex,NumericHandle,User_Info);
111 inline int umfpack_get_determinant(std::complex<double> *Mx,
double *Ex,
void *NumericHandle,
double User_Info [UMFPACK_INFO])
113 double& mx_real = numext::real_ref(*Mx);
114 return umfpack_zi_get_determinant(&mx_real,0,Ex,NumericHandle,User_Info);
133 template<
typename _MatrixType>
138 using Base::m_isInitialized;
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;
152 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
153 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
161 : m_dummy(0,0), mp_matrix(m_dummy)
166 template<
typename InputMatrixType>
167 explicit UmfPackLU(
const InputMatrixType& matrix)
176 if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
177 if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
180 inline Index rows()
const {
return mp_matrix.rows(); }
181 inline Index cols()
const {
return mp_matrix.cols(); }
190 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
194 inline const LUMatrixType& matrixL()
const 196 if (m_extractedDataAreDirty) extractData();
200 inline const LUMatrixType& matrixU()
const 202 if (m_extractedDataAreDirty) extractData();
206 inline const IntColVectorType& permutationP()
const 208 if (m_extractedDataAreDirty) extractData();
212 inline const IntRowVectorType& permutationQ()
const 214 if (m_extractedDataAreDirty) extractData();
222 template<
typename InputMatrixType>
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();
238 template<
typename InputMatrixType>
241 if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
242 if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
244 grab(matrix.derived());
246 analyzePattern_impl();
256 eigen_assert(m_numeric &&
"UmfPackLU: you must first call factorize()");
257 return m_fact_errorCode;
288 template<
typename InputMatrixType>
291 eigen_assert(m_analysisIsOk &&
"UmfPackLU: you must first call analyzePattern()");
293 umfpack_free_numeric(&m_numeric,Scalar());
295 grab(matrix.derived());
301 template<
typename BDerived,
typename XDerived>
304 Scalar determinant()
const;
306 void extractData()
const;
313 m_isInitialized =
false;
316 m_extractedDataAreDirty =
true;
319 void analyzePattern_impl()
321 umfpack_defaults(m_control.
data(), Scalar());
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);
328 m_isInitialized =
true;
330 m_analysisIsOk =
true;
331 m_factorizationIsOk =
false;
332 m_extractedDataAreDirty =
true;
335 void factorize_impl()
337 m_fact_errorCode = umfpack_numeric(mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
338 m_symbolic, &m_numeric, m_control.
data(), 0);
341 m_factorizationIsOk =
true;
342 m_extractedDataAreDirty =
true;
345 template<
typename MatrixDerived>
348 mp_matrix.~UmfpackMatrixRef();
349 ::new (&mp_matrix) UmfpackMatrixRef(A.
derived());
352 void grab(
const UmfpackMatrixRef &A)
354 if(&(A.derived()) != &mp_matrix)
356 mp_matrix.~UmfpackMatrixRef();
357 ::new (&mp_matrix) UmfpackMatrixRef(A);
362 mutable LUMatrixType m_l;
363 int m_fact_errorCode;
364 UmfpackControl m_control;
366 mutable LUMatrixType m_u;
367 mutable IntColVectorType m_p;
368 mutable IntRowVectorType m_q;
370 UmfpackMatrixType m_dummy;
371 UmfpackMatrixRef mp_matrix;
377 int m_factorizationIsOk;
379 mutable bool m_extractedDataAreDirty;
386 template<
typename MatrixType>
389 if (m_extractedDataAreDirty)
392 int lnz, unz, rows, cols, nz_udiag;
393 umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());
396 m_l.
resize(rows,(std::min)(rows,cols));
397 m_l.resizeNonZeros(lnz);
399 m_u.
resize((std::min)(rows,cols),cols);
400 m_u.resizeNonZeros(unz);
408 m_p.
data(), m_q.
data(), 0, 0, 0, m_numeric);
410 m_extractedDataAreDirty =
false;
414 template<
typename MatrixType>
418 umfpack_get_determinant(&det, 0, m_numeric, 0);
422 template<
typename MatrixType>
423 template<
typename BDerived,
typename XDerived>
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");
434 if(x.innerStride()!=1)
437 x_ptr = x_tmp.
data();
439 for (
int j=0; j<rhsCols; ++j)
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)
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