|
Reference documentation for deal.II version 9.2.0
|
\(\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\}}\)
Go to the documentation of this file.
16 #ifndef dealii_sundials_arkode_h
17 #define dealii_sundials_arkode_h
23 #ifdef DEAL_II_WITH_SUNDIALS
29 # ifdef DEAL_II_WITH_PETSC
36 # include <arkode/arkode.h>
37 # include <arkode/arkode_impl.h>
38 # include <nvector/nvector_serial.h>
39 # ifdef DEAL_II_WITH_MPI
40 # include <nvector/nvector_parallel.h>
42 # include <boost/signals2.hpp>
44 # include <sundials/sundials_math.h>
45 # include <sundials/sundials_types.h>
53 # define AssertARKode(code) Assert(code >= 0, ExcARKodeError(code))
312 template <
typename VectorType = Vector<
double>>
514 const MPI_Comm mpi_comm = MPI_COMM_WORLD);
674 std::function<
int(
const int convfail,
679 bool & j_is_current)>
726 std::function<
int(
const double t,
805 std::function<void(
const double t,
807 const unsigned int step_number)>
842 <<
"One of the SUNDIALS ARKode internal functions "
843 <<
" returned a negative error code: " << arg1
844 <<
". Please consult SUNDIALS manual.");
854 <<
"Please provide an implementation for the function \""
897 # ifdef DEAL_II_WITH_PETSC
898 # ifdef PETSC_USE_COMPLEX
900 "Sundials does not support complex scalar types, "
901 "but PETSc is configured to use a complex scalar type!");
905 "Sundials does not support complex scalar types, "
906 "but PETSc is configured to use a complex scalar type!");
907 # endif // PETSC_USE_COMPLEX
908 # endif // DEAL_II_WITH_PETSC
long double gamma(const unsigned int n)
unsigned int maximum_non_linear_iterations
void reset(const double t, const double h, const VectorType &y)
double absolute_tolerance
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
unsigned int maximum_order
std::function< int(const int convfail, const double t, const double gamma, const VectorType &ypred, const VectorType &fpred, bool &j_is_current)> setup_jacobian
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)
unsigned int solve_ode(VectorType &solution)
double relative_tolerance
std::function< void(const double t, const VectorType &sol, const unsigned int step_number)> output_step
#define DEAL_II_NAMESPACE_OPEN
DeclException1(ExcARKodeError, int,<< "One of the SUNDIALS ARKode internal functions "<< " returned a negative error code: "<< arg1<< ". Please consult SUNDIALS manual.")
std::function< VectorType &()> get_local_tolerances
void add_parameters(ParameterHandler &prm)
bool implicit_function_is_time_independent
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 &explicit_f)> explicit_function
std::function< void(VectorType &)> reinit_vector
GrowingVectorMemory< VectorType > mem
bool implicit_function_is_linear
ARKode(const AdditionalData &data=AdditionalData(), const MPI_Comm mpi_comm=MPI_COMM_WORLD)
void set_functions_to_trigger_an_assert()
std::function< int(const VectorType &rhs, VectorType &dst)> solve_mass_system
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)
std::function< int(const double t)> setup_mass
std::function< bool(const double t, VectorType &sol)> solver_should_restart
#define DEAL_II_NAMESPACE_CLOSE
std::function< int(const double t, const VectorType &y, VectorType &res)> implicit_function