deal.II version GIT relicensing-2289-g1e5549a87a 2024-12-21 21:30:00+00:00
\(\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.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2017 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15
16
17#include <deal.II/base/config.h>
18
20
21#ifdef DEAL_II_WITH_SUNDIALS
22
25
29# ifdef DEAL_II_WITH_TRILINOS
32# endif
33# ifdef DEAL_II_WITH_PETSC
36# endif
37
41
42# include <arkode/arkode_arkstep.h>
43# include <sunlinsol/sunlinsol_spgmr.h>
44# include <sunnonlinsol/sunnonlinsol_fixedpoint.h>
45
46# include <iostream>
47
49
50namespace SUNDIALS
51{
52 template <typename VectorType>
54 : ARKode(data, MPI_COMM_SELF)
55 {}
56
57
58 template <typename VectorType>
60 const MPI_Comm mpi_comm)
61 : data(data)
62 , arkode_mem(nullptr)
64 , arkode_ctx(nullptr)
65# endif
66 , mpi_communicator(mpi_comm)
67 , last_end_time(data.initial_time)
68 , pending_exception(nullptr)
69 {
71
72 // SUNDIALS will always duplicate communicators if we provide them. This
73 // can cause problems if SUNDIALS is configured with MPI and we pass along
74 // MPI_COMM_SELF in a serial application as MPI won't be
75 // initialized. Hence, work around that by just not providing a
76 // communicator in that case.
77# if DEAL_II_SUNDIALS_VERSION_GTE(7, 0, 0)
78 const int status =
79 SUNContext_Create(mpi_communicator == MPI_COMM_SELF ? SUN_COMM_NULL :
81 &arkode_ctx);
82 (void)status;
83 AssertARKode(status);
84# elif DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
85 const int status =
86 SUNContext_Create(mpi_communicator == MPI_COMM_SELF ? nullptr :
88 &arkode_ctx);
89 (void)status;
90 AssertARKode(status);
91# endif
92 }
93
94
95
96 template <typename VectorType>
98 {
99 ARKStepFree(&arkode_mem);
100
101# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
102 const int status = SUNContext_Free(&arkode_ctx);
103 (void)status;
104 AssertARKode(status);
105# endif
106
107 Assert(pending_exception == nullptr, ExcInternalError());
108 }
109
110
111
112 template <typename VectorType>
113 unsigned int
114 ARKode<VectorType>::solve_ode(VectorType &solution)
115 {
116 DiscreteTime time(data.initial_time,
117 data.final_time,
118 data.initial_step_size);
119
120 return do_evolve_time(solution, time, /* force_solver_restart = */ true);
121 }
122
123
124
125 template <typename VectorType>
126 unsigned int
128 const double intermediate_time,
129 const bool reset_solver)
130 {
132 intermediate_time > last_end_time,
134 "The requested intermediate time is smaller than the last requested "
135 "intermediate time."));
136
137 const bool do_reset = reset_solver || arkode_mem == nullptr;
138 DiscreteTime time(last_end_time, intermediate_time, data.initial_step_size);
139 return do_evolve_time(solution, time, do_reset);
140 }
141
142
143
144 template <typename VectorType>
145 unsigned int
147 DiscreteTime &time,
148 const bool do_reset)
149 {
150 if (do_reset)
151 {
152 reset(time.get_current_time(), time.get_next_step_size(), solution);
153 if (output_step)
154 output_step(time.get_current_time(),
155 solution,
156 time.get_step_number());
157 }
158 else
159 {
160 // If we don't do a full reset then we still need to fix the end time.
161 // In SUNDIALS 6 and later, SUNDIALS will not do timesteps if the
162 // current time is past the set end point (i.e., ARKStepEvolve will
163 // return ARK_TSTOP_RETURN).
164 const int status = ARKStepSetStopTime(arkode_mem, time.get_end_time());
165 (void)status;
166 AssertARKode(status);
167 }
168
169 auto solution_nvector = internal::make_nvector_view(solution
171 ,
172 arkode_ctx
173# endif
174 );
175
176 while (!time.is_at_end())
177 {
178 time.set_desired_next_step_size(data.output_period);
179
180 // Having set up all of the ancillary things, finally call the main
181 // ARKode function. Once we return, check what happened:
182 // - If we have a pending recoverable exception, ignore it if SUNDIAL's
183 // return code was zero -- in that case, SUNDIALS managed to indeed
184 // recover and we no longer need the exception
185 // - If we have any other exception, rethrow it
186 // - If no exception, test that SUNDIALS really did successfully return
187 Assert(pending_exception == nullptr, ExcInternalError());
188 double actual_next_time;
189 const auto status = ARKStepEvolve(arkode_mem,
190 time.get_next_time(),
191 solution_nvector,
192 &actual_next_time,
193 ARK_NORMAL);
194 if (pending_exception)
195 {
196 try
197 {
198 std::rethrow_exception(pending_exception);
199 }
200 catch (const RecoverableUserCallbackError &exc)
201 {
202 pending_exception = nullptr;
203 if (status == 0)
204 /* just eat the exception */;
205 else
206 throw;
207 }
208 catch (...)
209 {
210 pending_exception = nullptr;
211 throw;
212 }
213 }
214
215 AssertARKode(status);
216
217 // Then reflect this time advancement in our own DiscreteTime object:
218 time.set_next_step_size(actual_next_time - time.get_current_time());
219 time.advance_time();
220
221 // Finally check whether resets or output calls are desired at this
222 // time:
223 while (solver_should_restart(time.get_current_time(), solution))
224 reset(time.get_current_time(),
226 solution);
227
228 if (output_step)
229 output_step(time.get_current_time(),
230 solution,
231 time.get_step_number());
232 }
233 last_end_time = time.get_current_time();
234
235 long int n_steps;
236 const auto status = ARKStepGetNumSteps(arkode_mem, &n_steps);
237 (void)status;
238 AssertARKode(status);
239
240 return n_steps;
241 }
242
243
244
245 template <typename VectorType>
246 void
247 ARKode<VectorType>::reset(const double current_time,
248 const double current_time_step,
249 const VectorType &solution)
250 {
251 last_end_time = current_time;
252 int status;
253 (void)status;
254
255# if DEAL_II_SUNDIALS_VERSION_GTE(7, 0, 0)
256 status = SUNContext_Free(&arkode_ctx);
257 AssertARKode(status);
258
259 // Same comment applies as in class constructor:
260 status =
261 SUNContext_Create(mpi_communicator == MPI_COMM_SELF ? SUN_COMM_NULL :
262 mpi_communicator,
263 &arkode_ctx);
264 AssertARKode(status);
265# elif DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
266 status = SUNContext_Free(&arkode_ctx);
267 AssertARKode(status);
268
269 // Same comment applies as in class constructor:
270 status =
271 SUNContext_Create(mpi_communicator == MPI_COMM_SELF ? nullptr :
272 &mpi_communicator,
273 &arkode_ctx);
274 AssertARKode(status);
275# endif
276
277 if (arkode_mem)
278 {
279 ARKStepFree(&arkode_mem);
280 // Initialization is version-dependent: do that in a moment
281 }
282
283 // just a view on the memory in solution, all write operations on yy by
284 // ARKODE will automatically be mirrored to solution
285 auto initial_condition_nvector = internal::make_nvector_view(solution
287 ,
288 arkode_ctx
289# endif
290 );
291
292 Assert(explicit_function || implicit_function,
293 ExcFunctionNotProvided("explicit_function || implicit_function"));
294
295 auto explicit_function_callback = [](SUNDIALS::realtype tt,
296 N_Vector yy,
297 N_Vector yp,
298 void *user_data) -> int {
299 Assert(user_data != nullptr, ExcInternalError());
300 ARKode<VectorType> &solver =
301 *static_cast<ARKode<VectorType> *>(user_data);
302
303 auto *src_yy = internal::unwrap_nvector_const<VectorType>(yy);
304 auto *dst_yp = internal::unwrap_nvector<VectorType>(yp);
305
307 solver.explicit_function,
308 solver.pending_exception,
309 tt,
310 *src_yy,
311 *dst_yp);
312 };
313
314
315 auto implicit_function_callback = [](SUNDIALS::realtype tt,
316 N_Vector yy,
317 N_Vector yp,
318 void *user_data) -> int {
319 Assert(user_data != nullptr, ExcInternalError());
320 ARKode<VectorType> &solver =
321 *static_cast<ARKode<VectorType> *>(user_data);
322
323 auto *src_yy = internal::unwrap_nvector_const<VectorType>(yy);
324 auto *dst_yp = internal::unwrap_nvector<VectorType>(yp);
325
327 solver.implicit_function,
328 solver.pending_exception,
329 tt,
330 *src_yy,
331 *dst_yp);
332 };
333
334 arkode_mem = ARKStepCreate(explicit_function ? explicit_function_callback :
335 ARKRhsFn(nullptr),
336 implicit_function ? implicit_function_callback :
337 ARKRhsFn(nullptr),
338 current_time,
339 initial_condition_nvector
340# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
341 ,
342 arkode_ctx
343# endif
344 );
345 Assert(arkode_mem != nullptr, ExcInternalError());
346
347 if (get_local_tolerances)
348 {
349 const auto abs_tols = internal::make_nvector_view(get_local_tolerances()
351 ,
352 arkode_ctx
353# endif
354 );
355 status =
356 ARKStepSVtolerances(arkode_mem, data.relative_tolerance, abs_tols);
357 AssertARKode(status);
358 }
359 else
360 {
361 status = ARKStepSStolerances(arkode_mem,
362 data.relative_tolerance,
363 data.absolute_tolerance);
364 AssertARKode(status);
365 }
366
367 // for version 4.0.0 this call must be made before solver settings
368 status = ARKStepSetUserData(arkode_mem, this);
369 AssertARKode(status);
370
371 setup_system_solver(solution);
372
373 setup_mass_solver(solution);
374
375 status = ARKStepSetInitStep(arkode_mem, current_time_step);
376 AssertARKode(status);
377
378 status = ARKStepSetStopTime(arkode_mem, data.final_time);
379 AssertARKode(status);
380
381 status = ARKStepSetOrder(arkode_mem, data.maximum_order);
382 AssertARKode(status);
383
384 if (custom_setup)
385 custom_setup(arkode_mem);
386 }
387
388
389
390 template <typename VectorType>
391 void
392 ARKode<VectorType>::setup_system_solver(const VectorType &solution)
393 {
394 // no point in setting up a solver when there is no linear system
395 if (!implicit_function)
396 return;
397
398 int status;
399 (void)status;
400 // Initialize solver
401 // currently only iterative linear solvers are supported
402 if (jacobian_times_vector)
403 {
404 SUNLinearSolver sun_linear_solver;
405 if (solve_linearized_system)
406 {
407 linear_solver =
408 std::make_unique<internal::LinearSolverWrapper<VectorType>>(
409 solve_linearized_system,
410 pending_exception
411# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
412 ,
413 arkode_ctx
414# endif
415 );
416 sun_linear_solver = *linear_solver;
417 }
418 else
419 {
420 // use default solver from SUNDIALS
421 // TODO give user options
422 auto y_template = internal::make_nvector_view(solution
424 ,
425 arkode_ctx
426# endif
427 );
428# if DEAL_II_SUNDIALS_VERSION_LT(6, 0, 0)
429 sun_linear_solver =
430 SUNLinSol_SPGMR(y_template,
431 PREC_NONE,
432 0 /*krylov subvectors, 0 uses default*/);
433# else
434 sun_linear_solver =
435 SUNLinSol_SPGMR(y_template,
436 SUN_PREC_NONE,
437 0 /*krylov subvectors, 0 uses default*/,
438 arkode_ctx);
439# endif
440 }
441 status = ARKStepSetLinearSolver(arkode_mem, sun_linear_solver, nullptr);
442 AssertARKode(status);
443
444 auto jacobian_times_vector_callback = [](N_Vector v,
445 N_Vector Jv,
447 N_Vector y,
448 N_Vector fy,
449 void *user_data,
450 N_Vector) -> int {
451 Assert(user_data != nullptr, ExcInternalError());
452 ARKode<VectorType> &solver =
453 *static_cast<ARKode<VectorType> *>(user_data);
454
455 auto *src_v = internal::unwrap_nvector_const<VectorType>(v);
456 auto *src_y = internal::unwrap_nvector_const<VectorType>(y);
457 auto *src_fy = internal::unwrap_nvector_const<VectorType>(fy);
458
459 auto *dst_Jv = internal::unwrap_nvector<VectorType>(Jv);
460
463 solver.pending_exception,
464 *src_v,
465 *dst_Jv,
466 t,
467 *src_y,
468 *src_fy);
469 };
470
471 auto jacobian_times_vector_setup_callback = [](SUNDIALS::realtype t,
472 N_Vector y,
473 N_Vector fy,
474 void *user_data) -> int {
475 Assert(user_data != nullptr, ExcInternalError());
476 ARKode<VectorType> &solver =
477 *static_cast<ARKode<VectorType> *>(user_data);
478
479 auto *src_y = internal::unwrap_nvector_const<VectorType>(y);
480 auto *src_fy = internal::unwrap_nvector_const<VectorType>(fy);
481
484 solver.pending_exception,
485 t,
486 *src_y,
487 *src_fy);
488 };
489 status = ARKStepSetJacTimes(arkode_mem,
490 jacobian_times_setup ?
491 jacobian_times_vector_setup_callback :
492 ARKLsJacTimesSetupFn(nullptr),
493 jacobian_times_vector_callback);
494 AssertARKode(status);
495 if (jacobian_preconditioner_solve)
496 {
497 auto solve_with_jacobian_callback = [](SUNDIALS::realtype t,
498 N_Vector y,
499 N_Vector fy,
500 N_Vector r,
501 N_Vector z,
502 SUNDIALS::realtype gamma,
503 SUNDIALS::realtype delta,
504 int lr,
505 void *user_data) -> int {
506 Assert(user_data != nullptr, ExcInternalError());
507 ARKode<VectorType> &solver =
508 *static_cast<ARKode<VectorType> *>(user_data);
509
510 auto *src_y = internal::unwrap_nvector_const<VectorType>(y);
511 auto *src_fy = internal::unwrap_nvector_const<VectorType>(fy);
512 auto *src_r = internal::unwrap_nvector_const<VectorType>(r);
513
514 auto *dst_z = internal::unwrap_nvector<VectorType>(z);
515
518 solver.pending_exception,
519 t,
520 *src_y,
521 *src_fy,
522 *src_r,
523 *dst_z,
524 gamma,
525 delta,
526 lr);
527 };
528
529 auto jacobian_solver_setup_callback =
531 N_Vector y,
532 N_Vector fy,
534 SUNDIALS::booltype *jcurPtr,
535 SUNDIALS::realtype gamma,
536 void *user_data) -> int {
537 Assert(user_data != nullptr, ExcInternalError());
538 ARKode<VectorType> &solver =
539 *static_cast<ARKode<VectorType> *>(user_data);
540
541 auto *src_y = internal::unwrap_nvector_const<VectorType>(y);
542 auto *src_fy = internal::unwrap_nvector_const<VectorType>(fy);
543
546 solver.pending_exception,
547 t,
548 *src_y,
549 *src_fy,
550 jok,
551 *jcurPtr,
552 gamma);
553 };
554
555 status = ARKStepSetPreconditioner(arkode_mem,
556 jacobian_preconditioner_setup ?
557 jacobian_solver_setup_callback :
558 ARKLsPrecSetupFn(nullptr),
559 solve_with_jacobian_callback);
560 AssertARKode(status);
561 }
562 if (data.implicit_function_is_linear)
563 {
564 status = ARKStepSetLinear(
565 arkode_mem, data.implicit_function_is_time_independent ? 0 : 1);
566 AssertARKode(status);
567 }
568 }
569 else
570 {
571 auto y_template = internal::make_nvector_view(solution
573 ,
574 arkode_ctx
575# endif
576 );
577
578# if DEAL_II_SUNDIALS_VERSION_LT(6, 0, 0)
579 SUNNonlinearSolver fixed_point_solver =
580 SUNNonlinSol_FixedPoint(y_template,
581 data.anderson_acceleration_subspace);
582# else
583 SUNNonlinearSolver fixed_point_solver =
584 SUNNonlinSol_FixedPoint(y_template,
585 data.anderson_acceleration_subspace,
586 arkode_ctx);
587# endif
588
589 status = ARKStepSetNonlinearSolver(arkode_mem, fixed_point_solver);
590 AssertARKode(status);
591 }
592
593 status =
594 ARKStepSetMaxNonlinIters(arkode_mem, data.maximum_non_linear_iterations);
595 AssertARKode(status);
596 }
597
598
599
600 template <typename VectorType>
601 void
602 ARKode<VectorType>::setup_mass_solver(const VectorType &solution)
603 {
604 int status;
605 (void)status;
606
607 if (mass_times_vector)
608 {
609 SUNLinearSolver sun_mass_linear_solver;
610 if (solve_mass)
611 {
612 mass_solver =
613 std::make_unique<internal::LinearSolverWrapper<VectorType>>(
614 solve_mass,
615
616 pending_exception
617# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
618 ,
619 arkode_ctx
620# endif
621 );
622 sun_mass_linear_solver = *mass_solver;
623 }
624 else
625 {
626 auto y_template = internal::make_nvector_view(solution
628 ,
629 arkode_ctx
630# endif
631 );
632# if DEAL_II_SUNDIALS_VERSION_LT(6, 0, 0)
633 sun_mass_linear_solver =
634 SUNLinSol_SPGMR(y_template,
635 PREC_NONE,
636 0 /*krylov subvectors, 0 uses default*/);
637# else
638 sun_mass_linear_solver =
639 SUNLinSol_SPGMR(y_template,
640 SUN_PREC_NONE,
641 0 /*krylov subvectors, 0 uses default*/,
642 arkode_ctx);
643# endif
644 }
645
646 SUNDIALS::booltype mass_time_dependent =
647 data.mass_is_time_independent ? SUNFALSE : SUNTRUE;
648
649 status = ARKStepSetMassLinearSolver(arkode_mem,
650 sun_mass_linear_solver,
651 nullptr,
652 mass_time_dependent);
653 AssertARKode(status);
654
655 auto mass_matrix_times_vector_setup_callback =
656 [](SUNDIALS::realtype t, void *mtimes_data) -> int {
657 Assert(mtimes_data != nullptr, ExcInternalError());
658 ARKode<VectorType> &solver =
659 *static_cast<ARKode<VectorType> *>(mtimes_data);
660
662 solver.mass_times_setup, solver.pending_exception, t);
663 };
664
665 auto mass_matrix_times_vector_callback = [](N_Vector v,
666 N_Vector Mv,
668 void *mtimes_data) -> int {
669 Assert(mtimes_data != nullptr, ExcInternalError());
670 ARKode<VectorType> &solver =
671 *static_cast<ARKode<VectorType> *>(mtimes_data);
672
673 auto *src_v = internal::unwrap_nvector_const<VectorType>(v);
674 auto *dst_Mv = internal::unwrap_nvector<VectorType>(Mv);
675
677 solver.mass_times_vector,
678 solver.pending_exception,
679 t,
680 *src_v,
681 *dst_Mv);
682 };
683
684 status = ARKStepSetMassTimes(arkode_mem,
685 mass_times_setup ?
686 mass_matrix_times_vector_setup_callback :
687 ARKLsMassTimesSetupFn(nullptr),
688 mass_matrix_times_vector_callback,
689 this);
690 AssertARKode(status);
691
692 if (mass_preconditioner_solve)
693 {
694 auto mass_matrix_solver_setup_callback =
695 [](SUNDIALS::realtype t, void *user_data) -> int {
696 Assert(user_data != nullptr, ExcInternalError());
697 ARKode<VectorType> &solver =
698 *static_cast<ARKode<VectorType> *>(user_data);
699
702 };
703
704 auto solve_with_mass_matrix_callback = [](SUNDIALS::realtype t,
705 N_Vector r,
706 N_Vector z,
707 SUNDIALS::realtype delta,
708 int lr,
709 void *user_data) -> int {
710 Assert(user_data != nullptr, ExcInternalError());
711 ARKode<VectorType> &solver =
712 *static_cast<ARKode<VectorType> *>(user_data);
713
714 auto *src_r = internal::unwrap_nvector_const<VectorType>(r);
715 auto *dst_z = internal::unwrap_nvector<VectorType>(z);
716
719 solver.pending_exception,
720 t,
721 *src_r,
722 *dst_z,
723 delta,
724 lr);
725 };
726
727 status =
728 ARKStepSetMassPreconditioner(arkode_mem,
729 mass_preconditioner_setup ?
730 mass_matrix_solver_setup_callback :
731 ARKLsMassPrecSetupFn(nullptr),
732 solve_with_mass_matrix_callback);
733 AssertARKode(status);
734 }
735 }
736 }
737
738
739
740 template <typename VectorType>
741 void
743 {
744 solver_should_restart = [](const double, VectorType &) -> bool {
745 return false;
746 };
747 }
748
749
750
751 template <typename VectorType>
752 void *
754 {
755 return arkode_mem;
756 }
757
758 template class ARKode<Vector<double>>;
759 template class ARKode<BlockVector<double>>;
760
763
764# ifdef DEAL_II_WITH_MPI
765
766# ifdef DEAL_II_WITH_TRILINOS
769
770# endif // DEAL_II_WITH_TRILINOS
771
772# ifdef DEAL_II_WITH_PETSC
773# ifndef PETSC_USE_COMPLEX
776# endif // PETSC_USE_COMPLEX
777# endif // DEAL_II_WITH_PETSC
778
779# endif // DEAL_II_WITH_MPI
780
781} // namespace SUNDIALS
782
784
785#endif
#define AssertARKode(code)
Definition arkode.h:60
double get_next_step_size() const
void set_desired_next_step_size(const double time_step_size)
unsigned int get_step_number() const
bool is_at_end() const
double get_current_time() const
void set_next_step_size(const double time_step_size)
double get_next_time() const
double get_previous_step_size() const
void advance_time()
double get_end_time() const
void setup_mass_solver(const VectorType &solution)
Definition arkode.cc:602
std::function< void(const double t)> mass_preconditioner_setup
Definition arkode.h:942
std::function< void(const double t, const VectorType &y, VectorType &explicit_f)> explicit_function
Definition arkode.h:593
std::exception_ptr pending_exception
Definition arkode.h:1096
std::function< void(const double t, const VectorType &y, const VectorType &fy, const VectorType &r, VectorType &z, const double gamma, const double tol, const int lr)> jacobian_preconditioner_solve
Definition arkode.h:832
ARKode(const AdditionalData &data=AdditionalData())
Definition arkode.cc:53
void * get_arkode_memory() const
Definition arkode.cc:753
std::function< void(const double t, const VectorType &y, VectorType &res)> implicit_function
Definition arkode.h:612
std::function< void(const double t, const VectorType &r, VectorType &z, const double tol, const int lr)> mass_preconditioner_solve
Definition arkode.h:916
MPI_Comm mpi_communicator
Definition arkode.h:1081
void reset(const double t, const double h, const VectorType &y)
Definition arkode.cc:247
SUNContext arkode_ctx
Definition arkode.h:1074
std::function< void(const double t, const VectorType &y, const VectorType &fy, const int jok, int &jcur, const double gamma)> jacobian_preconditioner_setup
Definition arkode.h:882
unsigned int do_evolve_time(VectorType &solution, ::DiscreteTime &time, const bool do_reset)
Definition arkode.cc:146
std::function< void(const double t)> mass_times_setup
Definition arkode.h:668
void set_functions_to_trigger_an_assert()
Definition arkode.cc:742
std::function< void(const double t, const VectorType &v, VectorType &Mv)> mass_times_vector
Definition arkode.h:632
std::function< void(const double t, const VectorType &y, const VectorType &fy)> jacobian_times_setup
Definition arkode.h:739
unsigned int solve_ode(VectorType &solution)
Definition arkode.cc:114
void setup_system_solver(const VectorType &solution)
Definition arkode.cc:392
unsigned int solve_ode_incrementally(VectorType &solution, const double intermediate_time, const bool reset_solver=false)
Definition arkode.cc:127
std::function< void(const VectorType &v, VectorType &Jv, const double t, const VectorType &y, const VectorType &fy)> jacobian_times_vector
Definition arkode.h:701
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_SUNDIALS_VERSION_GTE(major, minor, patch)
Definition config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
static ::ExceptionBase & ExcFunctionNotProvided(std::string arg1)
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & RecoverableUserCallbackError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
std::vector< index_type > data
Definition mpi.cc:735
int call_and_possibly_capture_exception(const F &f, std::exception_ptr &eptr, Args &&...args)
Definition utilities.h:45
sunbooleantype booltype
::sunrealtype realtype