16#ifndef dealii_sundials_ida_h
17#define dealii_sundials_ida_h
22#ifdef DEAL_II_WITH_SUNDIALS
23# if DEAL_II_SUNDIALS_VERSION_LT(4, 0, 0)
29# ifdef DEAL_II_WITH_PETSC
36# ifdef DEAL_II_SUNDIALS_WITH_IDAS
37# include <idas/idas.h>
42# include <sundials/sundials_config.h>
43# if DEAL_II_SUNDIALS_VERSION_LT(3, 0, 0)
44# include <ida/ida_spbcgs.h>
45# include <ida/ida_spgmr.h>
46# include <ida/ida_sptfqmr.h>
48# include <boost/signals2.hpp>
50# include <nvector/nvector_serial.h>
51# include <sundials/sundials_math.h>
52# include <sundials/sundials_types.h>
60# define AssertIDA(code) Assert(code >= 0, ExcIDAError(code))
233 template <
typename VectorType = Vector<
double>>
252 enum InitialConditionCorrection
305 const double initial_time = 0.0,
306 const double final_time = 1.0,
307 const double initial_step_size = 1
e-2,
308 const double output_period = 1
e-1,
310 const double minimum_step_size = 1
e-6,
311 const unsigned int maximum_order = 5,
312 const unsigned int maximum_non_linear_iterations = 10,
314 const double absolute_tolerance = 1
e-6,
315 const double relative_tolerance = 1
e-5,
316 const bool ignore_algebraic_terms_for_errors =
true,
318 const InitialConditionCorrection &ic_type = use_y_diff,
319 const InitialConditionCorrection &reset_type = use_y_diff,
320 const unsigned int maximum_non_linear_iterations_ic = 5)
321 : initial_time(initial_time)
322 , final_time(final_time)
323 , initial_step_size(initial_step_size)
324 , minimum_step_size(minimum_step_size)
325 , absolute_tolerance(absolute_tolerance)
326 , relative_tolerance(relative_tolerance)
327 , maximum_order(maximum_order)
328 , output_period(output_period)
329 , ignore_algebraic_terms_for_errors(ignore_algebraic_terms_for_errors)
331 , reset_type(reset_type)
332 , maximum_non_linear_iterations_ic(maximum_non_linear_iterations_ic)
333 , maximum_non_linear_iterations(maximum_non_linear_iterations)
383 prm.
add_parameter(
"Time interval between each output", output_period);
390 maximum_non_linear_iterations);
394 prm.
add_parameter(
"Absolute error tolerance", absolute_tolerance);
395 prm.
add_parameter(
"Relative error tolerance", relative_tolerance);
397 "Ignore algebraic terms for error computations",
398 ignore_algebraic_terms_for_errors,
399 "Indicate whether or not to suppress algebraic variables "
400 "in the local error test.");
404 static std::string ic_type_str =
"use_y_diff";
406 "Correction type at initial time",
408 "This is one of the following three options for the "
409 "initial condition calculation. \n"
410 " none: do not try to make initial conditions consistent. \n"
411 " use_y_diff: compute the algebraic components of y and differential\n"
412 " components of y_dot, given the differential components of y. \n"
413 " This option requires that the user specifies differential and \n"
414 " algebraic components in the function get_differential_components.\n"
415 " use_y_dot: compute all components of y, given y_dot.",
417 prm.
add_action(
"Correction type at initial time",
418 [&](
const std::string &value) {
419 if (value ==
"use_y_diff")
420 ic_type = use_y_diff;
421 else if (value ==
"use_y_dot")
423 else if (value ==
"none")
429 static std::string reset_type_str =
"use_y_diff";
431 "Correction type after restart",
433 "This is one of the following three options for the "
434 "initial condition calculation. \n"
435 " none: do not try to make initial conditions consistent. \n"
436 " use_y_diff: compute the algebraic components of y and differential\n"
437 " components of y_dot, given the differential components of y. \n"
438 " This option requires that the user specifies differential and \n"
439 " algebraic components in the function get_differential_components.\n"
440 " use_y_dot: compute all components of y, given y_dot.",
442 prm.
add_action(
"Correction type after restart",
443 [&](
const std::string &value) {
444 if (value ==
"use_y_diff")
445 reset_type = use_y_diff;
446 else if (value ==
"use_y_dot")
447 reset_type = use_y_dot;
448 else if (value ==
"none")
454 maximum_non_linear_iterations_ic);
471 double initial_step_size;
476 double minimum_step_size;
481 double absolute_tolerance;
486 double relative_tolerance;
491 unsigned int maximum_order;
496 double output_period;
501 bool ignore_algebraic_terms_for_errors;
517 InitialConditionCorrection ic_type;
529 InitialConditionCorrection reset_type;
534 unsigned maximum_non_linear_iterations_ic;
539 unsigned int maximum_non_linear_iterations;
583 IDA(
const AdditionalData &data = AdditionalData(),
584 const MPI_Comm & mpi_comm = MPI_COMM_WORLD);
596 solve_dae(VectorType &solution, VectorType &solution_dot);
620 reset(
const double t,
const double h, VectorType &y, VectorType &yp);
625 std::function<void(VectorType &)> reinit_vector;
637 std::function<
int(
const double t,
639 const VectorType &y_dot,
673 std::function<
int(
const double t,
675 const VectorType &y_dot,
707 std::function<
int(
const VectorType &rhs, VectorType &dst)>
708 solve_jacobian_system;
724 std::function<void(
const double t,
725 const VectorType & sol,
726 const VectorType & sol_dot,
727 const unsigned int step_number)>
746 std::function<
bool(
const double t, VectorType &sol, VectorType &sol_dot)>
747 solver_should_restart;
763 std::function<
IndexSet()> differential_components;
771 std::function<VectorType &()> get_local_tolerances;
778 <<
"One of the SUNDIALS IDA internal functions "
779 <<
" returned a negative error code: " << arg1
780 <<
". Please consult SUNDIALS manual.");
790 <<
"Please provide an implementation for the function \""
799 set_functions_to_trigger_an_assert();
843# ifdef DEAL_II_WITH_PETSC
844# ifdef PETSC_USE_COMPLEX
845 static_assert(!std::is_same<VectorType, PETScWrappers::MPI::Vector>::value,
846 "Sundials does not support complex scalar types, "
847 "but PETSc is configured to use a complex scalar type!");
850 !std::is_same<VectorType, PETScWrappers::MPI::BlockVector>::value,
851 "Sundials does not support complex scalar types, "
852 "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)
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcInternalError()
#define AssertThrow(cond, exc)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
DeclException1(ExcARKodeError, int,<< "One of the SUNDIALS ARKode internal functions "<< " returned a negative error code: "<< arg1<< ". Please consult SUNDIALS manual.")