]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add a KINSOL test.
authorWolfgang Bangerth <bangerth@colostate.edu>
Wed, 30 Apr 2025 19:58:16 +0000 (13:58 -0600)
committerWolfgang Bangerth <bangerth@colostate.edu>
Wed, 30 Apr 2025 22:15:43 +0000 (16:15 -0600)
tests/sundials/kinsol_06.cc
tests/sundials/kinsol_06_v2.cc [new file with mode: 0644]
tests/sundials/kinsol_06_v2.output [new file with mode: 0644]
tests/sundials/kinsol_06_v2_input.prm [new file with mode: 0644]

index 1d7f7e4e32c17caa96dd2b00ee3a3446ed39f370..5dcca728ae526584cbff7073465f80c45b8e02b0 100644 (file)
@@ -14,7 +14,6 @@
 
 #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>
diff --git a/tests/sundials/kinsol_06_v2.cc b/tests/sundials/kinsol_06_v2.cc
new file mode 100644 (file)
index 0000000..a384d44
--- /dev/null
@@ -0,0 +1,101 @@
+// ------------------------------------------------------------------------
+//
+// SPDX-License-Identifier: LGPL-2.1-or-later
+// Copyright (C) 2017 - 2023 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// Part of the source code is dual licensed under Apache-2.0 WITH
+// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
+// governing the source code and code contributions can be found in
+// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
+//
+// ------------------------------------------------------------------------
+
+
+#include <deal.II/base/parameter_handler.h>
+
+#include <deal.II/lac/vector.h>
+
+#include <deal.II/sundials/kinsol.h>
+
+#include "../tests.h"
+
+
+// Like the kinsol_06 test, but let the residual callback fail
+// consistently. KINSOL will just give up after a few tries, and we
+// should get a catchable exception.
+//
+// This testcase is a variation of the kinsol_06 test, modified by
+// Simon Wiesheier, and posted on the mailing list. Then further
+// adapted by Wolfgang Bangerth.
+
+
+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_06_v2_input.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) {
+    // Count how many times this function has been called. Let it call
+    // *every* time after the first attempt:
+    static int count = 0;
+    deallog << "Computing residual for the " << count + 1
+            << "th time, at u=" << u[0] << std::endl;
+    if ((u[0] < -10) || (u[0] > 20) || (count > 0))
+      {
+        deallog << "Reporting recoverable failure." << std::endl;
+        throw RecoverableUserCallbackError();
+      }
+    count++;
+
+
+    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);
+    }
+  catch (const std::exception &e)
+    {
+      deallog << "KINSOL threw an exception with the following message:"
+              << std::endl
+              << e.what() << std::endl;
+    }
+}
diff --git a/tests/sundials/kinsol_06_v2.output b/tests/sundials/kinsol_06_v2.output
new file mode 100644 (file)
index 0000000..4553ed2
--- /dev/null
@@ -0,0 +1,21 @@
+
+DEAL::Computing residual for the 1th time, at u=10.0000
+DEAL::Setting up Jacobian system at u=10.0000
+DEAL::Computing residual for the 2th time, at u=10.0000
+DEAL::Reporting recoverable failure.
+DEAL::KINSOL threw an exception with the following message:
+DEAL::
+--------------------------------------------------------
+An error occurred in line <0> of file <> in function
+    
+The violated condition was: 
+    
+Additional information: 
+    A user call-back function encountered a recoverable error, but the
+    underlying library that called the call-back did not manage to recover
+    from the error and aborted its operation.
+    
+    See the glossary entry on user call-back functions for more
+    information.
+--------------------------------------------------------
+
diff --git a/tests/sundials/kinsol_06_v2_input.prm b/tests/sundials/kinsol_06_v2_input.prm
new file mode 100644 (file)
index 0000000..fd9f762
--- /dev/null
@@ -0,0 +1,19 @@
+set Function norm stopping tolerance       = 0.000000
+set Maximum number of nonlinear iterations = 5
+set Scaled step stopping tolerance         = 0.000000
+set Solution strategy                      = newton
+
+subsection Fixed point and Picard parameters
+  set Anderson acceleration subspace size = 5
+end
+
+subsection Linesearch parameters
+  set Maximum number of beta-condition failures = 0
+end
+
+subsection Newton parameters
+  set Maximum allowable scaled length of the Newton step = 0.000000
+  set Maximum iterations without matrix setup            = 0
+  set No initial matrix setup                            = false
+  set Relative error for different quotient computation  = 0.000000
+end

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.