10 #ifndef EIGEN_UMFPACKSUPPORT_H 11 #define EIGEN_UMFPACKSUPPORT_H 20 { umfpack_di_free_numeric(Numeric); *Numeric = 0; }
23 { umfpack_zi_free_numeric(Numeric); *Numeric = 0; }
26 { umfpack_di_free_symbolic(Symbolic); *Symbolic = 0; }
29 { umfpack_zi_free_symbolic(Symbolic); *Symbolic = 0; }
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);
39 const int Ap[],
const int Ai[],
const std::complex<double> Ax[],
void **Symbolic,
40 const double Control [UMFPACK_CONTROL],
double Info [UMFPACK_INFO])
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])
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);
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)
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);
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])
107 return umfpack_zi_get_determinant(&mx_real,0,Ex,NumericHandle,User_Info);
123 template<
typename _MatrixType>
128 typedef typename MatrixType::Scalar
Scalar;
130 typedef typename MatrixType::Index
Index;
205 template<
typename Rhs>
210 &&
"UmfPackLU::solve(): invalid number of rows of the right hand side matrix b");
218 template<
typename Rhs>
223 &&
"UmfPackLU::solve(): invalid number of rows of the right hand side matrix b");
274 #ifndef EIGEN_PARSED_BY_DOXYGEN 276 template<
typename BDerived,
typename XDerived>
301 if( ((MatrixType::Flags&
RowMajorBit)==RowMajorBit) ||
sizeof(
typename MatrixType::Index)!=
sizeof(
int) || !mat.isCompressed() )
320 mutable IntColVectorType
m_p;
321 mutable IntRowVectorType
m_q;
341 template<
typename MatrixType>
369 template<
typename MatrixType>
377 template<
typename MatrixType>
378 template<
typename BDerived,
typename XDerived>
381 const int rhsCols = b.cols();
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)
391 &x.
col(j).coeffRef(0), &b.const_cast_derived().
col(j).coeffRef(0),
m_numeric, 0, 0);
402 template<
typename _MatrixType,
typename Rhs>
409 template<typename Dest>
void evalTo(Dest& dst)
const 411 dec()._solve(rhs(),dst);
415 template<
typename _MatrixType,
typename Rhs>
422 template<typename Dest>
void evalTo(Dest& dst)
const 424 this->defaultEvalTo(dst);
432 #endif // EIGEN_UMFPACKSUPPORT_H Matrix< int, MatrixType::RowsAtCompileTime, 1 > IntColVectorType
int umfpack_get_lunz(int *lnz, int *unz, int *n_row, int *n_col, int *nz_udiag, void *Numeric, double)
A sparse LU factorization and solver based on UmfPack.
UmfPackLU< _MatrixType > Dec
void compute(const MatrixType &matrix)
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])
void resizeNonZeros(Index size)
const internal::solve_retval< UmfPackLU, Rhs > solve(const MatrixBase< Rhs > &b) const
Matrix< Scalar, Dynamic, 1 > Vector
const IntColVectorType & permutationP() const
UmfPackLU< _MatrixType > Dec
int umfpack_get_determinant(double *Mx, double *Ex, void *NumericHandle, double User_Info[UMFPACK_INFO])
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])
const unsigned int RowMajorBit
MatrixType::RealScalar RealScalar
const Index * outerIndexPtr() const
const IntRowVectorType & permutationQ() const
SparseMatrix< Scalar > LUMatrixType
void umfpack_free_numeric(void **Numeric, double)
Base class of any sparse matrices or sparse expressions.
void grapInput(const MatrixType &mat)
EIGEN_STRONG_INLINE void resize(Index nbRows, Index nbCols)
void resize(Index rows, Index cols)
const LUMatrixType & matrixL() const
TFSIMD_FORCE_INLINE const tfScalar & x() const
const Scalar * m_valuePtr
const int * m_outerIndexPtr
#define EIGEN_MAKE_SPARSE_SOLVE_HELPERS(DecompositionType, Rhs)
ComputationInfo info() const
Reports whether previous computation was successful.
const LUMatrixType & matrixU() const
EIGEN_STRONG_INLINE const Scalar * data() const
Scalar determinant() const
UmfpackMatrixType m_copyMatrix
bool m_extractedDataAreDirty
Matrix< int, 1, MatrixType::ColsAtCompileTime > IntRowVectorType
const Derived & derived() const
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)
MatrixType::Scalar Scalar
void factorize(const MatrixType &matrix)
const int * m_innerIndexPtr
const Index * innerIndexPtr() const
#define EIGEN_MAKE_SOLVE_HELPERS(DecompositionType, Rhs)
internal::add_const_on_value_type< EIGEN_MATHFUNC_RETVAL(real_ref, Scalar) >::type real_ref(const Scalar &x)
const internal::sparse_solve_retval< UmfPackLU, Rhs > solve(const SparseMatrixBase< Rhs > &b) const
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])
const Scalar * valuePtr() const
bool _solve(const MatrixBase< BDerived > &b, MatrixBase< XDerived > &x) const
Base class for all dense matrices, vectors, and expressions.
void umfpack_free_symbolic(void **Symbolic, double)
SparseMatrix< Scalar, ColMajor, int > UmfpackMatrixType
UmfPackLU(const MatrixType &matrix)
void analyzePattern(const MatrixType &matrix)