aboutsummaryrefslogtreecommitdiffhomepage
path: root/unsupported
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2015-02-10 18:57:41 +0100
committerGravatar Gael Guennebaud <g.gael@free.fr>2015-02-10 18:57:41 +0100
commitc6e8caf0900ae303e9e7399bed00af705015ff17 (patch)
tree384831cf695f94eb5fca5744a1b05a9cf8930e85 /unsupported
parentd10d6a40dda3fb5ac9f401b8e6d9cede3f3ca34a (diff)
Allows Lower|Upper as a template argument of CG and MINRES: in this case the full matrix will be considered.
Diffstat (limited to 'unsupported')
-rw-r--r--unsupported/Eigen/src/IterativeSolvers/MINRES.h7
-rw-r--r--unsupported/test/minres.cpp2
2 files changed, 8 insertions, 1 deletions
diff --git a/unsupported/Eigen/src/IterativeSolvers/MINRES.h b/unsupported/Eigen/src/IterativeSolvers/MINRES.h
index a34902001..65cffc255 100644
--- a/unsupported/Eigen/src/IterativeSolvers/MINRES.h
+++ b/unsupported/Eigen/src/IterativeSolvers/MINRES.h
@@ -250,6 +250,11 @@ namespace Eigen {
template<typename Rhs,typename Dest>
void _solve_with_guess_impl(const Rhs& b, Dest& x) const
{
+ typedef typename internal::conditional<UpLo==(Lower|Upper),
+ Ref<const MatrixType>&,
+ SparseSelfAdjointView<const Ref<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);
}
diff --git a/unsupported/test/minres.cpp b/unsupported/test/minres.cpp
index 81b762c37..f6e526bbd 100644
--- a/unsupported/test/minres.cpp
+++ b/unsupported/test/minres.cpp
@@ -21,6 +21,7 @@ template<typename T> void test_minres_T()
// 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, DiagonalPreconditioner<T> > minres_colmajor_uplo_diag;
// call tests for SPD matrix
CALL_SUBTEST( check_sparse_spd_solving(minres_colmajor_lower_I) );
@@ -28,6 +29,7 @@ template<typename T> void test_minres_T()
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