10#ifndef EIGEN_UMFPACKSUPPORT_H
11#define EIGEN_UMFPACKSUPPORT_H
20inline void umfpack_defaults(
double control[UMFPACK_CONTROL],
double)
21{ umfpack_di_defaults(control); }
23inline void umfpack_defaults(
double control[UMFPACK_CONTROL], std::complex<double>)
24{ umfpack_zi_defaults(control); }
26inline void umfpack_free_numeric(
void **Numeric,
double)
27{ umfpack_di_free_numeric(Numeric); *Numeric = 0; }
29inline void umfpack_free_numeric(
void **Numeric, std::complex<double>)
30{ umfpack_zi_free_numeric(Numeric); *Numeric = 0; }
32inline void umfpack_free_symbolic(
void **Symbolic,
double)
33{ umfpack_di_free_symbolic(Symbolic); *Symbolic = 0; }
35inline void umfpack_free_symbolic(
void **Symbolic, std::complex<double>)
36{ umfpack_zi_free_symbolic(Symbolic); *Symbolic = 0; }
38inline 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);
45inline 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);
52inline 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);
59inline 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);
66inline 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);
73inline 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);
80inline 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);
85inline 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);
90inline 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);
96inline 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);
106inline 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);
111inline 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);
131template<
typename _MatrixType>
136 using Base::m_isInitialized;
138 using Base::_solve_impl;
139 typedef _MatrixType MatrixType;
140 typedef typename MatrixType::Scalar Scalar;
141 typedef typename MatrixType::RealScalar RealScalar;
142 typedef typename MatrixType::StorageIndex StorageIndex;
150 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
151 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
159 : m_dummy(0,0), mp_matrix(m_dummy)
164 template<
typename InputMatrixType>
174 if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
175 if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
178 inline Index rows()
const {
return mp_matrix.rows(); }
179 inline Index cols()
const {
return mp_matrix.cols(); }
188 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
192 inline const LUMatrixType& matrixL()
const
194 if (m_extractedDataAreDirty) extractData();
198 inline const LUMatrixType& matrixU()
const
200 if (m_extractedDataAreDirty) extractData();
204 inline const IntColVectorType& permutationP()
const
206 if (m_extractedDataAreDirty) extractData();
210 inline const IntRowVectorType& permutationQ()
const
212 if (m_extractedDataAreDirty) extractData();
220 template<
typename InputMatrixType>
223 if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
224 if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
225 grab(matrix.derived());
226 analyzePattern_impl();
236 template<
typename InputMatrixType>
239 if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
240 if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
242 grab(matrix.derived());
244 analyzePattern_impl();
254 eigen_assert(m_numeric &&
"UmfPackLU: you must first call factorize()");
255 return m_fact_errorCode;
286 template<
typename InputMatrixType>
289 eigen_assert(m_analysisIsOk &&
"UmfPackLU: you must first call analyzePattern()");
291 umfpack_free_numeric(&m_numeric,Scalar());
293 grab(matrix.derived());
299 template<
typename BDerived,
typename XDerived>
302 Scalar determinant()
const;
304 void extractData()
const;
311 m_isInitialized =
false;
314 m_extractedDataAreDirty =
true;
317 void analyzePattern_impl()
319 umfpack_defaults(m_control.
data(), Scalar());
321 errorCode = umfpack_symbolic(internal::convert_index<int>(mp_matrix.rows()),
322 internal::convert_index<int>(mp_matrix.cols()),
323 mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
324 &m_symbolic, m_control.
data(), 0);
326 m_isInitialized =
true;
328 m_analysisIsOk =
true;
329 m_factorizationIsOk =
false;
330 m_extractedDataAreDirty =
true;
333 void factorize_impl()
335 m_fact_errorCode = umfpack_numeric(mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
336 m_symbolic, &m_numeric, m_control.
data(), 0);
339 m_factorizationIsOk =
true;
340 m_extractedDataAreDirty =
true;
343 template<
typename MatrixDerived>
344 void grab(
const EigenBase<MatrixDerived> &A)
346 mp_matrix.~UmfpackMatrixRef();
347 ::new (&mp_matrix) UmfpackMatrixRef(A.derived());
350 void grab(
const UmfpackMatrixRef &A)
352 if(&(A.derived()) != &mp_matrix)
354 mp_matrix.~UmfpackMatrixRef();
355 ::new (&mp_matrix) UmfpackMatrixRef(A);
360 mutable LUMatrixType m_l;
361 int m_fact_errorCode;
362 UmfpackControl m_control;
364 mutable LUMatrixType m_u;
365 mutable IntColVectorType m_p;
366 mutable IntRowVectorType m_q;
368 UmfpackMatrixType m_dummy;
369 UmfpackMatrixRef mp_matrix;
375 int m_factorizationIsOk;
377 mutable bool m_extractedDataAreDirty;
380 UmfPackLU(UmfPackLU& ) { }
384template<
typename MatrixType>
385void UmfPackLU<MatrixType>::extractData()
const
387 if (m_extractedDataAreDirty)
390 int lnz, unz, rows, cols, nz_udiag;
391 umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());
394 m_l.resize(rows,(std::min)(rows,cols));
395 m_l.resizeNonZeros(lnz);
397 m_u.resize((std::min)(rows,cols),cols);
398 m_u.resizeNonZeros(unz);
404 umfpack_get_numeric(m_l.outerIndexPtr(), m_l.innerIndexPtr(), m_l.valuePtr(),
405 m_u.outerIndexPtr(), m_u.innerIndexPtr(), m_u.valuePtr(),
406 m_p.data(), m_q.data(), 0, 0, 0, m_numeric);
408 m_extractedDataAreDirty =
false;
412template<
typename MatrixType>
413typename UmfPackLU<MatrixType>::Scalar UmfPackLU<MatrixType>::determinant()
const
416 umfpack_get_determinant(&det, 0, m_numeric, 0);
420template<
typename MatrixType>
421template<
typename BDerived,
typename XDerived>
422bool UmfPackLU<MatrixType>::_solve_impl(
const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x)
const
424 Index rhsCols = b.cols();
425 eigen_assert((BDerived::Flags&
RowMajorBit)==0 &&
"UmfPackLU backend does not support non col-major rhs yet");
426 eigen_assert((XDerived::Flags&
RowMajorBit)==0 &&
"UmfPackLU backend does not support non col-major result yet");
427 eigen_assert(b.derived().data() != x.derived().data() &&
" Umfpack does not support inplace solve");
432 if(x.innerStride()!=1)
434 x_tmp.resize(x.rows());
435 x_ptr = x_tmp.data();
437 for (
int j=0; j<rhsCols; ++j)
439 if(x.innerStride()==1)
440 x_ptr = &x.col(j).coeffRef(0);
441 errorCode = umfpack_solve(UMFPACK_A,
442 mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
443 x_ptr, &b.const_cast_derived().col(j).coeffRef(0), m_numeric, m_control.data(), 0);
444 if(x.innerStride()!=1)
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar * data() const
Definition PlainObjectBase.h:228
Pseudo expression representing a solving operation.
Definition Solve.h:63
A base class for sparse solvers.
Definition SparseSolverBase.h:54
A sparse LU factorization and solver based on UmfPack.
Definition UmfPackSupport.h:133
void compute(const InputMatrixType &matrix)
Computes the sparse Cholesky decomposition of matrix Note that the matrix should be column-major,...
Definition UmfPackSupport.h:221
void factorize(const InputMatrixType &matrix)
Performs a numeric decomposition of matrix.
Definition UmfPackSupport.h:287
const UmfpackControl & umfpackControl() const
Provides access to the control settings array used by UmfPack.
Definition UmfPackSupport.h:264
ComputationInfo info() const
Reports whether previous computation was successful.
Definition UmfPackSupport.h:186
UmfpackControl & umfpackControl()
Provides access to the control settings array used by UmfPack.
Definition UmfPackSupport.h:275
int umfpackFactorizeReturncode() const
Provides the return status code returned by UmfPack during the numeric factorization.
Definition UmfPackSupport.h:252
void analyzePattern(const InputMatrixType &matrix)
Performs a symbolic decomposition on the sparcity of matrix.
Definition UmfPackSupport.h:237
ComputationInfo
Enum for reporting the status of a computation.
Definition Constants.h:430
@ NumericalIssue
The provided data did not satisfy the prerequisites.
Definition Constants.h:434
@ InvalidInput
The inputs are invalid, or the algorithm has been improperly called.
Definition Constants.h:439
@ Success
Computation was successful.
Definition Constants.h:432
const unsigned int RowMajorBit
for a matrix, this means that the storage order is row-major.
Definition Constants.h:61
Definition inference.c:32