mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-04-29 23:34:12 +08:00
Merged in joaoruileal/eigen (pull request PR-276)
Minor improvements to Umfpack support
This commit is contained in:
commit
c6882a72ed
@ -10,19 +10,37 @@
|
|||||||
#ifndef EIGEN_UMFPACKSUPPORT_H
|
#ifndef EIGEN_UMFPACKSUPPORT_H
|
||||||
#define EIGEN_UMFPACKSUPPORT_H
|
#define EIGEN_UMFPACKSUPPORT_H
|
||||||
|
|
||||||
namespace Eigen {
|
namespace Eigen {
|
||||||
|
|
||||||
/* TODO extract L, extract U, compute det, etc... */
|
/* TODO extract L, extract U, compute det, etc... */
|
||||||
|
|
||||||
// generic double/complex<double> wrapper functions:
|
// generic double/complex<double> wrapper functions:
|
||||||
|
|
||||||
|
|
||||||
inline void umfpack_defaults(double control[UMFPACK_CONTROL], double)
|
inline void umfpack_defaults(double control[UMFPACK_CONTROL], double)
|
||||||
{ umfpack_di_defaults(control); }
|
{ umfpack_di_defaults(control); }
|
||||||
|
|
||||||
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)
|
||||||
@ -215,7 +234,7 @@ class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
|
|||||||
return m_q;
|
return m_q;
|
||||||
}
|
}
|
||||||
|
|
||||||
/** Computes the sparse Cholesky decomposition of \a matrix
|
/** Computes the sparse Cholesky decomposition of \a matrix
|
||||||
* Note that the matrix should be column-major, and in compressed format for best performance.
|
* Note that the matrix should be column-major, and in compressed format for best performance.
|
||||||
* \sa SparseMatrix::makeCompressed().
|
* \sa SparseMatrix::makeCompressed().
|
||||||
*/
|
*/
|
||||||
@ -240,7 +259,7 @@ class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
|
|||||||
{
|
{
|
||||||
if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
|
if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
|
||||||
if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
|
if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
|
||||||
|
|
||||||
grab(matrix.derived());
|
grab(matrix.derived());
|
||||||
|
|
||||||
analyzePattern_impl();
|
analyzePattern_impl();
|
||||||
@ -267,7 +286,7 @@ class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
|
|||||||
{
|
{
|
||||||
return m_control;
|
return m_control;
|
||||||
}
|
}
|
||||||
|
|
||||||
/** Provides access to the control settings array used by UmfPack.
|
/** Provides access to the control settings array used by UmfPack.
|
||||||
*
|
*
|
||||||
* If this array contains NaN's, the default values are used.
|
* If this array contains NaN's, the default values are used.
|
||||||
@ -278,7 +297,7 @@ class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
|
|||||||
{
|
{
|
||||||
return m_control;
|
return m_control;
|
||||||
}
|
}
|
||||||
|
|
||||||
/** Performs a numeric decomposition of \a matrix
|
/** Performs a numeric decomposition of \a matrix
|
||||||
*
|
*
|
||||||
* The given matrix must has the same sparcity than the matrix on which the pattern anylysis has been performed.
|
* The given matrix must has the same sparcity than the matrix on which the pattern anylysis has been performed.
|
||||||
@ -293,10 +312,38 @@ class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
|
|||||||
umfpack_free_numeric(&m_numeric,Scalar());
|
umfpack_free_numeric(&m_numeric,Scalar());
|
||||||
|
|
||||||
grab(matrix.derived());
|
grab(matrix.derived());
|
||||||
|
|
||||||
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,41 +361,42 @@ 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;
|
internal::convert_index<int>(mp_matrix.cols()),
|
||||||
errorCode = umfpack_symbolic(internal::convert_index<int>(mp_matrix.rows()),
|
mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
|
||||||
internal::convert_index<int>(mp_matrix.cols()),
|
&m_symbolic, m_control.data(), m_umfpackInfo.data());
|
||||||
mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
|
|
||||||
&m_symbolic, m_control.data(), 0);
|
|
||||||
|
|
||||||
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;
|
||||||
}
|
}
|
||||||
|
|
||||||
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;
|
||||||
m_extractedDataAreDirty = true;
|
m_extractedDataAreDirty = true;
|
||||||
}
|
}
|
||||||
|
|
||||||
template<typename MatrixDerived>
|
template<typename MatrixDerived>
|
||||||
void grab(const EigenBase<MatrixDerived> &A)
|
void grab(const EigenBase<MatrixDerived> &A)
|
||||||
{
|
{
|
||||||
mp_matrix.~UmfpackMatrixRef();
|
mp_matrix.~UmfpackMatrixRef();
|
||||||
::new (&mp_matrix) UmfpackMatrixRef(A.derived());
|
::new (&mp_matrix) UmfpackMatrixRef(A.derived());
|
||||||
}
|
}
|
||||||
|
|
||||||
void grab(const UmfpackMatrixRef &A)
|
void grab(const UmfpackMatrixRef &A)
|
||||||
{
|
{
|
||||||
if(&(A.derived()) != &mp_matrix)
|
if(&(A.derived()) != &mp_matrix)
|
||||||
@ -357,19 +405,20 @@ class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
|
|||||||
::new (&mp_matrix) UmfpackMatrixRef(A);
|
::new (&mp_matrix) UmfpackMatrixRef(A);
|
||||||
}
|
}
|
||||||
}
|
}
|
||||||
|
|
||||||
// cached data to reduce reallocation, etc.
|
// cached data to reduce reallocation, etc.
|
||||||
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;
|
||||||
mutable IntRowVectorType m_q;
|
mutable IntRowVectorType m_q;
|
||||||
|
|
||||||
UmfpackMatrixType m_dummy;
|
UmfpackMatrixType m_dummy;
|
||||||
UmfpackMatrixRef mp_matrix;
|
UmfpackMatrixRef mp_matrix;
|
||||||
|
|
||||||
void* m_numeric;
|
void* m_numeric;
|
||||||
void* m_symbolic;
|
void* m_symbolic;
|
||||||
|
|
||||||
@ -377,7 +426,7 @@ class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
|
|||||||
int m_factorizationIsOk;
|
int m_factorizationIsOk;
|
||||||
int m_analysisIsOk;
|
int m_analysisIsOk;
|
||||||
mutable bool m_extractedDataAreDirty;
|
mutable bool m_extractedDataAreDirty;
|
||||||
|
|
||||||
private:
|
private:
|
||||||
UmfPackLU(const UmfPackLU& ) { }
|
UmfPackLU(const UmfPackLU& ) { }
|
||||||
};
|
};
|
||||||
@ -427,7 +476,7 @@ bool UmfPackLU<MatrixType>::_solve_impl(const MatrixBase<BDerived> &b, MatrixBas
|
|||||||
eigen_assert((BDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major rhs yet");
|
eigen_assert((BDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major rhs yet");
|
||||||
eigen_assert((XDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major result yet");
|
eigen_assert((XDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major result yet");
|
||||||
eigen_assert(b.derived().data() != x.derived().data() && " Umfpack does not support inplace solve");
|
eigen_assert(b.derived().data() != x.derived().data() && " Umfpack does not support inplace solve");
|
||||||
|
|
||||||
int errorCode;
|
int errorCode;
|
||||||
Scalar* x_ptr = 0;
|
Scalar* x_ptr = 0;
|
||||||
Matrix<Scalar,Dynamic,1> x_tmp;
|
Matrix<Scalar,Dynamic,1> x_tmp;
|
||||||
@ -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