16 #ifndef dealii_sundials_arkode_h 17 #define dealii_sundials_arkode_h 19 #include <deal.II/base/config.h> 21 #include <deal.II/base/mpi.h> 23 #ifdef DEAL_II_WITH_SUNDIALS 25 # include <deal.II/base/conditional_ostream.h> 26 # include <deal.II/base/exceptions.h> 27 # include <deal.II/base/logstream.h> 28 # include <deal.II/base/mpi.h> 29 # include <deal.II/base/parameter_handler.h> 30 # ifdef DEAL_II_WITH_PETSC 31 # include <deal.II/lac/petsc_block_vector.h> 32 # include <deal.II/lac/petsc_vector.h> 34 # include <deal.II/lac/vector.h> 35 # include <deal.II/lac/vector_memory.h> 37 # include <arkode/arkode.h> 38 # include <arkode/arkode_impl.h> 39 # include <nvector/nvector_serial.h> 40 # ifdef DEAL_II_WITH_MPI 41 # include <nvector/nvector_parallel.h> 43 # include <boost/signals2.hpp> 45 # include <sundials/sundials_math.h> 46 # include <sundials/sundials_types.h> 51 DEAL_II_NAMESPACE_OPEN
54 # define AssertARKode(code) Assert(code >= 0, ExcARKodeError(code)) 313 template <
typename VectorType = Vector<
double>>
515 const MPI_Comm mpi_comm = MPI_COMM_WORLD);
544 reset(
const double t,
const double h,
const VectorType &y);
569 int(
const double t,
const VectorType &y, VectorType &explicit_f)>
588 std::function<int(const double t, const VectorType &y, VectorType &res)>
675 std::function<int(
const int convfail,
678 const VectorType &ypred,
679 const VectorType &fpred,
680 bool & j_is_current)>
727 std::function<int(
const double t,
729 const VectorType &ycur,
730 const VectorType &fcur,
731 const VectorType &rhs,
789 std::function<int(const VectorType &rhs, VectorType &dst)>
806 std::function<void(
const double t,
807 const VectorType & sol,
808 const unsigned int step_number)>
843 <<
"One of the SUNDIALS ARKode internal functions " 844 <<
" returned a negative error code: " << arg1
845 <<
". Please consult SUNDIALS manual.");
855 <<
"Please provide an implementation for the function \"" 898 # ifdef DEAL_II_WITH_PETSC 899 # ifdef PETSC_USE_COMPLEX 900 static_assert(!std::is_same<VectorType, PETScWrappers::MPI::Vector>::value,
901 "Sundials does not support complex scalar types, " 902 "but PETSc is configured to use a complex scalar type!");
905 !std::is_same<VectorType, PETScWrappers::MPI::BlockVector>::value,
906 "Sundials does not support complex scalar types, " 907 "but PETSc is configured to use a complex scalar type!");
908 # endif // PETSC_USE_COMPLEX 909 # endif // DEAL_II_WITH_PETSC 914 DEAL_II_NAMESPACE_CLOSE
void reset(const double t, const double h, const VectorType &y)
bool implicit_function_is_time_independent
GrowingVectorMemory< VectorType > mem
std::function< int(const VectorType &rhs, VectorType &dst)> solve_mass_system
double relative_tolerance
std::function< void(const double t, const VectorType &sol, const unsigned int step_number)> output_step
std::function< VectorType &()> get_local_tolerances
void add_parameter(const std::string &entry, ParameterType ¶meter, const std::string &documentation="", const Patterns::PatternBase &pattern= *Patterns::Tools::Convert< ParameterType >::to_pattern())
unsigned int maximum_non_linear_iterations
std::function< int(const double t)> setup_mass
std::function< bool(const double t, VectorType &sol)> solver_should_restart
void enter_subsection(const std::string &subsection)
std::function< int(const double t, const double gamma, const VectorType &ycur, const VectorType &fcur, const VectorType &rhs, VectorType &dst)> solve_jacobian_system
std::function< int(const double t, const VectorType &y, VectorType &res)> implicit_function
bool implicit_function_is_linear
DeclException1(ExcARKodeError, int,<< "One of the SUNDIALS ARKode internal functions "<< " returned a negative error code: "<< arg1<< ". Please consult SUNDIALS manual.")
unsigned int maximum_order
void add_parameters(ParameterHandler &prm)
std::function< int(const double t, const VectorType &y, VectorType &explicit_f)> explicit_function
std::function< void(VectorType &)> reinit_vector
unsigned int solve_ode(VectorType &solution)
std::function< int(const int convfail, const double t, const double gamma, const VectorType &ypred, const VectorType &fpred, bool &j_is_current)> setup_jacobian
void set_functions_to_trigger_an_assert()
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 double absolute_tolerance=1e-6, const double relative_tolerance=1e-5)
ARKode(const AdditionalData &data=AdditionalData(), const MPI_Comm mpi_comm=MPI_COMM_WORLD)
double absolute_tolerance