Reference documentation for deal.II version 9.2.0
\(\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\}}\)
kinsol.h
Go to the documentation of this file.
1 //-----------------------------------------------------------
2 //
3 // Copyright (C) 2017 - 2020 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 #ifndef dealii_sundials_kinsol_h
17 #define dealii_sundials_kinsol_h
18 
19 #include <deal.II/base/config.h>
20 #ifdef DEAL_II_WITH_SUNDIALS
21 
23 # include <deal.II/base/exceptions.h>
24 # include <deal.II/base/logstream.h>
25 # include <deal.II/base/mpi.h>
27 
28 # include <deal.II/lac/vector.h>
30 
31 # include <boost/signals2.hpp>
32 
33 # include <kinsol/kinsol.h>
34 # include <kinsol/kinsol_impl.h>
35 # include <nvector/nvector_serial.h>
36 # include <sundials/sundials_math.h>
37 # include <sundials/sundials_types.h>
38 
39 # include <memory>
40 
41 
43 
44 // Shorthand notation for KINSOL error codes.
45 # define AssertKINSOL(code) Assert(code >= 0, ExcKINSOLError(code))
46 
47 namespace SUNDIALS
48 {
197  template <typename VectorType = Vector<double>>
198  class KINSOL
199  {
200  public:
205  {
206  public:
212  {
216  newton = KIN_NONE,
220  linesearch = KIN_LINESEARCH,
224  fixed_point = KIN_FP,
228  picard = KIN_PICARD,
229  };
230 
260  // Global parameters
262  const unsigned int maximum_non_linear_iterations = 200,
263  const double function_tolerance = 0.0,
264  const double step_tolerance = 0.0,
265  const bool no_init_setup = false,
266  const unsigned int maximum_setup_calls = 0,
267  const double maximum_newton_step = 0.0,
268  const double dq_relative_error = 0.0,
269  const unsigned int maximum_beta_failures = 0,
270  const unsigned int anderson_subspace_size = 0)
271  : strategy(strategy)
281  {}
282 
321  void
323  {
324  static std::string strategy_str("newton");
325  prm.add_parameter("Solution strategy",
326  strategy_str,
327  "Choose among newton|linesearch|fixed_point|picard",
329  "newton|linesearch|fixed_point|picard"));
330  prm.add_action("Solution strategy", [&](const std::string &value) {
331  if (value == "newton")
332  strategy = newton;
333  else if (value == "linesearch")
335  else if (value == "fixed_point")
337  else if (value == "picard")
338  strategy = picard;
339  else
340  Assert(false, ExcInternalError());
341  });
342  prm.add_parameter("Maximum number of nonlinear iterations",
344  prm.add_parameter("Function norm stopping tolerance",
346  prm.add_parameter("Scaled step stopping tolerance", step_tolerance);
347 
348  prm.enter_subsection("Newton parameters");
349  prm.add_parameter("No initial matrix setup", no_init_setup);
350  prm.add_parameter("Maximum iterations without matrix setup",
352  prm.add_parameter("Maximum allowable scaled length of the Newton step",
354  prm.add_parameter("Relative error for different quotient computation",
356  prm.leave_subsection();
357 
358  prm.enter_subsection("Linesearch parameters");
359  prm.add_parameter("Maximum number of beta-condition failures",
361  prm.leave_subsection();
362 
363 
364  prm.enter_subsection("Fixed point and Picard parameters");
365  prm.add_parameter("Anderson acceleration subspace size",
367  prm.leave_subsection();
368  }
369 
378 
383 
391 
399 
409 
416  unsigned int maximum_setup_calls;
417 
424 
433 
439  unsigned int maximum_beta_failures;
440 
448  };
449 
459  const MPI_Comm mpi_comm = MPI_COMM_WORLD);
460 
464  ~KINSOL();
465 
471  unsigned int
472  solve(VectorType &initial_guess_and_solution);
473 
478  std::function<void(VectorType &)> reinit_vector;
479 
492  std::function<int(const VectorType &src, VectorType &dst)> residual;
493 
507  std::function<int(const VectorType &src, VectorType &dst)>
509 
510 
551  std::function<int(const VectorType &current_u, const VectorType &current_f)>
553 
591  std::function<int(const VectorType &ycur,
592  const VectorType &fcur,
593  const VectorType &rhs,
594  VectorType & dst)>
596 
603  std::function<VectorType &()> get_solution_scaling;
604 
612  std::function<VectorType &()> get_function_scaling;
613 
618  int,
619  << "One of the SUNDIALS KINSOL internal functions "
620  << " returned a negative error code: " << arg1
621  << ". Please consult SUNDIALS manual.");
622 
623 
624  private:
630  std::string,
631  << "Please provide an implementation for the function \""
632  << arg1 << "\"");
633 
639  void
641 
646 
650  void *kinsol_mem;
651 
655  N_Vector solution;
656 
660  N_Vector u_scale;
661 
665  N_Vector f_scale;
666 
671 
676  };
677 
678 } // namespace SUNDIALS
679 
680 
682 #endif
683 
684 
685 #endif
SUNDIALS::KINSOL::setup_jacobian
std::function< int(const VectorType &current_u, const VectorType &current_f)> setup_jacobian
Definition: kinsol.h:552
SUNDIALS::KINSOL::iteration_function
std::function< int(const VectorType &src, VectorType &dst)> iteration_function
Definition: kinsol.h:508
SUNDIALS::KINSOL::AdditionalData::function_tolerance
double function_tolerance
Definition: kinsol.h:390
SUNDIALS::KINSOL::AdditionalData::maximum_beta_failures
unsigned int maximum_beta_failures
Definition: kinsol.h:439
SUNDIALS::KINSOL::data
AdditionalData data
Definition: kinsol.h:645
SUNDIALS::KINSOL
Definition: kinsol.h:198
SUNDIALS::KINSOL::AdditionalData::anderson_subspace_size
unsigned int anderson_subspace_size
Definition: kinsol.h:447
SUNDIALS::KINSOL::u_scale
N_Vector u_scale
Definition: kinsol.h:660
SUNDIALS::KINSOL::f_scale
N_Vector f_scale
Definition: kinsol.h:665
SUNDIALS::KINSOL::AdditionalData::newton
@ newton
Definition: kinsol.h:216
VectorType
MPI_Comm
SUNDIALS::KINSOL::solution
N_Vector solution
Definition: kinsol.h:655
SUNDIALS
Definition: arkode.h:58
SUNDIALS::KINSOL::AdditionalData
Definition: kinsol.h:204
GrowingVectorMemory
Definition: vector_memory.h:320
conditional_ostream.h
SUNDIALS::KINSOL::AdditionalData::linesearch
@ linesearch
Definition: kinsol.h:220
SUNDIALS::KINSOL::solve_jacobian_system
std::function< int(const VectorType &ycur, const VectorType &fcur, const VectorType &rhs, VectorType &dst)> solve_jacobian_system
Definition: kinsol.h:595
SUNDIALS::KINSOL::set_functions_to_trigger_an_assert
void set_functions_to_trigger_an_assert()
Definition: kinsol.cc:344
SUNDIALS::KINSOL::ExcFunctionNotProvided
static ::ExceptionBase & ExcFunctionNotProvided(std::string arg1)
SUNDIALS::KINSOL::kinsol_mem
void * kinsol_mem
Definition: kinsol.h:650
mpi.h
SUNDIALS::KINSOL::AdditionalData::strategy
SolutionStrategy strategy
Definition: kinsol.h:377
SUNDIALS::KINSOL::get_function_scaling
std::function< VectorType &()> get_function_scaling
Definition: kinsol.h:612
SUNDIALS::KINSOL::communicator
MPI_Comm communicator
Definition: kinsol.h:670
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
SUNDIALS::KINSOL::AdditionalData::fixed_point
@ fixed_point
Definition: kinsol.h:224
SUNDIALS::KINSOL::AdditionalData::dq_relative_error
double dq_relative_error
Definition: kinsol.h:432
SUNDIALS::KINSOL::ExcKINSOLError
static ::ExceptionBase & ExcKINSOLError(int arg1)
SUNDIALS::KINSOL::~KINSOL
~KINSOL()
Definition: kinsol.cc:170
DeclException1
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:518
SUNDIALS::KINSOL::AdditionalData::picard
@ picard
Definition: kinsol.h:228
parameter_handler.h
ParameterHandler::add_action
void add_action(const std::string &entry, const std::function< void(const std::string &value)> &action)
Definition: parameter_handler.cc:830
exceptions.h
SUNDIALS::KINSOL::mem
GrowingVectorMemory< VectorType > mem
Definition: kinsol.h:675
SUNDIALS::KINSOL::AdditionalData::AdditionalData
AdditionalData(const SolutionStrategy &strategy=linesearch, const unsigned int maximum_non_linear_iterations=200, const double function_tolerance=0.0, const double step_tolerance=0.0, const bool no_init_setup=false, const unsigned int maximum_setup_calls=0, const double maximum_newton_step=0.0, const double dq_relative_error=0.0, const unsigned int maximum_beta_failures=0, const unsigned int anderson_subspace_size=0)
Definition: kinsol.h:259
SUNDIALS::KINSOL::get_solution_scaling
std::function< VectorType &()> get_solution_scaling
Definition: kinsol.h:603
value
static const bool value
Definition: dof_tools_constraints.cc:433
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
SUNDIALS::KINSOL::solve
unsigned int solve(VectorType &initial_guess_and_solution)
Definition: kinsol.cc:188
SUNDIALS::KINSOL::AdditionalData::SolutionStrategy
SolutionStrategy
Definition: kinsol.h:211
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
SUNDIALS::KINSOL::AdditionalData::step_tolerance
double step_tolerance
Definition: kinsol.h:398
vector.h
SUNDIALS::KINSOL::AdditionalData::no_init_setup
bool no_init_setup
Definition: kinsol.h:408
ParameterHandler::add_parameter
void add_parameter(const std::string &entry, ParameterType &parameter, const std::string &documentation="", const Patterns::PatternBase &pattern= *Patterns::Tools::Convert< ParameterType >::to_pattern(), const bool has_to_be_set=false)
Definition: parameter_handler.h:2328
SUNDIALS::KINSOL::AdditionalData::maximum_newton_step
double maximum_newton_step
Definition: kinsol.h:423
vector_memory.h
SUNDIALS::KINSOL::KINSOL
KINSOL(const AdditionalData &data=AdditionalData(), const MPI_Comm mpi_comm=MPI_COMM_WORLD)
Definition: kinsol.cc:153
ParameterHandler
Definition: parameter_handler.h:845
config.h
ParameterHandler::enter_subsection
void enter_subsection(const std::string &subsection)
Definition: parameter_handler.cc:927
SUNDIALS::KINSOL::reinit_vector
std::function< void(VectorType &)> reinit_vector
Definition: kinsol.h:478
Patterns::Selection
Definition: patterns.h:383
SUNDIALS::KINSOL::AdditionalData::maximum_non_linear_iterations
unsigned int maximum_non_linear_iterations
Definition: kinsol.h:382
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
ParameterHandler::leave_subsection
void leave_subsection()
Definition: parameter_handler.cc:941
logstream.h
SUNDIALS::KINSOL::AdditionalData::add_parameters
void add_parameters(ParameterHandler &prm)
Definition: kinsol.h:322
SUNDIALS::KINSOL::residual
std::function< int(const VectorType &src, VectorType &dst)> residual
Definition: kinsol.h:492
int
SUNDIALS::KINSOL::AdditionalData::maximum_setup_calls
unsigned int maximum_setup_calls
Definition: kinsol.h:416