17#ifndef dealii_sundials_ida_h
18#define dealii_sundials_ida_h
23#ifdef DEAL_II_WITH_SUNDIALS
28# ifdef DEAL_II_WITH_PETSC
35# ifdef DEAL_II_SUNDIALS_WITH_IDAS
36# include <idas/idas.h>
43# include <boost/signals2.hpp>
45# include <nvector/nvector_serial.h>
46# include <sundials/sundials_config.h>
47# include <sundials/sundials_math.h>
48# include <sundials/sundials_types.h>
56# define AssertIDA(code) Assert(code >= 0, ExcIDAError(code))
233 template <
typename VectorType = Vector<
double>>
402 "Ignore algebraic terms for error computations",
404 "Indicate whether or not to suppress algebraic variables "
405 "in the local error test.");
409 static std::string ic_type_str =
"use_y_diff";
411 "Correction type at initial time",
413 "This is one of the following three options for the "
414 "initial condition calculation. \n"
415 " none: do not try to make initial conditions consistent. \n"
416 " use_y_diff: compute the algebraic components of y and differential\n"
417 " components of y_dot, given the differential components of y. \n"
418 " This option requires that the user specifies differential and \n"
419 " algebraic components in the function get_differential_components.\n"
420 " use_y_dot: compute all components of y, given y_dot.",
422 prm.
add_action(
"Correction type at initial time",
423 [&](
const std::string &value) {
424 if (value ==
"use_y_diff")
426 else if (value ==
"use_y_dot")
428 else if (value ==
"none")
434 static std::string reset_type_str =
"use_y_diff";
436 "Correction type after restart",
438 "This is one of the following three options for the "
439 "initial condition calculation. \n"
440 " none: do not try to make initial conditions consistent. \n"
441 " use_y_diff: compute the algebraic components of y and differential\n"
442 " components of y_dot, given the differential components of y. \n"
443 " This option requires that the user specifies differential and \n"
444 " algebraic components in the function get_differential_components.\n"
445 " use_y_dot: compute all components of y, given y_dot.",
447 prm.
add_action(
"Correction type after restart",
448 [&](
const std::string &value) {
449 if (value ==
"use_y_diff")
451 else if (value ==
"use_y_dot")
453 else if (value ==
"none")
461 "Factor to use when converting from the integrator tolerance to the linear solver tolerance",
619 solve_dae(VectorType &solution, VectorType &solution_dot);
643 reset(
const double t,
const double h, VectorType &y, VectorType &yp);
660 std::function<
int(
const double t,
662 const VectorType &y_dot,
699 std::function<
int(
const double t,
701 const VectorType &y_dot,
737 std::function<
int(
const VectorType &rhs, VectorType &dst)>
782 int(
const VectorType &rhs, VectorType &dst,
const double tolerance)>
799 std::function<void(
const double t,
800 const VectorType & sol,
801 const VectorType & sol_dot,
802 const unsigned int step_number)>
821 std::function<
bool(
const double t, VectorType &sol, VectorType &sol_dot)>
853 <<
"One of the SUNDIALS IDA internal functions "
854 <<
" returned a negative error code: " << arg1
855 <<
". Please consult SUNDIALS manual.");
865 <<
"Please provide an implementation for the function \""
886# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
904# ifdef DEAL_II_WITH_PETSC
905# ifdef PETSC_USE_COMPLEX
906 static_assert(!std::is_same<VectorType, PETScWrappers::MPI::Vector>::value,
907 "Sundials does not support complex scalar types, "
908 "but PETSc is configured to use a complex scalar type!");
911 !std::is_same<VectorType, PETScWrappers::MPI::BlockVector>::value,
912 "Sundials does not support complex scalar types, "
913 "but PETSc is configured to use a complex scalar type!");
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)
void add_action(const std::string &entry, const std::function< void(const std::string &value)> &action)
void enter_subsection(const std::string &subsection)
double relative_tolerance
InitialConditionCorrection ic_type
bool ignore_algebraic_terms_for_errors
void add_parameters(ParameterHandler &prm)
AdditionalData(const double initial_time=0.0, const double final_time=1.0, const double initial_step_size=1e-2, const double output_period=1e-1, const double minimum_step_size=1e-6, const unsigned int maximum_order=5, const unsigned int maximum_non_linear_iterations=10, const double ls_norm_factor=0, const double absolute_tolerance=1e-6, const double relative_tolerance=1e-5, const bool ignore_algebraic_terms_for_errors=true, const InitialConditionCorrection &ic_type=use_y_diff, const InitialConditionCorrection &reset_type=use_y_diff, const unsigned int maximum_non_linear_iterations_ic=5)
unsigned maximum_non_linear_iterations_ic
InitialConditionCorrection reset_type
unsigned int maximum_order
InitialConditionCorrection
unsigned int maximum_non_linear_iterations
double absolute_tolerance
DEAL_II_DEPRECATED std::function< int(const VectorType &rhs, VectorType &dst)> solve_jacobian_system
unsigned int solve_dae(VectorType &solution, VectorType &solution_dot)
MPI_Comm mpi_communicator
void set_functions_to_trigger_an_assert()
const AdditionalData data
void reset(const double t, const double h, VectorType &y, VectorType &yp)
std::function< IndexSet()> differential_components
DeclException1(ExcIDAError, int,<< "One of the SUNDIALS IDA internal functions "<< " returned a negative error code: "<< arg1<< ". Please consult SUNDIALS manual.")
std::function< int(const double t, const VectorType &y, const VectorType &y_dot, VectorType &res)> residual
std::function< VectorType &()> get_local_tolerances
std::function< void(const double t, const VectorType &sol, const VectorType &sol_dot, const unsigned int step_number)> output_step
std::function< bool(const double t, VectorType &sol, VectorType &sol_dot)> solver_should_restart
std::function< int(const double t, const VectorType &y, const VectorType &y_dot, const double alpha)> setup_jacobian
DeclException1(ExcFunctionNotProvided, std::string,<< "Please provide an implementation for the function \""<< arg1<< "\"")
GrowingVectorMemory< VectorType > mem
std::function< void(VectorType &)> reinit_vector
std::function< int(const VectorType &rhs, VectorType &dst, const double tolerance)> solve_with_jacobian
#define DEAL_II_DEPRECATED
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcInternalError()
#define AssertThrow(cond, exc)