]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add a test with an irrecoverable failure of KINSOL. 15189/head
authorWolfgang Bangerth <bangerth@colostate.edu>
Tue, 9 May 2023 16:11:03 +0000 (10:11 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Tue, 9 May 2023 20:16:39 +0000 (14:16 -0600)
tests/sundials/kinsol_07.cc [new file with mode: 0644]
tests/sundials/kinsol_07.output [new file with mode: 0644]

diff --git a/tests/sundials/kinsol_07.cc b/tests/sundials/kinsol_07.cc
new file mode 100644 (file)
index 0000000..b45c69d
--- /dev/null
@@ -0,0 +1,99 @@
+//-----------------------------------------------------------
+//
+//    Copyright (C) 2017 - 2021 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.
+//
+//-----------------------------------------------------------
+
+#include <deal.II/base/parameter_handler.h>
+
+#include <deal.II/lac/full_matrix.h>
+#include <deal.II/lac/vector.h>
+
+#include <deal.II/sundials/kinsol.h>
+
+#include "../tests.h"
+
+// This is a variation of the kinsol_06 test, which throws
+// "recoverable" errors and forces KINSOL to backtrack. This test
+// instead throws an irrecoverable error and forces KINSOL to
+// abort. The test checks that the thrown exception is re-thrown after
+// KINSOL returns.
+
+
+int
+main()
+{
+  initlog();
+
+  using VectorType = Vector<double>;
+
+  // Size of the problem
+  const unsigned int N = 1;
+
+  SUNDIALS::KINSOL<VectorType>::AdditionalData data;
+  ParameterHandler                             prm;
+  data.add_parameters(prm);
+
+  std::ifstream ifile(SOURCE_DIR "/kinsol_newton.prm");
+  prm.parse_input(ifile);
+
+  SUNDIALS::KINSOL<VectorType> kinsol(data);
+
+  kinsol.reinit_vector = [N](VectorType &v) { v.reinit(N); };
+
+  kinsol.residual = [&](const VectorType &u, VectorType &F) {
+    deallog << "Computing residual at " << u[0] << std::endl;
+
+    if ((u[0] < -10) || (u[0] > 20))
+      {
+        deallog << "Reporting irrecoverable failure." << std::endl;
+        throw ExcMessage("Irrecoverable failure.");
+      }
+
+
+    F.reinit(u);
+    F[0] = std::atan(u[0]) - 0.5;
+  };
+
+  double J_inverse;
+
+  kinsol.setup_jacobian = [&J_inverse](const VectorType &u,
+                                       const VectorType &F) {
+    deallog << "Setting up Jacobian system at u=" << u[0] << std::endl;
+
+    const double J = 1. / (1 + u[0] * u[0]);
+    J_inverse      = 1. / J;
+  };
+
+
+  kinsol.solve_with_jacobian = [&](const VectorType &rhs,
+                                   VectorType &      dst,
+                                   double) { dst[0] = J_inverse * rhs[0]; };
+
+  VectorType v(N);
+  v[0] = 10;
+
+  try
+    {
+      auto niter = kinsol.solve(v);
+
+      deallog << "Found solution at " << std::flush;
+      v.print(deallog.get_file_stream());
+      deallog << "Converged in " << niter << " iterations." << std::endl;
+    }
+
+  catch (const std::exception &exc)
+    {
+      deallog << "Caught an irrecoverable exception in a callback:" << std::endl
+              << exc.what() << std::endl;
+    }
+}
diff --git a/tests/sundials/kinsol_07.output b/tests/sundials/kinsol_07.output
new file mode 100644 (file)
index 0000000..ba030d4
--- /dev/null
@@ -0,0 +1,17 @@
+
+DEAL::Computing residual at 10.0000
+DEAL::Setting up Jacobian system at u=10.0000
+DEAL::Computing residual at 10.0000
+DEAL::Computing residual at -88.0839
+DEAL::Reporting irrecoverable failure.
+DEAL::Caught an irrecoverable exception in a callback:
+DEAL::
+--------------------------------------------------------
+An error occurred in line <0> of file <> in function
+    
+The violated condition was: 
+    
+Additional information: 
+    Irrecoverable failure.
+--------------------------------------------------------
+

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.