21 #ifdef DEAL_II_WITH_SUNDIALS 26 # ifdef DEAL_II_WITH_TRILINOS 30 # ifdef DEAL_II_WITH_PETSC 37 # include <arkode/arkode_impl.h> 38 # include <sundials/sundials_config.h> 54 template <
typename VectorType>
56 t_arkode_explicit_function(realtype tt,
61 ARKode<VectorType> &solver =
62 *
static_cast<ARKode<VectorType> *
>(user_data);
66 solver.reinit_vector(*src_yy);
69 solver.reinit_vector(*dst_yp);
73 int err = solver.explicit_function(tt, *src_yy, *dst_yp);
82 template <
typename VectorType>
84 t_arkode_implicit_function(realtype tt,
89 ARKode<VectorType> &solver =
90 *
static_cast<ARKode<VectorType> *
>(user_data);
94 solver.reinit_vector(*src_yy);
97 solver.reinit_vector(*dst_yp);
101 int err = solver.implicit_function(tt, *src_yy, *dst_yp);
110 template <
typename VectorType>
112 t_arkode_setup_jacobian(ARKodeMem arkode_mem,
116 booleantype *jcurPtr,
121 ARKode<VectorType> &solver =
122 *
static_cast<ARKode<VectorType> *
>(arkode_mem->ark_user_data);
126 solver.reinit_vector(*src_ypred);
129 solver.reinit_vector(*src_fpred);
131 copy(*src_ypred, ypred);
132 copy(*src_fpred, fpred);
135 bool jcurPtr_tmp =
false;
136 int err = solver.setup_jacobian(convfail,
138 arkode_mem->ark_gamma,
142 # if DEAL_II_SUNDIALS_VERSION_GTE(2, 0, 0) 143 *jcurPtr = jcurPtr_tmp ? SUNTRUE : SUNFALSE;
145 *jcurPtr = jcurPtr_tmp ? TRUE : FALSE;
153 template <
typename VectorType>
155 t_arkode_solve_jacobian(ARKodeMem arkode_mem,
163 ARKode<VectorType> &solver =
164 *
static_cast<ARKode<VectorType> *
>(arkode_mem->ark_user_data);
168 solver.reinit_vector(*src);
171 solver.reinit_vector(*src_ycur);
174 solver.reinit_vector(*src_fcur);
177 solver.reinit_vector(*dst);
180 copy(*src_ycur, ycur);
181 copy(*src_fcur, fcur);
183 int err = solver.solve_jacobian_system(arkode_mem->ark_tn,
184 arkode_mem->ark_gamma,
196 template <
typename VectorType>
198 t_arkode_setup_mass(ARKodeMem arkode_mem, N_Vector, N_Vector, N_Vector)
200 ARKode<VectorType> &solver =
201 *
static_cast<ARKode<VectorType> *
>(arkode_mem->ark_user_data);
202 int err = solver.setup_mass(arkode_mem->ark_tn);
208 template <
typename VectorType>
210 t_arkode_solve_mass(ARKodeMem arkode_mem,
219 ARKode<VectorType> &solver =
220 *
static_cast<ARKode<VectorType> *
>(arkode_mem->ark_user_data);
224 solver.reinit_vector(*src);
227 solver.reinit_vector(*dst);
231 int err = solver.solve_mass_system(*src, *dst);
238 template <
typename VectorType>
240 const MPI_Comm mpi_comm)
242 , arkode_mem(nullptr)
252 template <
typename VectorType>
257 # ifdef DEAL_II_WITH_MPI 269 template <
typename VectorType>
273 unsigned int system_size = solution.size();
277 unsigned int step_number = 0;
285 # ifdef DEAL_II_WITH_MPI 288 const IndexSet is = solution.locally_owned_elements();
289 const std::size_t local_system_size = is.
n_elements();
291 yy = N_VNew_Parallel(
communicator, local_system_size, system_size);
294 N_VNew_Parallel(
communicator, local_system_size, system_size);
301 "Trying to use a serial code with a parallel vector."));
302 yy = N_VNew_Serial(system_size);
326 reset(t, h, solution);
335 # ifdef DEAL_II_WITH_MPI 338 N_VDestroy_Parallel(
yy);
344 N_VDestroy_Serial(
yy);
351 template <
typename VectorType>
354 const double current_time_step,
357 unsigned int system_size;
367 # ifdef DEAL_II_WITH_MPI 370 N_VDestroy_Parallel(
yy);
376 N_VDestroy_Serial(
yy);
383 system_size = solution.size();
384 # ifdef DEAL_II_WITH_MPI 387 const IndexSet is = solution.locally_owned_elements();
388 const std::size_t local_system_size = is.
n_elements();
390 yy = N_VNew_Parallel(
communicator, local_system_size, system_size);
393 N_VNew_Parallel(
communicator, local_system_size, system_size);
398 yy = N_VNew_Serial(system_size);
405 ExcFunctionNotProvided(
"explicit_function || implicit_function"));
430 status = ARKodeSetInitStep(
arkode_mem, current_time_step);
444 auto ARKode_mem =
static_cast<ARKodeMem
>(
arkode_mem);
452 status = ARKodeSetLinear(
458 ARKode_mem->ark_lsolve = t_arkode_solve_jacobian<VectorType>;
461 ARKode_mem->ark_lsetup = t_arkode_setup_jacobian<VectorType>;
462 # if DEAL_II_SUNDIALS_VERSION_LT(3, 0, 0) 463 ARKode_mem->ark_setupNonNull =
true;
477 ARKode_mem->ark_msolve = t_arkode_solve_mass<VectorType>;
481 ARKode_mem->ark_msetup = t_arkode_setup_mass<VectorType>;
482 # if DEAL_II_SUNDIALS_VERSION_LT(3, 0, 0) 483 ARKode_mem->ark_MassSetupNonNull =
true;
492 template <
typename VectorType>
497 AssertThrow(
false, ExcFunctionNotProvided(
"reinit_vector"));
508 # ifdef DEAL_II_WITH_MPI 510 # ifdef DEAL_II_WITH_TRILINOS 513 # endif // DEAL_II_WITH_TRILINOS 515 # ifdef DEAL_II_WITH_PETSC 516 # ifndef PETSC_USE_COMPLEX 519 # endif // PETSC_USE_COMPLEX 520 # endif // DEAL_II_WITH_PETSC 522 # endif // DEAL_II_WITH_MPI void reset(const double t, const double h, const VectorType &y)
bool implicit_function_is_time_independent
#define AssertNothrow(cond, exc)
const auto & SundialsARKode
#define DEAL_II_SUNDIALS_VERSION_LT(major, minor, patch)
void copy(TrilinosWrappers::MPI::Vector &dst, const N_Vector &src)
std::function< int(const VectorType &rhs, VectorType &dst)> solve_mass_system
double relative_tolerance
#define AssertThrow(cond, exc)
std::function< void(const double t, const VectorType &sol, const unsigned int step_number)> output_step
std::function< VectorType &()> get_local_tolerances
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
std::function< int(const double t, const double gamma, const VectorType &ycur, const VectorType &fcur, const VectorType &rhs, VectorType &dst)> solve_jacobian_system
#define Assert(cond, exc)
#define DEAL_II_NAMESPACE_CLOSE
std::function< int(const double t, const VectorType &y, VectorType &res)> implicit_function
bool implicit_function_is_linear
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
unsigned int maximum_order
MPI_Comm duplicate_communicator(const MPI_Comm &mpi_communicator)
#define DEAL_II_NAMESPACE_OPEN
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
#define AssertARKode(code)
void set_functions_to_trigger_an_assert()
ARKode(const AdditionalData &data=AdditionalData(), const MPI_Comm mpi_comm=MPI_COMM_WORLD)
size_type n_elements() const
void copy(const T *begin, const T *end, U *dest)
static ::ExceptionBase & ExcInternalError()
double absolute_tolerance