mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-04-30 07:44:10 +08:00
Merged in joaoruileal/eigen (pull request PR-276)
Minor improvements to Umfpack support
This commit is contained in:
commit
c6882a72ed
@ -23,6 +23,24 @@ inline void umfpack_defaults(double control[UMFPACK_CONTROL], double)
|
|||||||
inline void umfpack_defaults(double control[UMFPACK_CONTROL], std::complex<double>)
|
inline void umfpack_defaults(double control[UMFPACK_CONTROL], std::complex<double>)
|
||||||
{ umfpack_zi_defaults(control); }
|
{ umfpack_zi_defaults(control); }
|
||||||
|
|
||||||
|
inline void umfpack_report_info(double control[UMFPACK_CONTROL], double info[UMFPACK_INFO], double)
|
||||||
|
{ umfpack_di_report_info(control, info);}
|
||||||
|
|
||||||
|
inline void umfpack_report_info(double control[UMFPACK_CONTROL], double info[UMFPACK_INFO], std::complex<double>)
|
||||||
|
{ umfpack_zi_report_info(control, info);}
|
||||||
|
|
||||||
|
inline void umfpack_report_status(double control[UMFPACK_CONTROL], int status, double)
|
||||||
|
{ umfpack_di_report_status(control, status);}
|
||||||
|
|
||||||
|
inline void umfpack_report_status(double control[UMFPACK_CONTROL], int status, std::complex<double>)
|
||||||
|
{ umfpack_zi_report_status(control, status);}
|
||||||
|
|
||||||
|
inline void umfpack_report_control(double control[UMFPACK_CONTROL], double)
|
||||||
|
{ umfpack_di_report_control(control);}
|
||||||
|
|
||||||
|
inline void umfpack_report_control(double control[UMFPACK_CONTROL], std::complex<double>)
|
||||||
|
{ umfpack_zi_report_control(control);}
|
||||||
|
|
||||||
inline void umfpack_free_numeric(void **Numeric, double)
|
inline void umfpack_free_numeric(void **Numeric, double)
|
||||||
{ umfpack_di_free_numeric(Numeric); *Numeric = 0; }
|
{ umfpack_di_free_numeric(Numeric); *Numeric = 0; }
|
||||||
|
|
||||||
@ -156,6 +174,7 @@ class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
|
|||||||
public:
|
public:
|
||||||
|
|
||||||
typedef Array<double, UMFPACK_CONTROL, 1> UmfpackControl;
|
typedef Array<double, UMFPACK_CONTROL, 1> UmfpackControl;
|
||||||
|
typedef Array<double, UMFPACK_INFO, 1> UmfpackInfo;
|
||||||
|
|
||||||
UmfPackLU()
|
UmfPackLU()
|
||||||
: m_dummy(0,0), mp_matrix(m_dummy)
|
: m_dummy(0,0), mp_matrix(m_dummy)
|
||||||
@ -297,6 +316,34 @@ class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
|
|||||||
factorize_impl();
|
factorize_impl();
|
||||||
}
|
}
|
||||||
|
|
||||||
|
/** Prints the current UmfPack control settings.
|
||||||
|
*
|
||||||
|
* \sa umfpackControl()
|
||||||
|
*/
|
||||||
|
void printUmfpackControl()
|
||||||
|
{
|
||||||
|
umfpack_report_control(m_control.data(), Scalar());
|
||||||
|
}
|
||||||
|
|
||||||
|
/** Prints statistics collected by UmfPack.
|
||||||
|
*
|
||||||
|
* \sa analyzePattern(), compute()
|
||||||
|
*/
|
||||||
|
void printUmfpackInfo()
|
||||||
|
{
|
||||||
|
eigen_assert(m_analysisIsOk && "UmfPackLU: you must first call analyzePattern()");
|
||||||
|
umfpack_report_info(m_control.data(), m_umfpackInfo.data(), Scalar());
|
||||||
|
}
|
||||||
|
|
||||||
|
/** Prints the status of the previous factorization operation performed by UmfPack (symbolic or numerical factorization).
|
||||||
|
*
|
||||||
|
* \sa analyzePattern(), compute()
|
||||||
|
*/
|
||||||
|
void printUmfpackStatus() {
|
||||||
|
eigen_assert(m_analysisIsOk && "UmfPackLU: you must first call analyzePattern()");
|
||||||
|
umfpack_report_status(m_control.data(), m_fact_errorCode, Scalar());
|
||||||
|
}
|
||||||
|
|
||||||
/** \internal */
|
/** \internal */
|
||||||
template<typename BDerived,typename XDerived>
|
template<typename BDerived,typename XDerived>
|
||||||
bool _solve_impl(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const;
|
bool _solve_impl(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const;
|
||||||
@ -314,19 +361,19 @@ class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
|
|||||||
m_numeric = 0;
|
m_numeric = 0;
|
||||||
m_symbolic = 0;
|
m_symbolic = 0;
|
||||||
m_extractedDataAreDirty = true;
|
m_extractedDataAreDirty = true;
|
||||||
|
|
||||||
|
umfpack_defaults(m_control.data(), Scalar());
|
||||||
}
|
}
|
||||||
|
|
||||||
void analyzePattern_impl()
|
void analyzePattern_impl()
|
||||||
{
|
{
|
||||||
umfpack_defaults(m_control.data(), Scalar());
|
m_fact_errorCode = umfpack_symbolic(internal::convert_index<int>(mp_matrix.rows()),
|
||||||
int errorCode = 0;
|
|
||||||
errorCode = umfpack_symbolic(internal::convert_index<int>(mp_matrix.rows()),
|
|
||||||
internal::convert_index<int>(mp_matrix.cols()),
|
internal::convert_index<int>(mp_matrix.cols()),
|
||||||
mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
|
mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
|
||||||
&m_symbolic, m_control.data(), 0);
|
&m_symbolic, m_control.data(), m_umfpackInfo.data());
|
||||||
|
|
||||||
m_isInitialized = true;
|
m_isInitialized = true;
|
||||||
m_info = errorCode ? InvalidInput : Success;
|
m_info = m_fact_errorCode ? InvalidInput : Success;
|
||||||
m_analysisIsOk = true;
|
m_analysisIsOk = true;
|
||||||
m_factorizationIsOk = false;
|
m_factorizationIsOk = false;
|
||||||
m_extractedDataAreDirty = true;
|
m_extractedDataAreDirty = true;
|
||||||
@ -334,8 +381,9 @@ class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
|
|||||||
|
|
||||||
void factorize_impl()
|
void factorize_impl()
|
||||||
{
|
{
|
||||||
|
|
||||||
m_fact_errorCode = umfpack_numeric(mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
|
m_fact_errorCode = umfpack_numeric(mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
|
||||||
m_symbolic, &m_numeric, m_control.data(), 0);
|
m_symbolic, &m_numeric, m_control.data(), m_umfpackInfo.data());
|
||||||
|
|
||||||
m_info = m_fact_errorCode == UMFPACK_OK ? Success : NumericalIssue;
|
m_info = m_fact_errorCode == UMFPACK_OK ? Success : NumericalIssue;
|
||||||
m_factorizationIsOk = true;
|
m_factorizationIsOk = true;
|
||||||
@ -362,6 +410,7 @@ class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
|
|||||||
mutable LUMatrixType m_l;
|
mutable LUMatrixType m_l;
|
||||||
int m_fact_errorCode;
|
int m_fact_errorCode;
|
||||||
UmfpackControl m_control;
|
UmfpackControl m_control;
|
||||||
|
UmfpackInfo m_umfpackInfo;
|
||||||
|
|
||||||
mutable LUMatrixType m_u;
|
mutable LUMatrixType m_u;
|
||||||
mutable IntColVectorType m_p;
|
mutable IntColVectorType m_p;
|
||||||
@ -442,7 +491,7 @@ bool UmfPackLU<MatrixType>::_solve_impl(const MatrixBase<BDerived> &b, MatrixBas
|
|||||||
x_ptr = &x.col(j).coeffRef(0);
|
x_ptr = &x.col(j).coeffRef(0);
|
||||||
errorCode = umfpack_solve(UMFPACK_A,
|
errorCode = umfpack_solve(UMFPACK_A,
|
||||||
mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
|
mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
|
||||||
x_ptr, &b.const_cast_derived().col(j).coeffRef(0), m_numeric, m_control.data(), 0);
|
x_ptr, &b.const_cast_derived().col(j).coeffRef(0), m_numeric, m_control.data(), m_umfpackInfo.data());
|
||||||
if(x.innerStride()!=1)
|
if(x.innerStride()!=1)
|
||||||
x.col(j) = x_tmp;
|
x.col(j) = x_tmp;
|
||||||
if (errorCode!=0)
|
if (errorCode!=0)
|
||||||
|
Loading…
x
Reference in New Issue
Block a user