|
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\}}\)
Go to the documentation of this file.
16 #ifndef dealii_sundials_kinsol_h
17 #define dealii_sundials_kinsol_h
20 #ifdef DEAL_II_WITH_SUNDIALS
31 # include <boost/signals2.hpp>
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>
45 # define AssertKINSOL(code) Assert(code >= 0, ExcKINSOLError(code))
197 template <
typename VectorType = Vector<
double>>
324 static std::string strategy_str(
"newton");
327 "Choose among newton|linesearch|fixed_point|picard",
329 "newton|linesearch|fixed_point|picard"));
331 if (
value ==
"newton")
333 else if (
value ==
"linesearch")
335 else if (
value ==
"fixed_point")
337 else if (
value ==
"picard")
352 prm.
add_parameter(
"Maximum allowable scaled length of the Newton step",
354 prm.
add_parameter(
"Relative error for different quotient computation",
359 prm.
add_parameter(
"Maximum number of beta-condition failures",
459 const MPI_Comm mpi_comm = MPI_COMM_WORLD);
619 <<
"One of the SUNDIALS KINSOL internal functions "
620 <<
" returned a negative error code: " << arg1
621 <<
". Please consult SUNDIALS manual.");
631 <<
"Please provide an implementation for the function \""
std::function< int(const VectorType ¤t_u, const VectorType ¤t_f)> setup_jacobian
std::function< int(const VectorType &src, VectorType &dst)> iteration_function
double function_tolerance
unsigned int maximum_beta_failures
unsigned int anderson_subspace_size
std::function< int(const VectorType &ycur, const VectorType &fcur, const VectorType &rhs, VectorType &dst)> solve_jacobian_system
void set_functions_to_trigger_an_assert()
static ::ExceptionBase & ExcFunctionNotProvided(std::string arg1)
SolutionStrategy strategy
std::function< VectorType &()> get_function_scaling
#define DEAL_II_NAMESPACE_OPEN
static ::ExceptionBase & ExcKINSOLError(int arg1)
#define DeclException1(Exception1, type1, outsequence)
void add_action(const std::string &entry, const std::function< void(const std::string &value)> &action)
GrowingVectorMemory< VectorType > mem
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)
std::function< VectorType &()> get_solution_scaling
static ::ExceptionBase & ExcInternalError()
unsigned int solve(VectorType &initial_guess_and_solution)
#define Assert(cond, exc)
void add_parameter(const std::string &entry, ParameterType ¶meter, const std::string &documentation="", const Patterns::PatternBase &pattern= *Patterns::Tools::Convert< ParameterType >::to_pattern(), const bool has_to_be_set=false)
double maximum_newton_step
KINSOL(const AdditionalData &data=AdditionalData(), const MPI_Comm mpi_comm=MPI_COMM_WORLD)
void enter_subsection(const std::string &subsection)
std::function< void(VectorType &)> reinit_vector
unsigned int maximum_non_linear_iterations
#define DEAL_II_NAMESPACE_CLOSE
void add_parameters(ParameterHandler &prm)
std::function< int(const VectorType &src, VectorType &dst)> residual
unsigned int maximum_setup_calls