]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add test case
authorMartin Kronbichler <martin.kronbichler@uni-a.de>
Fri, 6 May 2022 12:02:27 +0000 (14:02 +0200)
committerMartin Kronbichler <martin.kronbichler@uni-a.de>
Mon, 9 May 2022 09:03:19 +0000 (11:03 +0200)
tests/lac/solver_cg_flexible.cc [new file with mode: 0644]
tests/lac/solver_cg_flexible.output [new file with mode: 0644]

diff --git a/tests/lac/solver_cg_flexible.cc b/tests/lac/solver_cg_flexible.cc
new file mode 100644 (file)
index 0000000..681a001
--- /dev/null
@@ -0,0 +1,101 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 1998 - 2020 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.
+//
+// ---------------------------------------------------------------------
+
+// Check flexible variant of SolverCG by verifying that it converges more
+// quickly in case of a variable preconditioner
+
+
+#include <deal.II/lac/diagonal_matrix.h>
+#include <deal.II/lac/precondition.h>
+#include <deal.II/lac/solver_cg.h>
+#include <deal.II/lac/vector.h>
+
+#include "../tests.h"
+
+
+// Class to represent an approximate inverse with variable number of
+// iterations to stress the flexible CG implementation
+class IterativeInverse
+{
+public:
+  IterativeInverse(const DiagonalMatrix<Vector<double>> &matrix)
+    : matrix(matrix)
+    , apply_count(0)
+  {}
+
+  void
+  vmult(Vector<double> &dst, const Vector<double> &src) const
+  {
+    IterationNumberControl control(apply_count % 2 ? 4 : 3, 1e-12);
+    SolverCG<>             inner_solver(control);
+
+    dst = 0;
+    inner_solver.solve(matrix, dst, src, PreconditionIdentity());
+
+    ++apply_count;
+  }
+
+private:
+  const DiagonalMatrix<Vector<double>> &matrix;
+  mutable unsigned int                  apply_count;
+};
+
+
+
+SolverControl::State
+monitor_norm(const unsigned int iteration,
+             const double       check_value,
+             const Vector<double> &)
+{
+  deallog << "   CG estimated residual at iteration " << iteration << ": "
+          << check_value << std::endl;
+  return SolverControl::success;
+}
+
+
+int
+main()
+{
+  initlog();
+
+  // Create diagonal matrix with entries between 1 and 30
+  DiagonalMatrix<Vector<double>> matrix;
+  matrix.get_vector().reinit(30);
+  for (unsigned int i = 0; i < matrix.m(); ++i)
+    matrix.get_vector()[i] = i + 1.0;
+
+  Vector<double> rhs(matrix.m()), sol(matrix.m());
+  rhs = 1.;
+
+  IterativeInverse variable_preconditioner(matrix);
+
+  // Try to solve with standard CG method - note that we solve with a coarse
+  // tolerance to avoid trouble with roundoff
+  {
+    SolverControl control(30, 1e-4);
+    SolverCG<>    solver(control);
+    solver.connect(&monitor_norm);
+    solver.solve(matrix, sol, rhs, variable_preconditioner);
+  }
+
+  // And now use the flexible CG method
+  {
+    SolverControl      control(30, 1e-4);
+    SolverFlexibleCG<> solver(control);
+    solver.connect(&monitor_norm);
+    sol = 0;
+    solver.solve(matrix, sol, rhs, variable_preconditioner);
+  }
+}
diff --git a/tests/lac/solver_cg_flexible.output b/tests/lac/solver_cg_flexible.output
new file mode 100644 (file)
index 0000000..0daa023
--- /dev/null
@@ -0,0 +1,91 @@
+
+DEAL:cg::Starting value 5.47723
+DEAL:cg::   CG estimated residual at iteration 0: 5.47723
+DEAL:cg:cg::Starting value 5.47723
+DEAL:cg:cg::Convergence step 3 value 1.69418
+DEAL:cg::   CG estimated residual at iteration 1: 1.69418
+DEAL:cg:cg::Starting value 1.69418
+DEAL:cg:cg::Convergence step 4 value 0.640168
+DEAL:cg::   CG estimated residual at iteration 2: 1.16797
+DEAL:cg:cg::Starting value 1.16797
+DEAL:cg:cg::Convergence step 3 value 0.340792
+DEAL:cg::   CG estimated residual at iteration 3: 0.581169
+DEAL:cg:cg::Starting value 0.581169
+DEAL:cg:cg::Convergence step 4 value 0.0744913
+DEAL:cg::   CG estimated residual at iteration 4: 0.326845
+DEAL:cg:cg::Starting value 0.326845
+DEAL:cg:cg::Convergence step 3 value 0.0750283
+DEAL:cg::   CG estimated residual at iteration 5: 0.110157
+DEAL:cg:cg::Starting value 0.110157
+DEAL:cg:cg::Convergence step 4 value 0.0248289
+DEAL:cg::   CG estimated residual at iteration 6: 0.0580666
+DEAL:cg:cg::Starting value 0.0580666
+DEAL:cg:cg::Convergence step 3 value 0.0299569
+DEAL:cg::   CG estimated residual at iteration 7: 0.0311519
+DEAL:cg:cg::Starting value 0.0311519
+DEAL:cg:cg::Convergence step 4 value 0.00678179
+DEAL:cg::   CG estimated residual at iteration 8: 0.0196211
+DEAL:cg:cg::Starting value 0.0196211
+DEAL:cg:cg::Convergence step 3 value 0.00345458
+DEAL:cg::   CG estimated residual at iteration 9: 0.00776886
+DEAL:cg:cg::Starting value 0.00776886
+DEAL:cg:cg::Convergence step 4 value 0.00116611
+DEAL:cg::   CG estimated residual at iteration 10: 0.00285322
+DEAL:cg:cg::Starting value 0.00285322
+DEAL:cg:cg::Convergence step 3 value 0.00157219
+DEAL:cg::   CG estimated residual at iteration 11: 0.00166226
+DEAL:cg:cg::Starting value 0.00166226
+DEAL:cg:cg::Convergence step 4 value 0.000451853
+DEAL:cg::   CG estimated residual at iteration 12: 0.00112345
+DEAL:cg:cg::Starting value 0.00112345
+DEAL:cg:cg::Convergence step 3 value 0.000250859
+DEAL:cg::   CG estimated residual at iteration 13: 0.000392171
+DEAL:cg:cg::Starting value 0.000392171
+DEAL:cg:cg::Convergence step 4 value 5.59076e-05
+DEAL:cg::   CG estimated residual at iteration 14: 0.000235130
+DEAL:cg:cg::Starting value 0.000235130
+DEAL:cg:cg::Convergence step 3 value 6.46587e-05
+DEAL:cg::   CG estimated residual at iteration 15: 0.000105880
+DEAL:cg:cg::Starting value 0.000105880
+DEAL:cg:cg::Convergence step 4 value 2.74842e-05
+DEAL:cg::Convergence step 16 value 4.78582e-05
+DEAL:cg::   CG estimated residual at iteration 16: 4.78582e-05
+DEAL:cg::Starting value 5.47723
+DEAL:cg::   CG estimated residual at iteration 0: 5.47723
+DEAL:cg:cg::Starting value 5.47723
+DEAL:cg:cg::Convergence step 3 value 1.69418
+DEAL:cg::   CG estimated residual at iteration 1: 1.69418
+DEAL:cg:cg::Starting value 1.69418
+DEAL:cg:cg::Convergence step 4 value 0.640168
+DEAL:cg::   CG estimated residual at iteration 2: 1.16797
+DEAL:cg:cg::Starting value 1.16797
+DEAL:cg:cg::Convergence step 3 value 0.340792
+DEAL:cg::   CG estimated residual at iteration 3: 0.399820
+DEAL:cg:cg::Starting value 0.399820
+DEAL:cg:cg::Convergence step 4 value 0.0868174
+DEAL:cg::   CG estimated residual at iteration 4: 0.200425
+DEAL:cg:cg::Starting value 0.200425
+DEAL:cg:cg::Convergence step 3 value 0.0482953
+DEAL:cg::   CG estimated residual at iteration 5: 0.0551085
+DEAL:cg:cg::Starting value 0.0551085
+DEAL:cg:cg::Convergence step 4 value 0.0149675
+DEAL:cg::   CG estimated residual at iteration 6: 0.0295425
+DEAL:cg:cg::Starting value 0.0295425
+DEAL:cg:cg::Convergence step 3 value 0.00812170
+DEAL:cg::   CG estimated residual at iteration 7: 0.00923062
+DEAL:cg:cg::Starting value 0.00923062
+DEAL:cg:cg::Convergence step 4 value 0.00218874
+DEAL:cg::   CG estimated residual at iteration 8: 0.00475654
+DEAL:cg:cg::Starting value 0.00475654
+DEAL:cg:cg::Convergence step 3 value 0.00103994
+DEAL:cg::   CG estimated residual at iteration 9: 0.00102251
+DEAL:cg:cg::Starting value 0.00102251
+DEAL:cg:cg::Convergence step 4 value 0.000231461
+DEAL:cg::   CG estimated residual at iteration 10: 0.000425978
+DEAL:cg:cg::Starting value 0.000425978
+DEAL:cg:cg::Convergence step 3 value 0.000140415
+DEAL:cg::   CG estimated residual at iteration 11: 0.000161655
+DEAL:cg:cg::Starting value 0.000161655
+DEAL:cg:cg::Convergence step 4 value 4.27639e-05
+DEAL:cg::Convergence step 12 value 7.45407e-05
+DEAL:cg::   CG estimated residual at iteration 12: 7.45407e-05

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.