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);
123 template<
typename _MatrixType>
127 typedef _MatrixType MatrixType;
128 typedef typename MatrixType::Scalar Scalar;
129 typedef typename MatrixType::RealScalar RealScalar;
130 typedef typename MatrixType::Index Index;
149 if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
150 if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
153 inline Index rows()
const {
return m_copyMatrix.
rows(); }
154 inline Index cols()
const {
return m_copyMatrix.
cols(); }
163 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
167 inline const LUMatrixType& matrixL()
const
169 if (m_extractedDataAreDirty) extractData();
173 inline const LUMatrixType& matrixU()
const
175 if (m_extractedDataAreDirty) extractData();
179 inline const IntColVectorType& permutationP()
const
181 if (m_extractedDataAreDirty) extractData();
185 inline const IntRowVectorType& permutationQ()
const
187 if (m_extractedDataAreDirty) extractData();
205 template<
typename Rhs>
208 eigen_assert(m_isInitialized &&
"UmfPackLU is not initialized.");
209 eigen_assert(rows()==b.rows()
210 &&
"UmfPackLU::solve(): invalid number of rows of the right hand side matrix b");
211 return internal::solve_retval<UmfPackLU, Rhs>(*
this, b.derived());
218 template<
typename Rhs>
221 eigen_assert(m_isInitialized &&
"UmfPackLU is not initialized.");
222 eigen_assert(rows()==b.
rows()
223 &&
"UmfPackLU::solve(): invalid number of rows of the right hand side matrix b");
224 return internal::sparse_solve_retval<UmfPackLU, Rhs>(*
this, b.
derived());
236 umfpack_free_symbolic(&m_symbolic,Scalar());
238 umfpack_free_numeric(&m_numeric,Scalar());
243 errorCode = umfpack_symbolic(matrix.rows(), matrix.cols(), m_outerIndexPtr, m_innerIndexPtr, m_valuePtr,
246 m_isInitialized =
true;
248 m_analysisIsOk =
true;
249 m_factorizationIsOk =
false;
260 eigen_assert(m_analysisIsOk &&
"UmfPackLU: you must first call analyzePattern()");
262 umfpack_free_numeric(&m_numeric,Scalar());
267 errorCode = umfpack_numeric(m_outerIndexPtr, m_innerIndexPtr, m_valuePtr,
268 m_symbolic, &m_numeric, 0, 0);
271 m_factorizationIsOk =
true;
274 #ifndef EIGEN_PARSED_BY_DOXYGEN
276 template<
typename BDerived,
typename XDerived>
280 Scalar determinant()
const;
282 void extractData()
const;
290 m_isInitialized =
false;
298 void grapInput(
const MatrixType& mat)
300 m_copyMatrix.
resize(mat.rows(), mat.cols());
301 if( ((MatrixType::Flags&
RowMajorBit)==RowMajorBit) ||
sizeof(
typename MatrixType::Index)!=
sizeof(
int) || !mat.isCompressed() )
307 m_valuePtr = m_copyMatrix.
valuePtr();
311 m_outerIndexPtr = mat.outerIndexPtr();
312 m_innerIndexPtr = mat.innerIndexPtr();
313 m_valuePtr = mat.valuePtr();
318 mutable LUMatrixType m_l;
319 mutable LUMatrixType m_u;
320 mutable IntColVectorType m_p;
321 mutable IntRowVectorType m_q;
323 UmfpackMatrixType m_copyMatrix;
324 const Scalar* m_valuePtr;
325 const int* m_outerIndexPtr;
326 const int* m_innerIndexPtr;
331 bool m_isInitialized;
332 int m_factorizationIsOk;
334 mutable bool m_extractedDataAreDirty;
337 UmfPackLU(UmfPackLU& ) { }
341 template<
typename MatrixType>
342 void UmfPackLU<MatrixType>::extractData()
const
344 if (m_extractedDataAreDirty)
347 int lnz, unz, rows, cols, nz_udiag;
348 umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());
351 m_l.resize(rows,(std::min)(rows,cols));
352 m_l.resizeNonZeros(lnz);
354 m_u.resize((std::min)(rows,cols),cols);
355 m_u.resizeNonZeros(unz);
361 umfpack_get_numeric(m_l.outerIndexPtr(), m_l.innerIndexPtr(), m_l.valuePtr(),
362 m_u.outerIndexPtr(), m_u.innerIndexPtr(), m_u.valuePtr(),
363 m_p.data(), m_q.data(), 0, 0, 0, m_numeric);
365 m_extractedDataAreDirty =
false;
369 template<
typename MatrixType>
370 typename UmfPackLU<MatrixType>::Scalar UmfPackLU<MatrixType>::determinant()
const
373 umfpack_get_determinant(&det, 0, m_numeric, 0);
377 template<
typename MatrixType>
378 template<
typename BDerived,
typename XDerived>
379 bool UmfPackLU<MatrixType>::_solve(
const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x)
const
381 const int rhsCols = b.cols();
382 eigen_assert((BDerived::Flags&RowMajorBit)==0 &&
"UmfPackLU backend does not support non col-major rhs yet");
383 eigen_assert((XDerived::Flags&RowMajorBit)==0 &&
"UmfPackLU backend does not support non col-major result yet");
384 eigen_assert(b.derived().data() != x.derived().data() &&
" Umfpack does not support inplace solve");
387 for (
int j=0; j<rhsCols; ++j)
389 errorCode = umfpack_solve(UMFPACK_A,
390 m_outerIndexPtr, m_innerIndexPtr, m_valuePtr,
391 &x.col(j).coeffRef(0), &b.const_cast_derived().col(j).coeffRef(0), m_numeric, 0, 0);
402 template<
typename _MatrixType,
typename Rhs>
403 struct solve_retval<UmfPackLU<_MatrixType>, Rhs>
404 : solve_retval_base<UmfPackLU<_MatrixType>, Rhs>
406 typedef UmfPackLU<_MatrixType> Dec;
407 EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
409 template<typename Dest>
void evalTo(Dest& dst)
const
411 dec()._solve(rhs(),dst);
415 template<
typename _MatrixType,
typename Rhs>
416 struct sparse_solve_retval<UmfPackLU<_MatrixType>, Rhs>
417 : sparse_solve_retval_base<UmfPackLU<_MatrixType>, Rhs>
419 typedef UmfPackLU<_MatrixType> Dec;
420 EIGEN_MAKE_SPARSE_SOLVE_HELPERS(Dec,Rhs)
422 template<typename Dest>
void evalTo(Dest& dst)
const
424 this->defaultEvalTo(dst);
432 #endif // EIGEN_UMFPACKSUPPORT_H