Reference documentation for deal.II version 9.4.1
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
kinsol.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2017 - 2022 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16
17#ifndef dealii_sundials_kinsol_h
18#define dealii_sundials_kinsol_h
19
20#include <deal.II/base/config.h>
21
22#ifdef DEAL_II_WITH_SUNDIALS
23
24
28# include <deal.II/base/mpi.h>
30
31# include <deal.II/lac/vector.h>
33
34# include <boost/signals2.hpp>
35
36# include <kinsol/kinsol.h>
37# include <nvector/nvector_serial.h>
38# include <sundials/sundials_math.h>
39# include <sundials/sundials_types.h>
40
41# include <memory>
42
43
45
46// Shorthand notation for KINSOL error codes.
47# define AssertKINSOL(code) Assert(code >= 0, ExcKINSOLError(code))
48
49namespace SUNDIALS
50{
180 template <typename VectorType = Vector<double>>
181 class KINSOL
182 {
183 public:
188 {
189 public:
195 {
199 newton = KIN_NONE,
203 linesearch = KIN_LINESEARCH,
207 fixed_point = KIN_FP,
211 picard = KIN_PICARD,
212 };
213
243 const unsigned int maximum_non_linear_iterations = 200,
244 const double function_tolerance = 0.0,
245 const double step_tolerance = 0.0,
246 const bool no_init_setup = false,
247 const unsigned int maximum_setup_calls = 0,
248 const double maximum_newton_step = 0.0,
249 const double dq_relative_error = 0.0,
250 const unsigned int maximum_beta_failures = 0,
251 const unsigned int anderson_subspace_size = 0);
252
291 void
293
302
307
315
323
333
343
350
359
366
374 };
375
386
394 KINSOL(const AdditionalData &data, const MPI_Comm &mpi_comm);
395
399 ~KINSOL();
400
406 unsigned int
407 solve(VectorType &initial_guess_and_solution);
408
416 std::function<void(VectorType &)> reinit_vector;
417
431 std::function<int(const VectorType &src, VectorType &dst)> residual;
432
447 std::function<int(const VectorType &src, VectorType &dst)>
449
495 std::function<int(const VectorType &current_u, const VectorType &current_f)>
497
556 std::function<int(const VectorType &ycur,
557 const VectorType &fcur,
558 const VectorType &rhs,
559 VectorType & dst)>
561
601 std::function<
602 int(const VectorType &rhs, VectorType &dst, const double tolerance)>
604
643 std::function<VectorType &()> get_solution_scaling;
644
659 std::function<VectorType &()> get_function_scaling;
660
665 int,
666 << "One of the SUNDIALS KINSOL internal functions "
667 << "returned a negative error code: " << arg1
668 << ". Please consult SUNDIALS manual.");
669
670
671 private:
677 std::string,
678 << "Please provide an implementation for the function \""
679 << arg1 << "\"");
680
686 void
688
693
698
703
704# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
708 SUNContext kinsol_ctx;
709# endif
710
714 N_Vector solution;
715
719 N_Vector u_scale;
720
724 N_Vector f_scale;
725
730 };
731
732} // namespace SUNDIALS
733
734
736
737#endif
738
739#endif
unsigned int maximum_setup_calls
Definition: kinsol.h:342
SolutionStrategy strategy
Definition: kinsol.h:301
unsigned int maximum_beta_failures
Definition: kinsol.h:365
void add_parameters(ParameterHandler &prm)
Definition: kinsol.cc:83
unsigned int maximum_non_linear_iterations
Definition: kinsol.h:306
unsigned int anderson_subspace_size
Definition: kinsol.h:373
std::function< VectorType &()> get_solution_scaling
Definition: kinsol.h:643
std::function< int(const VectorType &rhs, VectorType &dst, const double tolerance)> solve_with_jacobian
Definition: kinsol.h:603
std::function< int(const VectorType &current_u, const VectorType &current_f)> setup_jacobian
Definition: kinsol.h:496
void set_functions_to_trigger_an_assert()
Definition: kinsol.cc:624
MPI_Comm mpi_communicator
Definition: kinsol.h:697
N_Vector solution
Definition: kinsol.h:714
std::function< int(const VectorType &ycur, const VectorType &fcur, const VectorType &rhs, VectorType &dst)> solve_jacobian_system
Definition: kinsol.h:560
unsigned int solve(VectorType &initial_guess_and_solution)
Definition: kinsol.cc:351
SUNContext kinsol_ctx
Definition: kinsol.h:708
GrowingVectorMemory< VectorType > mem
Definition: kinsol.h:729
std::function< int(const VectorType &src, VectorType &dst)> residual
Definition: kinsol.h:431
N_Vector f_scale
Definition: kinsol.h:724
std::function< VectorType &()> get_function_scaling
Definition: kinsol.h:659
std::function< void(VectorType &)> reinit_vector
Definition: kinsol.h:416
AdditionalData data
Definition: kinsol.h:692
std::function< int(const VectorType &src, VectorType &dst)> iteration_function
Definition: kinsol.h:448
N_Vector u_scale
Definition: kinsol.h:719
void * kinsol_mem
Definition: kinsol.h:702
#define DEAL_II_DEPRECATED
Definition: config.h:164
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
static ::ExceptionBase & ExcKINSOLError(int arg1)
static ::ExceptionBase & ExcFunctionNotProvided(std::string arg1)
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:509