add analyzePattern/factorize API to iterative solvers and basic preconditioners

This commit is contained in:
Gael Guennebaud 2012-02-27 14:10:26 +01:00
parent 122f28626c
commit eb168ef8ed
3 changed files with 62 additions and 7 deletions

View File

@ -62,9 +62,15 @@ class DiagonalPreconditioner
Index rows() const { return m_invdiag.size(); }
Index cols() const { return m_invdiag.size(); }
template<typename MatrixType>
DiagonalPreconditioner& compute(const MatrixType& mat)
DiagonalPreconditioner& analyzePattern(const MatrixType& )
{
return *this;
}
template<typename MatrixType>
DiagonalPreconditioner& factorize(const MatrixType& mat)
{
m_invdiag.resize(mat.cols());
for(int j=0; j<mat.outerSize(); ++j)
@ -79,6 +85,12 @@ class DiagonalPreconditioner
m_isInitialized = true;
return *this;
}
template<typename MatrixType>
DiagonalPreconditioner& compute(const MatrixType& mat)
{
return factorize(mat);
}
template<typename Rhs, typename Dest>
void _solve(const Rhs& b, Dest& x) const
@ -130,6 +142,12 @@ class IdentityPreconditioner
template<typename MatrixType>
IdentityPreconditioner(const MatrixType& ) {}
template<typename MatrixType>
IdentityPreconditioner& analyzePattern(const MatrixType& ) { return *this; }
template<typename MatrixType>
IdentityPreconditioner& factorize(const MatrixType& ) { return *this; }
template<typename MatrixType>
IdentityPreconditioner& compute(const MatrixType& ) { return *this; }

View File

@ -97,10 +97,10 @@ class IncompleteLUT
void factorize(const MatrixType& amat);
/**
* Compute an incomplete LU factorization with dual threshold on the matrix mat
* No pivoting is done in this version
*
**/
* Compute an incomplete LU factorization with dual threshold on the matrix mat
* No pivoting is done in this version
*
**/
template<typename MatrixType>
IncompleteLUT<Scalar>& compute(const MatrixType& amat)
{

View File

@ -70,6 +70,39 @@ public:
}
~IterativeSolverBase() {}
/** Initializes the iterative solver for the sparcity pattern of the matrix \a A for further solving \c Ax=b problems.
*
* Currently, this function mostly call analyzePattern on the preconditioner. In the future
* we might, for instance, implement column reodering for faster matrix vector products.
*/
Derived& analyzePattern(const MatrixType& A)
{
m_preconditioner.analyzePattern(A);
m_isInitialized = true;
m_analysisIsOk = true;
m_info = Success;
return derived();
}
/** Initializes the iterative solver with the numerical values of the matrix \a A for further solving \c Ax=b problems.
*
* Currently, this function mostly call factorize on the preconditioner.
*
* \warning this class stores a reference to the matrix A as well as some
* precomputed values that depend on it. Therefore, if \a A is changed
* this class becomes invalid. Call compute() to update it with the new
* matrix A, or modify a copy of A.
*/
Derived& factorize(const MatrixType& A)
{
eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
mp_matrix = &A;
m_preconditioner.factorize(A);
m_factorizationIsOk = true;
m_info = Success;
return derived();
}
/** Initializes the iterative solver with the matrix \a A for further solving \c Ax=b problems.
*
@ -86,6 +119,8 @@ public:
mp_matrix = &A;
m_preconditioner.compute(A);
m_isInitialized = true;
m_analysisIsOk = true;
m_factorizationIsOk = true;
m_info = Success;
return derived();
}
@ -194,6 +229,8 @@ protected:
void init()
{
m_isInitialized = false;
m_analysisIsOk = false;
m_factorizationIsOk = false;
m_maxIterations = -1;
m_tolerance = NumTraits<Scalar>::epsilon();
}
@ -206,7 +243,7 @@ protected:
mutable RealScalar m_error;
mutable int m_iterations;
mutable ComputationInfo m_info;
mutable bool m_isInitialized;
mutable bool m_isInitialized, m_analysisIsOk, m_factorizationIsOk;
};
namespace internal {