From: Daniel Arndt Date: Fri, 15 Jun 2018 13:56:55 +0000 (+0200) Subject: Avoid crash for no iteration in GMRES eigenvalue approximation X-Git-Tag: v9.1.0-rc1~1030^2 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=7e4a64d167e0dc8f8fa9ef7b47c67d0c81f4bc76;p=dealii.git Avoid crash for no iteration in GMRES eigenvalue approximation --- diff --git a/doc/news/changes/minor/20180615DanielArndt b/doc/news/changes/minor/20180615DanielArndt new file mode 100644 index 0000000000..1185c39a79 --- /dev/null +++ b/doc/news/changes/minor/20180615DanielArndt @@ -0,0 +1,4 @@ +Fixed: The eigenvalue approximation in SolverGMRES failed if no iterations +were performed. Now, there is simply no output in this case. +
+(Daniel Arndt, 2018/06/15) diff --git a/include/deal.II/lac/solver_gmres.h b/include/deal.II/lac/solver_gmres.h index 8fb6ce7786..e45306b99b 100644 --- a/include/deal.II/lac/solver_gmres.h +++ b/include/deal.II/lac/solver_gmres.h @@ -743,8 +743,9 @@ SolverGMRES::compute_eigs_and_cond( const boost::signals2::signal &cond_signal) { // Avoid copying the Hessenberg matrix if it isn't needed. - if (!eigenvalues_signal.empty() || !hessenberg_signal.empty() || - !cond_signal.empty()) + if ((!eigenvalues_signal.empty() || !hessenberg_signal.empty() || + !cond_signal.empty()) && + dim > 0) { LAPACKFullMatrix mat(dim, dim); for (unsigned int i = 0; i < dim; ++i) diff --git a/tests/lac/gmres_eigenvalues_02.cc b/tests/lac/gmres_eigenvalues_02.cc new file mode 100644 index 0000000000..cd7eb3786f --- /dev/null +++ b/tests/lac/gmres_eigenvalues_02.cc @@ -0,0 +1,67 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2018 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + + +// test that GMRES eigenvalue approximation doesn't crash +// when no iterations are performed. + +#include +#include +#include +#include + +#include "../tests.h" + +template +void +test() +{ + const unsigned int n = 2; + ; + Vector rhs(n), sol(n); + rhs = 0.; + + LAPACKFullMatrix matrix(n, n); + + for (unsigned int i = 0; i < n; ++i) + matrix(i, i) = std::sqrt(i + 1); + + SolverControl control; + { + SolverGMRES> solver(control); + solver.connect_eigenvalues_slot( + [](const std::vector> &eigenvalues) { + deallog << "n_eigenvalues: " << eigenvalues.size(); + }); + solver.solve(matrix, sol, rhs, PreconditionIdentity()); + } + + { + SolverCG> solver(control); + solver.connect_eigenvalues_slot([](const std::vector &eigenvalues) { + deallog << "n_eigenvalues: " << eigenvalues.size(); + }); + solver.solve(matrix, sol, rhs, PreconditionIdentity()); + } + + deallog << "OK" << std::endl; +} + +int +main() +{ + initlog(); + test(); +} diff --git a/tests/lac/gmres_eigenvalues_02.with_lapack=true.output b/tests/lac/gmres_eigenvalues_02.with_lapack=true.output new file mode 100644 index 0000000000..19dbb7fac0 --- /dev/null +++ b/tests/lac/gmres_eigenvalues_02.with_lapack=true.output @@ -0,0 +1,6 @@ + +DEAL:GMRES::Starting value 0.00000 +DEAL:GMRES::Convergence step 0 value 0.00000 +DEAL:cg::Starting value 0.00000 +DEAL:cg::Convergence step 0 value 0.00000 +DEAL::OK