21 #ifdef DEAL_II_WITH_SUNDIALS 26 # ifdef DEAL_II_WITH_TRILINOS 30 # ifdef DEAL_II_WITH_PETSC 37 # ifdef DEAL_II_SUNDIALS_WITH_IDAS 38 # include <idas/idas_impl.h> 40 # include <ida/ida_impl.h> 54 template <
typename VectorType>
56 t_dae_residual(realtype tt,
62 IDA<VectorType> &solver = *
static_cast<IDA<VectorType> *
>(user_data);
66 solver.reinit_vector(*src_yy);
69 solver.reinit_vector(*src_yp);
72 solver.reinit_vector(*residual);
77 int err = solver.residual(tt, *src_yy, *src_yp, *residual);
86 template <
typename VectorType>
88 t_dae_lsetup(IDAMem IDA_mem,
100 IDA<VectorType> &solver =
101 *
static_cast<IDA<VectorType> *
>(IDA_mem->ida_user_data);
105 solver.reinit_vector(*src_yy);
108 solver.reinit_vector(*src_yp);
113 int err = solver.setup_jacobian(IDA_mem->ida_tn,
122 template <
typename VectorType>
124 t_dae_solve(IDAMem IDA_mem,
135 IDA<VectorType> &solver =
136 *
static_cast<IDA<VectorType> *
>(IDA_mem->ida_user_data);
140 solver.reinit_vector(*src);
143 solver.reinit_vector(*dst);
147 int err = solver.solve_jacobian_system(*src, *dst);
155 template <
typename VectorType>
170 template <
typename VectorType>
175 # ifdef DEAL_II_WITH_MPI 187 template <
typename VectorType>
191 unsigned int system_size = solution.size();
195 unsigned int step_number = 0;
203 # ifdef DEAL_II_WITH_MPI 206 const IndexSet is = solution.locally_owned_elements();
207 const std::size_t local_system_size = is.
n_elements();
209 yy = N_VNew_Parallel(
communicator, local_system_size, system_size);
211 yp = N_VNew_Parallel(
communicator, local_system_size, system_size);
216 N_VNew_Parallel(
communicator, local_system_size, system_size);
223 "Trying to use a serial code with a parallel vector."));
224 yy = N_VNew_Serial(system_size);
225 yp = N_VNew_Serial(system_size);
226 diff_id = N_VNew_Serial(system_size);
239 status = IDASolve(
ida_mem, next_time, &t, yy, yp, IDA_NORMAL);
242 status = IDAGetLastStep(
ida_mem, &h);
246 copy(solution_dot, yp);
249 reset(t, h, solution, solution_dot);
253 output_step(t, solution, solution_dot, step_number);
257 # ifdef DEAL_II_WITH_MPI 260 N_VDestroy_Parallel(yy);
261 N_VDestroy_Parallel(yp);
268 N_VDestroy_Serial(yy);
269 N_VDestroy_Serial(yp);
277 template <
typename VectorType>
280 const double current_time_step,
284 unsigned int system_size;
296 # ifdef DEAL_II_WITH_MPI 299 N_VDestroy_Parallel(yy);
300 N_VDestroy_Parallel(yp);
307 N_VDestroy_Serial(yy);
308 N_VDestroy_Serial(yp);
316 system_size = solution.size();
317 # ifdef DEAL_II_WITH_MPI 320 const IndexSet is = solution.locally_owned_elements();
321 const std::size_t local_system_size = is.
n_elements();
323 yy = N_VNew_Parallel(
communicator, local_system_size, system_size);
325 yp = N_VNew_Parallel(
communicator, local_system_size, system_size);
330 N_VNew_Parallel(
communicator, local_system_size, system_size);
335 yy = N_VNew_Serial(system_size);
336 yp = N_VNew_Serial(system_size);
337 diff_id = N_VNew_Serial(system_size);
342 copy(yp, solution_dot);
344 status = IDAInit(
ida_mem, t_dae_residual<VectorType>, current_time, yy, yp);
355 status = IDASStolerances(
ida_mem,
361 status = IDASetInitStep(
ida_mem, current_time_step);
364 status = IDASetUserData(
ida_mem,
this);
372 diff_comp_vector = 0.0;
374 for (
auto i = dc.begin(); i != dc.end(); ++i)
375 diff_comp_vector[*i] = 1.0;
393 auto IDA_mem =
static_cast<IDAMem
>(
ida_mem);
395 IDA_mem->ida_lsetup = t_dae_lsetup<VectorType>;
396 IDA_mem->ida_lsolve = t_dae_solve<VectorType>;
397 # if DEAL_II_SUNDIALS_VERSION_LT(3, 0, 0) 398 IDA_mem->ida_setupNonNull =
true;
418 IDACalcIC(
ida_mem, IDA_Y_INIT, current_time + current_time_step);
421 status = IDAGetConsistentIC(
ida_mem, yy, yp);
425 copy(solution_dot, yp);
430 IDACalcIC(
ida_mem, IDA_YA_YDP_INIT, current_time + current_time_step);
433 status = IDAGetConsistentIC(
ida_mem, yy, yp);
437 copy(solution_dot, yp);
441 template <
typename VectorType>
446 AssertThrow(
false, ExcFunctionNotProvided(
"reinit_vector"));
452 VectorType &) ->
int {
454 AssertThrow(
false, ExcFunctionNotProvided(
"residual"));
463 AssertThrow(
false, ExcFunctionNotProvided(
"setup_jacobian"));
469 AssertThrow(
false, ExcFunctionNotProvided(
"solve_jacobian_system"));
476 const unsigned int) {
return; };
479 [](
const double, VectorType &, VectorType &) ->
bool {
return false; };
485 const unsigned int size = v->size();
493 # ifdef DEAL_II_WITH_MPI 495 # ifdef DEAL_II_WITH_TRILINOS 498 # endif // DEAL_II_WITH_TRILINOS 500 # ifdef DEAL_II_WITH_PETSC 501 # ifndef PETSC_USE_COMPLEX 504 # endif // PETSC_USE_COMPLEX 505 # endif // DEAL_II_WITH_PETSC 507 # endif // DEAL_II_WITH_MPI 513 #endif // DEAL_II_WITH_SUNDIALS void set_functions_to_trigger_an_assert()
GrowingVectorMemory< VectorType > mem
std::function< void(VectorType &)> reinit_vector
#define AssertNothrow(cond, exc)
std::function< IndexSet()> differential_components
unsigned int maximum_order
void copy(TrilinosWrappers::MPI::Vector &dst, const N_Vector &src)
InitialConditionCorrection
unsigned maximum_non_linear_iterations_ic
InitialConditionCorrection reset_type
void reset(const double t, const double h, VectorType &y, VectorType &yp)
#define AssertThrow(cond, exc)
IDA(const AdditionalData &data=AdditionalData(), const MPI_Comm mpi_comm=MPI_COMM_WORLD)
std::function< VectorType &()> get_local_tolerances
#define Assert(cond, exc)
bool ignore_algebraic_terms_for_errors
#define DEAL_II_NAMESPACE_CLOSE
std::function< int(const double t, const VectorType &y, const VectorType &y_dot, VectorType &res)> residual
std::function< void(const double t, const VectorType &sol, const VectorType &sol_dot, const unsigned int step_number)> output_step
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
std::function< int(const VectorType &rhs, VectorType &dst)> solve_jacobian_system
IndexSet complete_index_set(const IndexSet::size_type N)
std::function< bool(const double t, VectorType &sol, VectorType &sol_dot)> solver_should_restart
MPI_Comm duplicate_communicator(const MPI_Comm &mpi_communicator)
unsigned int maximum_non_linear_iterations
#define DEAL_II_NAMESPACE_OPEN
unsigned int solve_dae(VectorType &solution, VectorType &solution_dot)
double relative_tolerance
double absolute_tolerance
size_type n_elements() const
void copy(const T *begin, const T *end, U *dest)
InitialConditionCorrection ic_type
static ::ExceptionBase & ExcInternalError()
std::function< int(const double t, const VectorType &y, const VectorType &y_dot, const double alpha)> setup_jacobian