10 #ifndef EIGEN_UMFPACKSUPPORT_H
11 #define EIGEN_UMFPACKSUPPORT_H
19 inline void umfpack_free_numeric(
void **Numeric,
double)
20 { umfpack_di_free_numeric(Numeric); *Numeric = 0; }
22 inline void umfpack_free_numeric(
void **Numeric, std::complex<double>)
23 { umfpack_zi_free_numeric(Numeric); *Numeric = 0; }
25 inline void umfpack_free_symbolic(
void **Symbolic,
double)
26 { umfpack_di_free_symbolic(Symbolic); *Symbolic = 0; }
28 inline void umfpack_free_symbolic(
void **Symbolic, std::complex<double>)
29 { umfpack_zi_free_symbolic(Symbolic); *Symbolic = 0; }
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])
35 return umfpack_di_symbolic(n_row,n_col,Ap,Ai,Ax,Symbolic,Control,Info);
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])
42 return umfpack_zi_symbolic(n_row,n_col,Ap,Ai,&numext::real_ref(Ax[0]),0,Symbolic,Control,Info);
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])
49 return umfpack_di_numeric(Ap,Ai,Ax,Symbolic,Numeric,Control,Info);
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])
56 return umfpack_zi_numeric(Ap,Ai,&numext::real_ref(Ax[0]),0,Symbolic,Numeric,Control,Info);
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])
63 return umfpack_di_solve(sys,Ap,Ai,Ax,X,B,Numeric,Control,Info);
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])
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);
73 inline int umfpack_get_lunz(
int *lnz,
int *unz,
int *n_row,
int *n_col,
int *nz_udiag,
void *Numeric,
double)
75 return umfpack_di_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
78 inline int umfpack_get_lunz(
int *lnz,
int *unz,
int *n_row,
int *n_col,
int *nz_udiag,
void *Numeric, std::complex<double>)
80 return umfpack_zi_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
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)
86 return umfpack_di_get_numeric(Lp,Lj,Lx,Up,Ui,Ux,P,Q,Dx,do_recip,Rs,Numeric);
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)
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);
99 inline int umfpack_get_determinant(
double *Mx,
double *Ex,
void *NumericHandle,
double User_Info [UMFPACK_INFO])
101 return umfpack_di_get_determinant(Mx,Ex,NumericHandle,User_Info);
104 inline int umfpack_get_determinant(std::complex<double> *Mx,
double *Ex,
void *NumericHandle,
double User_Info [UMFPACK_INFO])
106 double& mx_real = numext::real_ref(*Mx);
107 return umfpack_zi_get_determinant(&mx_real,0,Ex,NumericHandle,User_Info);
124 template<
typename _MatrixType>
129 using Base::m_isInitialized;
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;
146 : m_dummy(0,0), mp_matrix(m_dummy)
151 explicit UmfPackLU(
const MatrixType& matrix)
160 if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
161 if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
164 inline Index rows()
const {
return mp_matrix.rows(); }
165 inline Index cols()
const {
return mp_matrix.cols(); }
174 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
178 inline const LUMatrixType& matrixL()
const
180 if (m_extractedDataAreDirty) extractData();
184 inline const LUMatrixType& matrixU()
const
186 if (m_extractedDataAreDirty) extractData();
190 inline const IntColVectorType& permutationP()
const
192 if (m_extractedDataAreDirty) extractData();
196 inline const IntRowVectorType& permutationQ()
const
198 if (m_extractedDataAreDirty) extractData();
206 template<
typename InputMatrixType>
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();
222 template<
typename InputMatrixType>
225 if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
226 if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
228 grab(matrix.derived());
230 analyzePattern_impl();
239 template<
typename InputMatrixType>
242 eigen_assert(m_analysisIsOk &&
"UmfPackLU: you must first call analyzePattern()");
244 umfpack_free_numeric(&m_numeric,Scalar());
246 grab(matrix.derived());
252 template<
typename BDerived,
typename XDerived>
255 Scalar determinant()
const;
257 void extractData()
const;
264 m_isInitialized =
false;
267 m_extractedDataAreDirty =
true;
270 void analyzePattern_impl()
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(),
278 m_isInitialized =
true;
280 m_analysisIsOk =
true;
281 m_factorizationIsOk =
false;
282 m_extractedDataAreDirty =
true;
285 void factorize_impl()
288 errorCode = umfpack_numeric(mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
289 m_symbolic, &m_numeric, 0, 0);
292 m_factorizationIsOk =
true;
293 m_extractedDataAreDirty =
true;
296 template<
typename MatrixDerived>
297 void grab(
const EigenBase<MatrixDerived> &A)
299 mp_matrix.~UmfpackMatrixRef();
300 ::new (&mp_matrix) UmfpackMatrixRef(A.derived());
303 void grab(const UmfpackMatrixRef &A)
305 if(&(A.derived()) != &mp_matrix)
307 mp_matrix.~UmfpackMatrixRef();
308 ::new (&mp_matrix) UmfpackMatrixRef(A);
313 mutable LUMatrixType m_l;
314 mutable LUMatrixType m_u;
315 mutable IntColVectorType m_p;
316 mutable IntRowVectorType m_q;
318 UmfpackMatrixType m_dummy;
319 UmfpackMatrixRef mp_matrix;
325 int m_factorizationIsOk;
327 mutable
bool m_extractedDataAreDirty;
330 UmfPackLU(UmfPackLU& ) { }
334 template<
typename MatrixType>
335 void UmfPackLU<MatrixType>::extractData()
const
337 if (m_extractedDataAreDirty)
340 int lnz, unz, rows, cols, nz_udiag;
341 umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());
344 m_l.
resize(rows,(std::min)(rows,cols));
345 m_l.resizeNonZeros(lnz);
347 m_u.
resize((std::min)(rows,cols),cols);
348 m_u.resizeNonZeros(unz);
356 m_p.
data(), m_q.
data(), 0, 0, 0, m_numeric);
358 m_extractedDataAreDirty =
false;
362 template<
typename MatrixType>
363 typename UmfPackLU<MatrixType>::Scalar UmfPackLU<MatrixType>::determinant()
const
366 umfpack_get_determinant(&det, 0, m_numeric, 0);
370 template<
typename MatrixType>
371 template<
typename BDerived,
typename XDerived>
372 bool UmfPackLU<MatrixType>::_solve_impl(
const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x)
const
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");
381 Matrix<Scalar,Dynamic,1> x_tmp;
382 if(x.innerStride()!=1)
384 x_tmp.resize(x.rows());
385 x_ptr = x_tmp.data();
387 for (
int j=0; j<rhsCols; ++j)
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)
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
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