Reference documentation for deal.II version 9.4.1
\(\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\}}\)
Loading...
Searching...
No Matches
arkode.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2017 - 2022 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
17#ifndef dealii_sundials_arkode_h
18#define dealii_sundials_arkode_h
19
20#include <deal.II/base/config.h>
21
22#include <deal.II/base/mpi.h>
23
24#ifdef DEAL_II_WITH_SUNDIALS
25
30# ifdef DEAL_II_WITH_PETSC
33# endif
34# include <deal.II/lac/vector.h>
36
37# include <arkode/arkode.h>
38# include <nvector/nvector_serial.h>
39# ifdef DEAL_II_WITH_MPI
40# include <nvector/nvector_parallel.h>
41# endif
43
46
47# include <boost/signals2.hpp>
48
49# include <sundials/sundials_linearsolver.h>
50# include <sundials/sundials_math.h>
51# include <sundials/sundials_types.h>
52
53# include <memory>
54
55
57
58
59// Shorthand notation for ARKODE error codes.
60# define AssertARKode(code) Assert(code >= 0, ExcARKodeError(code))
61
65namespace SUNDIALS
66{
324 template <typename VectorType = Vector<double>>
325 class ARKode
326 {
327 public:
332 {
333 public:
365 // Initial parameters
366 const double initial_time = 0.0,
367 const double final_time = 1.0,
368 const double initial_step_size = 1e-2,
369 const double output_period = 1e-1,
370 // Running parameters
371 const double minimum_step_size = 1e-6,
372 const unsigned int maximum_order = 5,
373 const unsigned int maximum_non_linear_iterations = 10,
374 const bool implicit_function_is_linear = false,
375 const bool implicit_function_is_time_independent = false,
376 const bool mass_is_time_independent = false,
377 const int anderson_acceleration_subspace = 3,
378 // Error parameters
379 const double absolute_tolerance = 1e-6,
380 const double relative_tolerance = 1e-5);
381
396 void
398
403
408
413
418
423
428
432 unsigned int maximum_order;
433
439
445
450
456
462
468 };
469
480
488 ARKode(const AdditionalData &data, const MPI_Comm &mpi_comm);
489
493 ~ARKode();
494
502 unsigned int
503 solve_ode(VectorType &solution);
504
531 unsigned int
532 solve_ode_incrementally(VectorType & solution,
533 const double intermediate_time,
534 const bool reset_solver = false);
535
550 void
551 reset(const double t, const double h, const VectorType &y);
552
569 void *
570 get_arkode_memory() const;
571
581 std::function<void(VectorType &)> reinit_vector;
582
599 std::function<
600 int(const double t, const VectorType &y, VectorType &explicit_f)>
602
619 std::function<int(const double t, const VectorType &y, VectorType &res)>
621
639 std::function<int(double t, const VectorType &v, VectorType &Mv)>
641
676 std::function<int(const double t)> mass_times_setup;
677
704 std::function<int(const VectorType &v,
705 VectorType & Jv,
706 double t,
707 const VectorType &y,
708 const VectorType &fy)>
710
745 std::function<int(realtype t, const VectorType &y, const VectorType &fy)>
747
768
785
817 std::function<int(double t,
818 const VectorType &y,
819 const VectorType &fy,
820 const VectorType &r,
821 VectorType & z,
822 double gamma,
823 double tol,
824 int lr)>
826
869 std::function<int(double t,
870 const VectorType &y,
871 const VectorType &fy,
872 int jok,
873 int & jcur,
874 double gamma)>
876
904 std::function<
905 int(double t, const VectorType &r, VectorType &z, double tol, int lr)>
907
932 std::function<int(double t)> mass_preconditioner_setup;
933
948 std::function<void(const double t,
949 const VectorType & sol,
950 const unsigned int step_number)>
952
970 std::function<bool(const double t, VectorType &sol)> solver_should_restart;
971
978 std::function<VectorType &()> get_local_tolerances;
979
1004 std::function<void(void *arkode_mem)> custom_setup;
1005
1006 private:
1011 DeclException1(ExcFunctionNotProvided,
1012 std::string,
1013 << "Please provide an implementation for the function \""
1014 << arg1 << "\"");
1015
1019 int
1020 do_evolve_time(VectorType & solution,
1021 ::DiscreteTime &time,
1022 const bool do_reset);
1023
1030 void
1031 setup_system_solver(const VectorType &solution);
1032
1039 void
1040 setup_mass_solver(const VectorType &solution);
1041
1047 void
1049
1054
1059
1060# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
1064 SUNContext arkode_ctx;
1065# endif
1066
1072
1077
1078 std::unique_ptr<internal::LinearSolverWrapper<VectorType>> linear_solver;
1079 std::unique_ptr<internal::LinearSolverWrapper<VectorType>> mass_solver;
1080
1081# ifdef DEAL_II_WITH_PETSC
1082# ifdef PETSC_USE_COMPLEX
1083 static_assert(!std::is_same<VectorType, PETScWrappers::MPI::Vector>::value,
1084 "Sundials does not support complex scalar types, "
1085 "but PETSc is configured to use a complex scalar type!");
1086
1087 static_assert(
1088 !std::is_same<VectorType, PETScWrappers::MPI::BlockVector>::value,
1089 "Sundials does not support complex scalar types, "
1090 "but PETSc is configured to use a complex scalar type!");
1091# endif // PETSC_USE_COMPLEX
1092# endif // DEAL_II_WITH_PETSC
1093 };
1094
1095
1099 DeclException1(ExcARKodeError,
1100 int,
1101 << "One of the SUNDIALS ARKode internal functions "
1102 << " returned a negative error code: " << arg1
1103 << ". Please consult SUNDIALS manual.");
1104
1105
1106 template <typename VectorType>
1108 // Initial parameters
1109 const double initial_time,
1110 const double final_time,
1111 const double initial_step_size,
1112 const double output_period,
1113 // Running parameters
1114 const double minimum_step_size,
1115 const unsigned int maximum_order,
1116 const unsigned int maximum_non_linear_iterations,
1117 const bool implicit_function_is_linear,
1118 const bool implicit_function_is_time_independent,
1119 const bool mass_is_time_independent,
1120 const int anderson_acceleration_subspace,
1121 // Error parameters
1122 const double absolute_tolerance,
1123 const double relative_tolerance)
1124 : initial_time(initial_time)
1125 , final_time(final_time)
1126 , initial_step_size(initial_step_size)
1127 , minimum_step_size(minimum_step_size)
1128 , absolute_tolerance(absolute_tolerance)
1129 , relative_tolerance(relative_tolerance)
1130 , maximum_order(maximum_order)
1131 , output_period(output_period)
1132 , maximum_non_linear_iterations(maximum_non_linear_iterations)
1133 , implicit_function_is_linear(implicit_function_is_linear)
1134 , implicit_function_is_time_independent(
1135 implicit_function_is_time_independent)
1136 , mass_is_time_independent(mass_is_time_independent)
1137 , anderson_acceleration_subspace(anderson_acceleration_subspace)
1138 {}
1139
1140
1141
1142 template <typename VectorType>
1143 void
1145 {
1146 prm.add_parameter("Initial time", initial_time);
1147 prm.add_parameter("Final time", final_time);
1148 prm.add_parameter("Time interval between each output", output_period);
1149 prm.enter_subsection("Running parameters");
1150 prm.add_parameter("Initial step size", initial_step_size);
1151 prm.add_parameter("Minimum step size", minimum_step_size);
1152 prm.add_parameter("Maximum order of ARK", maximum_order);
1153 prm.add_parameter("Maximum number of nonlinear iterations",
1154 maximum_non_linear_iterations);
1155 prm.add_parameter("Implicit function is linear",
1156 implicit_function_is_linear);
1157 prm.add_parameter("Implicit function is time independent",
1158 implicit_function_is_time_independent);
1159 prm.add_parameter("Mass is time independent", mass_is_time_independent);
1160 prm.add_parameter("Anderson-acceleration subspace",
1161 anderson_acceleration_subspace);
1162 prm.leave_subsection();
1163 prm.enter_subsection("Error control");
1164 prm.add_parameter("Absolute error tolerance", absolute_tolerance);
1165 prm.add_parameter("Relative error tolerance", relative_tolerance);
1166 prm.leave_subsection();
1167 }
1168
1169} // namespace SUNDIALS
1170
1172#endif
1173
1174
1175#endif
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)
void enter_subsection(const std::string &subsection)
void add_parameters(ParameterHandler &prm)
Definition: arkode.h:1144
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)
Definition: arkode.h:1107
unsigned int maximum_non_linear_iterations
Definition: arkode.h:444
void setup_mass_solver(const VectorType &solution)
Definition: arkode.cc:606
std::function< int(const double t, const VectorType &y, VectorType &res)> implicit_function
Definition: arkode.h:620
void * arkode_mem
Definition: arkode.h:1058
std::function< int(double t, const VectorType &r, VectorType &z, double tol, int lr)> mass_preconditioner_solve
Definition: arkode.h:906
std::function< int(const VectorType &v, VectorType &Jv, double t, const VectorType &y, const VectorType &fy)> jacobian_times_vector
Definition: arkode.h:709
std::function< int(double t)> mass_preconditioner_setup
Definition: arkode.h:932
void * get_arkode_memory() const
Definition: arkode.cc:695
std::function< VectorType &()> get_local_tolerances
Definition: arkode.h:978
AdditionalData data
Definition: arkode.h:1053
std::function< int(realtype t, const VectorType &y, const VectorType &fy)> jacobian_times_setup
Definition: arkode.h:746
LinearSolveFunction< VectorType > solve_linearized_system
Definition: arkode.h:767
int do_evolve_time(VectorType &solution, ::DiscreteTime &time, const bool do_reset)
Definition: arkode.cc:337
MPI_Comm mpi_communicator
Definition: arkode.h:1071
std::function< int(const double t)> mass_times_setup
Definition: arkode.h:676
DeclException1(ExcFunctionNotProvided, std::string,<< "Please provide an implementation for the function \""<< arg1<< "\"")
void reset(const double t, const double h, const VectorType &y)
Definition: arkode.cc:405
SUNContext arkode_ctx
Definition: arkode.h:1064
std::unique_ptr< internal::LinearSolverWrapper< VectorType > > mass_solver
Definition: arkode.h:1079
std::function< int(double t, const VectorType &y, const VectorType &fy, int jok, int &jcur, double gamma)> jacobian_preconditioner_setup
Definition: arkode.h:875
std::function< int(const double t, const VectorType &y, VectorType &explicit_f)> explicit_function
Definition: arkode.h:601
std::function< bool(const double t, VectorType &sol)> solver_should_restart
Definition: arkode.h:970
std::function< int(double t, const VectorType &v, VectorType &Mv)> mass_times_vector
Definition: arkode.h:640
std::function< void(const double t, const VectorType &sol, const unsigned int step_number)> output_step
Definition: arkode.h:951
void set_functions_to_trigger_an_assert()
Definition: arkode.cc:680
std::function< int(double t, const VectorType &y, const VectorType &fy, const VectorType &r, VectorType &z, double gamma, double tol, int lr)> jacobian_preconditioner_solve
Definition: arkode.h:825
std::unique_ptr< internal::LinearSolverWrapper< VectorType > > linear_solver
Definition: arkode.h:1078
LinearSolveFunction< VectorType > solve_mass
Definition: arkode.h:784
unsigned int solve_ode(VectorType &solution)
Definition: arkode.cc:305
void setup_system_solver(const VectorType &solution)
Definition: arkode.cc:499
unsigned int solve_ode_incrementally(VectorType &solution, const double intermediate_time, const bool reset_solver=false)
Definition: arkode.cc:318
double last_end_time
Definition: arkode.h:1076
std::function< void(void *arkode_mem)> custom_setup
Definition: arkode.h:1004
DEAL_II_DEPRECATED std::function< void(VectorType &)> reinit_vector
Definition: arkode.h:581
#define DEAL_II_DEPRECATED
Definition: config.h:164
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:509
std::function< int(SundialsOperator< VectorType > &op, SundialsPreconditioner< VectorType > &prec, VectorType &x, const VectorType &b, double tol)> LinearSolveFunction