]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add tests checking result of SolverControl classes with deal.II solver
authorJean-Paul Pelteret <jppelteret@gmail.com>
Sun, 29 Jan 2017 10:20:09 +0000 (11:20 +0100)
committerJean-Paul Pelteret <jppelteret@gmail.com>
Mon, 30 Jan 2017 13:28:30 +0000 (14:28 +0100)
tests/lac/solver_control_01.cc [new file with mode: 0644]
tests/lac/solver_control_01.output [new file with mode: 0644]
tests/lac/solver_control_02.cc [new file with mode: 0644]
tests/lac/solver_control_02.output [new file with mode: 0644]
tests/lac/solver_control_03.cc [new file with mode: 0644]
tests/lac/solver_control_03.output [new file with mode: 0644]

diff --git a/tests/lac/solver_control_01.cc b/tests/lac/solver_control_01.cc
new file mode 100644 (file)
index 0000000..42f765f
--- /dev/null
@@ -0,0 +1,110 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// test SolverControl
+// This test is adapted from tests/trilinos/solver_03.cc
+
+
+#include "../tests.h"
+#include <deal.II/base/utilities.h>
+#include "../testmatrix.h"
+#include <cmath>
+#include <fstream>
+#include <iostream>
+#include <iomanip>
+#include <deal.II/lac/sparse_matrix.h>
+#include <deal.II/lac/vector.h>
+#include <deal.II/lac/vector_memory.h>
+#include <deal.II/lac/solver_control.h>
+#include <deal.II/lac/solver_richardson.h>
+#include <deal.II/lac/solver_relaxation.h>
+#include <deal.II/lac/precondition.h>
+
+template<typename MatrixType, typename VectorType, class PRECONDITION>
+void
+check_solve (SolverControl       &solver_control,
+             const MatrixType    &A,
+             VectorType          &u,
+             VectorType          &f,
+             const PRECONDITION  &P,
+             const bool           expected_result)
+{
+  SolverCG<VectorType> solver(solver_control);
+
+  u = 0.;
+  f = 1.;
+  bool success = false;
+  try
+    {
+      solver.solve(A,u,f,P);
+      deallog << "Success. ";
+      success = true;
+    }
+  catch (std::exception &e)
+    {
+      deallog << "Failure. ";
+    }
+
+  deallog << "Solver stopped after " << solver_control.last_step()
+          << " iterations" << std::endl;
+  Assert(success == expected_result, ExcMessage("Incorrect result."));
+}
+
+
+int main(int argc, char **argv)
+{
+  std::ofstream logfile("output");
+  logfile.precision(4);
+  deallog.attach(logfile);
+  deallog.threshold_double(1.e-10);
+
+  {
+    const unsigned int size = 32;
+    unsigned int dim = (size-1)*(size-1);
+
+    deallog << "Size " << size << " Unknowns " << dim << std::endl;
+
+    // Make matrix
+    FDMatrix testproblem(size, size);
+    SparsityPattern structure(dim, dim, 5);
+    testproblem.five_point_structure(structure);
+    structure.compress();
+    SparseMatrix<double>  A(structure);
+    testproblem.five_point(A);
+
+    Vector<double>  f(dim);
+    Vector<double>  u(dim);
+    f = 1.;
+
+    PreconditionJacobi<> preconditioner;
+    preconditioner.initialize(A);
+
+    deallog.push("Abs tol");
+    {
+      // Expects success
+      SolverControl solver_control(2000, 1.e-3);
+      check_solve (solver_control, A,u,f, preconditioner, true);
+    }
+    deallog.pop();
+    deallog.push("Iterations");
+    {
+      // Expects failure
+      SolverControl solver_control(20, 1.e-3);
+      check_solve (solver_control, A,u,f, preconditioner, false);
+    }
+    deallog.pop();
+  }
+}
diff --git a/tests/lac/solver_control_01.output b/tests/lac/solver_control_01.output
new file mode 100644 (file)
index 0000000..800f90b
--- /dev/null
@@ -0,0 +1,8 @@
+
+DEAL::Size 32 Unknowns 961
+DEAL:Abs tol:cg::Starting value 31.0000
+DEAL:Abs tol:cg::Convergence step 43 value 0.000818932
+DEAL:Abs tol::Success. Solver stopped after 43 iterations
+DEAL:Iterations:cg::Starting value 31.0000
+DEAL:Iterations:cg::Failure step 20 value 6.00198
+DEAL:Iterations::Failure. Solver stopped after 20 iterations
diff --git a/tests/lac/solver_control_02.cc b/tests/lac/solver_control_02.cc
new file mode 100644 (file)
index 0000000..9686acf
--- /dev/null
@@ -0,0 +1,118 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// test ReductionControl
+// This test is adapted from tests/trilinos/solver_03.cc
+
+
+#include "../tests.h"
+#include <deal.II/base/utilities.h>
+#include "../testmatrix.h"
+#include <cmath>
+#include <fstream>
+#include <iostream>
+#include <iomanip>
+#include <deal.II/lac/sparse_matrix.h>
+#include <deal.II/lac/vector.h>
+#include <deal.II/lac/vector_memory.h>
+#include <deal.II/lac/solver_control.h>
+#include <deal.II/lac/solver_richardson.h>
+#include <deal.II/lac/solver_relaxation.h>
+#include <deal.II/lac/precondition.h>
+
+template<typename MatrixType, typename VectorType, class PRECONDITION>
+void
+check_solve (SolverControl       &solver_control,
+             const MatrixType    &A,
+             VectorType          &u,
+             VectorType          &f,
+             const PRECONDITION  &P,
+             const bool           expected_result)
+{
+  SolverCG<VectorType> solver(solver_control);
+
+  u = 0.;
+  f = 1.;
+  bool success = false;
+  try
+    {
+      solver.solve(A,u,f,P);
+      deallog << "Success. ";
+      success = true;
+    }
+  catch (std::exception &e)
+    {
+      deallog << "Failure. ";
+    }
+
+  deallog << "Solver stopped after " << solver_control.last_step()
+          << " iterations" << std::endl;
+  Assert(success == expected_result, ExcMessage("Incorrect result."));
+}
+
+
+int main(int argc, char **argv)
+{
+  std::ofstream logfile("output");
+  logfile.precision(4);
+  deallog.attach(logfile);
+  deallog.threshold_double(1.e-10);
+
+
+  {
+    const unsigned int size = 32;
+    unsigned int dim = (size-1)*(size-1);
+
+    deallog << "Size " << size << " Unknowns " << dim << std::endl;
+
+    // Make matrix
+    FDMatrix testproblem(size, size);
+    SparsityPattern structure(dim, dim, 5);
+    testproblem.five_point_structure(structure);
+    structure.compress();
+    SparseMatrix<double>  A(structure);
+    testproblem.five_point(A);
+
+    Vector<double>  f(dim);
+    Vector<double>  u(dim);
+    f = 1.;
+
+    PreconditionJacobi<> preconditioner;
+    preconditioner.initialize(A);
+
+    deallog.push("Rel tol");
+    {
+      // Expects success
+      ReductionControl solver_control(2000, 1.e-30, 1e-6);
+      check_solve (solver_control, A,u,f, preconditioner, true);
+    }
+    deallog.pop();
+    deallog.push("Abs tol");
+    {
+      // Expects success
+      ReductionControl solver_control(2000, 1.e-3, 1e-9);
+      check_solve (solver_control, A,u,f, preconditioner, true);
+    }
+    deallog.pop();
+    deallog.push("Iterations");
+    {
+      // Expects failure
+      ReductionControl solver_control(20, 1.e-30, 1e-6);
+      check_solve (solver_control, A,u,f, preconditioner, false);
+    }
+    deallog.pop();
+  }
+}
diff --git a/tests/lac/solver_control_02.output b/tests/lac/solver_control_02.output
new file mode 100644 (file)
index 0000000..4ddd995
--- /dev/null
@@ -0,0 +1,11 @@
+
+DEAL::Size 32 Unknowns 961
+DEAL:Rel tol:cg::Starting value 31.0000
+DEAL:Rel tol:cg::Convergence step 50 value 2.11364e-05
+DEAL:Rel tol::Success. Solver stopped after 50 iterations
+DEAL:Abs tol:cg::Starting value 31.0000
+DEAL:Abs tol:cg::Convergence step 43 value 0.000818932
+DEAL:Abs tol::Success. Solver stopped after 43 iterations
+DEAL:Iterations:cg::Starting value 31.0000
+DEAL:Iterations:cg::Failure step 20 value 6.00198
+DEAL:Iterations::Failure. Solver stopped after 20 iterations
diff --git a/tests/lac/solver_control_03.cc b/tests/lac/solver_control_03.cc
new file mode 100644 (file)
index 0000000..85cb45d
--- /dev/null
@@ -0,0 +1,111 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2017 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 at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// test IterationNumberControl
+// This test is adapted from tests/trilinos/solver_03.cc
+
+
+#include "../tests.h"
+#include <deal.II/base/utilities.h>
+#include "../testmatrix.h"
+#include <cmath>
+#include <fstream>
+#include <iostream>
+#include <iomanip>
+#include <deal.II/lac/sparse_matrix.h>
+#include <deal.II/lac/vector.h>
+#include <deal.II/lac/vector_memory.h>
+#include <deal.II/lac/solver_control.h>
+#include <deal.II/lac/solver_richardson.h>
+#include <deal.II/lac/solver_relaxation.h>
+#include <deal.II/lac/precondition.h>
+
+template<typename MatrixType, typename VectorType, class PRECONDITION>
+void
+check_solve (SolverControl       &solver_control,
+             const MatrixType    &A,
+             VectorType          &u,
+             VectorType          &f,
+             const PRECONDITION  &P,
+             const bool           expected_result)
+{
+  SolverCG<VectorType> solver(solver_control);
+
+  u = 0.;
+  f = 1.;
+  bool success = false;
+  try
+    {
+      solver.solve(A,u,f,P);
+      deallog << "Success. ";
+      success = true;
+    }
+  catch (std::exception &e)
+    {
+      deallog << "Failure. ";
+    }
+
+  deallog << "Solver stopped after " << solver_control.last_step()
+          << " iterations" << std::endl;
+  Assert(success == expected_result, ExcMessage("Incorrect result."));
+}
+
+
+int main(int argc, char **argv)
+{
+  std::ofstream logfile("output");
+  logfile.precision(4);
+  deallog.attach(logfile);
+  deallog.threshold_double(1.e-10);
+
+
+  {
+    const unsigned int size = 32;
+    unsigned int dim = (size-1)*(size-1);
+
+    deallog << "Size " << size << " Unknowns " << dim << std::endl;
+
+    // Make matrix
+    FDMatrix testproblem(size, size);
+    SparsityPattern structure(dim, dim, 5);
+    testproblem.five_point_structure(structure);
+    structure.compress();
+    SparseMatrix<double>  A(structure);
+    testproblem.five_point(A);
+
+    Vector<double>  f(dim);
+    Vector<double>  u(dim);
+    f = 1.;
+
+    PreconditionJacobi<> preconditioner;
+    preconditioner.initialize(A);
+
+    deallog.push("Abs tol");
+    {
+      // Expects success
+      IterationNumberControl solver_control(2000, 1.e-3);
+      check_solve (solver_control, A,u,f, preconditioner, true);
+    }
+    deallog.pop();
+    deallog.push("Iterations");
+    {
+      // Expects success
+      IterationNumberControl solver_control(20, 1.e-3);
+      check_solve (solver_control, A,u,f, preconditioner, true);
+    }
+    deallog.pop();
+  }
+}
diff --git a/tests/lac/solver_control_03.output b/tests/lac/solver_control_03.output
new file mode 100644 (file)
index 0000000..eb675db
--- /dev/null
@@ -0,0 +1,8 @@
+
+DEAL::Size 32 Unknowns 961
+DEAL:Abs tol:cg::Starting value 31.0000
+DEAL:Abs tol:cg::Convergence step 43 value 0.000818932
+DEAL:Abs tol::Success. Solver stopped after 43 iterations
+DEAL:Iterations:cg::Starting value 31.0000
+DEAL:Iterations:cg::Convergence step 20 value 6.00198
+DEAL:Iterations::Success. Solver stopped after 20 iterations

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.