17#ifndef dealii_sundials_arkode_h
18#define dealii_sundials_arkode_h
23#ifdef DEAL_II_WITH_SUNDIALS
31# ifdef DEAL_II_WITH_PETSC
38# include <arkode/arkode.h>
39# include <nvector/nvector_serial.h>
40# ifdef DEAL_II_WITH_MPI
41# include <nvector/nvector_parallel.h>
48# include <boost/signals2.hpp>
50# include <sundials/sundials_linearsolver.h>
51# include <sundials/sundials_math.h>
52# include <sundials/sundials_types.h>
62# define AssertARKode(code) Assert(code >= 0, ExcARKodeError(code))
329 template <
typename VectorType = Vector<
double>>
538 const double intermediate_time,
539 const bool reset_solver =
false);
556 reset(
const double t,
const double h,
const VectorType &y);
594 void(
const double t,
const VectorType &y, VectorType &explicit_f)>
613 std::function<void(
const double t,
const VectorType &y, VectorType &res)>
633 std::function<void(
const double t,
const VectorType &v, VectorType &Mv)>
698 std::function<void(
const VectorType &v,
702 const VectorType &fy)>
740 void(
const double t,
const VectorType &y,
const VectorType &fy)>
826 std::function<void(
const double t,
828 const VectorType &fy,
878 std::function<void(
const double t,
880 const VectorType &fy,
913 std::function<void(
const double t,
960 std::function<void(
const double t,
961 const VectorType & sol,
962 const unsigned int step_number)>
1025 <<
"Please provide an implementation for the function \""
1034 const bool do_reset);
1072# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
1091 std::unique_ptr<internal::LinearSolverWrapper<VectorType>>
mass_solver;
1100# ifdef DEAL_II_WITH_PETSC
1101# ifdef PETSC_USE_COMPLEX
1102 static_assert(!std::is_same<VectorType, PETScWrappers::MPI::Vector>::value,
1103 "Sundials does not support complex scalar types, "
1104 "but PETSc is configured to use a complex scalar type!");
1107 !std::is_same<VectorType, PETScWrappers::MPI::BlockVector>::value,
1108 "Sundials does not support complex scalar types, "
1109 "but PETSc is configured to use a complex scalar type!");
1120 <<
"One of the SUNDIALS ARKode internal functions "
1121 <<
" returned a negative error code: " << arg1
1122 <<
". Please consult SUNDIALS manual.");
1125 template <
typename VectorType>
1128 const double initial_time,
1129 const double final_time,
1130 const double initial_step_size,
1131 const double output_period,
1133 const double minimum_step_size,
1134 const unsigned int maximum_order,
1135 const unsigned int maximum_non_linear_iterations,
1136 const bool implicit_function_is_linear,
1137 const bool implicit_function_is_time_independent,
1138 const bool mass_is_time_independent,
1139 const int anderson_acceleration_subspace,
1141 const double absolute_tolerance,
1142 const double relative_tolerance)
1143 : initial_time(initial_time)
1144 , final_time(final_time)
1145 , initial_step_size(initial_step_size)
1146 , minimum_step_size(minimum_step_size)
1147 , absolute_tolerance(absolute_tolerance)
1148 , relative_tolerance(relative_tolerance)
1149 , maximum_order(maximum_order)
1150 , output_period(output_period)
1151 , maximum_non_linear_iterations(maximum_non_linear_iterations)
1152 , implicit_function_is_linear(implicit_function_is_linear)
1153 , implicit_function_is_time_independent(
1154 implicit_function_is_time_independent)
1155 , mass_is_time_independent(mass_is_time_independent)
1156 , anderson_acceleration_subspace(anderson_acceleration_subspace)
1161 template <
typename VectorType>
1167 prm.
add_parameter(
"Time interval between each output", output_period);
1173 maximum_non_linear_iterations);
1175 implicit_function_is_linear);
1177 implicit_function_is_time_independent);
1178 prm.
add_parameter(
"Mass is time independent", mass_is_time_independent);
1180 anderson_acceleration_subspace);
1183 prm.
add_parameter(
"Absolute error tolerance", absolute_tolerance);
1184 prm.
add_parameter(
"Relative error tolerance", relative_tolerance);
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 enter_subsection(const std::string &subsection)
double absolute_tolerance
unsigned int maximum_order
bool implicit_function_is_linear
bool mass_is_time_independent
bool implicit_function_is_time_independent
double relative_tolerance
void add_parameters(ParameterHandler &prm)
int anderson_acceleration_subspace
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 bool implicit_function_is_linear=false, const bool implicit_function_is_time_independent=false, const bool mass_is_time_independent=false, const int anderson_acceleration_subspace=3, const double absolute_tolerance=1e-6, const double relative_tolerance=1e-5)
unsigned int maximum_non_linear_iterations
void setup_mass_solver(const VectorType &solution)
std::function< void(const double t)> mass_preconditioner_setup
std::function< void(const double t, const VectorType &y, VectorType &explicit_f)> explicit_function
std::exception_ptr pending_exception
std::function< void(const double t, const VectorType &y, const VectorType &fy, const int jok, int &jcur, const double gamma)> jacobian_preconditioner_setup
void * get_arkode_memory() const
std::function< VectorType &()> get_local_tolerances
std::function< void(const VectorType &v, VectorType &Jv, const double t, const VectorType &y, const VectorType &fy)> jacobian_times_vector
LinearSolveFunction< VectorType > solve_linearized_system
std::function< void(const double t, const VectorType &y, VectorType &res)> implicit_function
int do_evolve_time(VectorType &solution, ::DiscreteTime &time, const bool do_reset)
MPI_Comm mpi_communicator
void reset(const double t, const double h, const VectorType &y)
std::unique_ptr< internal::LinearSolverWrapper< VectorType > > mass_solver
std::function< void(const double t, const VectorType &y, const VectorType &fy, const VectorType &r, VectorType &z, const double gamma, const double tol, const int lr)> jacobian_preconditioner_solve
std::function< bool(const double t, VectorType &sol)> solver_should_restart
std::function< void(const double t)> mass_times_setup
std::function< void(const double t, const VectorType &sol, const unsigned int step_number)> output_step
void set_functions_to_trigger_an_assert()
std::unique_ptr< internal::LinearSolverWrapper< VectorType > > linear_solver
std::function< void(const double t, const VectorType &v, VectorType &Mv)> mass_times_vector
LinearSolveFunction< VectorType > solve_mass
std::function< void(const double t, const VectorType &r, VectorType &z, const double tol, const int lr)> mass_preconditioner_solve
std::function< void(const double t, const VectorType &y, const VectorType &fy)> jacobian_times_setup
unsigned int solve_ode(VectorType &solution)
void setup_system_solver(const VectorType &solution)
unsigned int solve_ode_incrementally(VectorType &solution, const double intermediate_time, const bool reset_solver=false)
std::function< void(void *arkode_mem)> custom_setup
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcARKodeError(int arg1)
static ::ExceptionBase & ExcFunctionNotProvided(std::string arg1)
#define DeclException1(Exception1, type1, outsequence)
std::function< void(SundialsOperator< VectorType > &op, SundialsPreconditioner< VectorType > &prec, VectorType &x, const VectorType &b, double tol)> LinearSolveFunction