Allows Lower|Upper as a template argument of CG and MINRES: in this case the full matrix will be considered.

This commit is contained in:
Gael Guennebaud 2015-02-10 18:57:41 +01:00
parent f9931a0392
commit 7b35b4cacc
5 changed files with 43 additions and 13 deletions

View File

@ -112,9 +112,9 @@ struct traits<ConjugateGradient<_MatrixType,_UpLo,_Preconditioner> >
* This class allows to solve for A.x = b sparse linear problems using a conjugate gradient algorithm.
* The sparse matrix A must be selfadjoint. The vectors x and b can be either dense or sparse.
*
* \tparam _MatrixType the type of the sparse matrix A, can be a dense or a sparse matrix.
* \tparam _UpLo the triangular part that will be used for the computations. It can be Lower
* or Upper. Default is Lower.
* \tparam _MatrixType the type of the matrix A, can be a dense or a sparse matrix.
* \tparam _UpLo the triangular part that will be used for the computations. It can be Lower,
* Upper, or Lower|Upper in which the full matrix entries will be considered. Default is Lower.
* \tparam _Preconditioner the type of the preconditioner. Default is DiagonalPreconditioner
*
* The maximal number of iterations and tolerance value can be controlled via the setMaxIterations()
@ -213,6 +213,10 @@ public:
template<typename Rhs,typename Dest>
void _solveWithGuess(const Rhs& b, Dest& x) const
{
typedef typename internal::conditional<UpLo==(Lower|Upper),
const MatrixType&,
SparseSelfAdjointView<const MatrixType, UpLo>
>::type MatrixWrapperType;
m_iterations = Base::maxIterations();
m_error = Base::m_tolerance;
@ -222,8 +226,7 @@ public:
m_error = Base::m_tolerance;
typename Dest::ColXpr xj(x,j);
internal::conjugate_gradient(mp_matrix->template selfadjointView<UpLo>(), b.col(j), xj,
Base::m_preconditioner, m_iterations, m_error);
internal::conjugate_gradient(MatrixWrapperType(*mp_matrix), b.col(j), xj, Base::m_preconditioner, m_iterations, m_error);
}
m_isInitialized = true;

View File

@ -12,13 +12,15 @@
template<typename T> void test_conjugate_gradient_T()
{
ConjugateGradient<SparseMatrix<T>, Lower> cg_colmajor_lower_diag;
ConjugateGradient<SparseMatrix<T>, Upper> cg_colmajor_upper_diag;
ConjugateGradient<SparseMatrix<T>, Lower > cg_colmajor_lower_diag;
ConjugateGradient<SparseMatrix<T>, Upper > cg_colmajor_upper_diag;
ConjugateGradient<SparseMatrix<T>, Lower|Upper> cg_colmajor_loup_diag;
ConjugateGradient<SparseMatrix<T>, Lower, IdentityPreconditioner> cg_colmajor_lower_I;
ConjugateGradient<SparseMatrix<T>, Upper, IdentityPreconditioner> cg_colmajor_upper_I;
CALL_SUBTEST( check_sparse_spd_solving(cg_colmajor_lower_diag) );
CALL_SUBTEST( check_sparse_spd_solving(cg_colmajor_upper_diag) );
CALL_SUBTEST( check_sparse_spd_solving(cg_colmajor_loup_diag) );
CALL_SUBTEST( check_sparse_spd_solving(cg_colmajor_lower_I) );
CALL_SUBTEST( check_sparse_spd_solving(cg_colmajor_upper_I) );
}

View File

@ -161,7 +161,10 @@ int generate_sparse_spd_problem(Solver& , typename Solver::MatrixType& A, typena
dA = dM * dM.adjoint();
halfA.resize(size,size);
halfA.template selfadjointView<Solver::UpLo>().rankUpdate(M);
if(Solver::UpLo==(Lower|Upper))
halfA = A;
else
halfA.template selfadjointView<Solver::UpLo>().rankUpdate(M);
return size;
}

View File

@ -250,6 +250,11 @@ namespace Eigen {
template<typename Rhs,typename Dest>
void _solveWithGuess(const Rhs& b, Dest& x) const
{
typedef typename internal::conditional<UpLo==(Lower|Upper),
const MatrixType&,
SparseSelfAdjointView<const MatrixType, UpLo>
>::type MatrixWrapperType;
m_iterations = Base::maxIterations();
m_error = Base::m_tolerance;
@ -259,7 +264,7 @@ namespace Eigen {
m_error = Base::m_tolerance;
typename Dest::ColXpr xj(x,j);
internal::minres(mp_matrix->template selfadjointView<UpLo>(), b.col(j), xj,
internal::minres(MatrixWrapperType(*mp_matrix), b.col(j), xj,
Base::m_preconditioner, m_iterations, m_error);
}

View File

@ -14,15 +14,32 @@
template<typename T> void test_minres_T()
{
MINRES<SparseMatrix<T>, Lower, DiagonalPreconditioner<T> > minres_colmajor_diag;
MINRES<SparseMatrix<T>, Lower, IdentityPreconditioner > minres_colmajor_I;
MINRES<SparseMatrix<T>, Lower|Upper, DiagonalPreconditioner<T> > minres_colmajor_diag;
MINRES<SparseMatrix<T>, Lower, IdentityPreconditioner > minres_colmajor_lower_I;
MINRES<SparseMatrix<T>, Upper, IdentityPreconditioner > minres_colmajor_upper_I;
// MINRES<SparseMatrix<T>, Lower, IncompleteLUT<T> > minres_colmajor_ilut;
//minres<SparseMatrix<T>, SSORPreconditioner<T> > minres_colmajor_ssor;
CALL_SUBTEST( check_sparse_square_solving(minres_colmajor_diag) );
CALL_SUBTEST( check_sparse_spd_solving(minres_colmajor_I) );
// CALL_SUBTEST( check_sparse_square_solving(minres_colmajor_diag) );
// CALL_SUBTEST( check_sparse_square_solving(minres_colmajor_ilut) );
//CALL_SUBTEST( check_sparse_square_solving(minres_colmajor_ssor) );
// Diagonal preconditioner
MINRES<SparseMatrix<T>, Lower, DiagonalPreconditioner<T> > minres_colmajor_lower_diag;
MINRES<SparseMatrix<T>, Upper, DiagonalPreconditioner<T> > minres_colmajor_upper_diag;
MINRES<SparseMatrix<T>, Upper|Lower, DiagonalPreconditioner<T> > minres_colmajor_uplo_diag;
// call tests for SPD matrix
CALL_SUBTEST( check_sparse_spd_solving(minres_colmajor_lower_I) );
CALL_SUBTEST( check_sparse_spd_solving(minres_colmajor_upper_I) );
CALL_SUBTEST( check_sparse_spd_solving(minres_colmajor_lower_diag) );
CALL_SUBTEST( check_sparse_spd_solving(minres_colmajor_upper_diag) );
// CALL_SUBTEST( check_sparse_spd_solving(minres_colmajor_uplo_diag) );
// TO DO: symmetric semi-definite matrix
// TO DO: symmetric indefinite matrix
}
void test_minres()