1#ifndef AMREX_SUNDIALS_INTEGRATOR_H
2#define AMREX_SUNDIALS_INTEGRATOR_H
6#include <AMReX_Config.H>
14#include <nvector/nvector_manyvector.h>
15#include <sunnonlinsol/sunnonlinsol_fixedpoint.h>
16#include <sunlinsol/sunlinsol_spgmr.h>
17#include <arkode/arkode_arkstep.h>
18#include <arkode/arkode_mristep.h>
25 std::function<int(amrex::Real, N_Vector, N_Vector,
void*)>
f;
29 std::function<int(amrex::Real, N_Vector, N_Vector,
void*)>
fi;
30 std::function<int(amrex::Real, N_Vector, N_Vector,
void*)>
fe;
33 std::function<int(amrex::Real, N_Vector, N_Vector,
void*)>
ff;
36 std::function<int(amrex::Real, N_Vector,
void*)>
post_stage;
37 std::function<int(amrex::Real, N_Vector,
void*)>
post_step;
44namespace SundialsUserFun {
45 static int f (amrex::Real t, N_Vector y_data, N_Vector y_rhs,
void *user_data) {
47 return udata->
f(t, y_data, y_rhs, user_data);
50 static int fi (amrex::Real t, N_Vector y_data, N_Vector y_rhs,
void *user_data) {
52 return udata->
fi(t, y_data, y_rhs, user_data);
55 static int fe (amrex::Real t, N_Vector y_data, N_Vector y_rhs,
void *user_data) {
57 return udata->
fe(t, y_data, y_rhs, user_data);
60 static int ff (amrex::Real t, N_Vector y_data, N_Vector y_rhs,
void *user_data) {
62 return udata->
ff(t, y_data, y_rhs, user_data);
65 static int post_stage (amrex::Real t, N_Vector y_data,
void *user_data) {
67 return udata->
post_stage(t, y_data, user_data);
70 static int post_step (amrex::Real t, N_Vector y_data,
void *user_data) {
72 return udata->
post_step(t, y_data, user_data);
80 static int post_fast_step (amrex::Real t, N_Vector y_data,
void *user_data) {
133 SUNLinearSolver
LS =
nullptr;
134 SUNNonlinearSolver
NLS =
nullptr;
161 if (
type ==
"ERK" ||
type ==
"DIRK" ||
type ==
"IMEX-RK") {
164 else if (
type ==
"EX-MRI" ||
type ==
"IM-MRI" ||
type ==
"IMEX-MRI") {
168 std::string msg(
"Unknown method type: ");
190 void SetupRK (amrex::Real time, N_Vector y_data)
200 if (
method !=
"DEFAULT") {
205 else if (
type ==
"DIRK") {
209 if (
method !=
"DEFAULT") {
214 else if (
type ==
"IMEX-RK") {
241 if (
type ==
"DIRK" ||
type ==
"IMEX-RK") {
247 NLS = SUNNonlinSol_FixedPoint(y_data, 0,
sunctx);
363 if (
type ==
"EX-MRI") {
369 else if (
type ==
"IM-MRI") {
375 else if (
type ==
"IMEX-MRI") {
383 if (
method !=
"DEFAULT") {
384 MRIStepCoupling MRIC = MRIStepCoupling_LoadTableByName(
method.c_str());
388 MRIStepCoupling_Free(MRIC);
406 if (
type ==
"IM-MRI" ||
type ==
"IMEX-MRI") {
412 NLS = SUNNonlinSol_FixedPoint(y_data, 0,
sunctx);
457 const int num_vecs = N_VGetNumSubvectors_ManyVector(y_data);
458 S_data.resize(num_vecs);
460 for(
int i = 0; i < num_vecs; i++)
472 auto get_length = [&](
int index) -> sunindextype {
473 auto* p_mf = &S_data[index];
474 return p_mf->nComp() * (p_mf->boxArray()).numPts();
477 sunindextype NV_len = S_data.
size();
478 N_Vector* NV_array =
new N_Vector[NV_len];
480 for (
int i = 0; i < NV_len; ++i) {
485 N_Vector y_data = N_VNew_ManyVector(NV_len, NV_array,
sunctx);
495 auto get_length = [&](
int index) -> sunindextype {
496 auto* p_mf = &S_data[index];
497 return p_mf->nComp() * (p_mf->boxArray()).numPts();
500 sunindextype NV_len = S_data.
size();
501 N_Vector* NV_array =
new N_Vector[NV_len];
503 for (
int i = 0; i < NV_len; ++i) {
519 N_Vector y_data = N_VNew_ManyVector(NV_len, NV_array,
sunctx);
574 void initialize (
const T& S_data,
const amrex::Real time = 0.0)
578#if defined(SUNDIALS_VERSION_MAJOR) && (SUNDIALS_VERSION_MAJOR < 7)
580 sunctx = ::sundials::Context(&mpi_comm);
582 sunctx = ::sundials::Context(
nullptr);
586 sunctx = ::sundials::Context(mpi_comm);
588 sunctx = ::sundials::Context(SUN_COMM_NULL);
593 udata.
f = [&](amrex::Real rhs_time, N_Vector y_data, N_Vector y_rhs,
607 udata.
fi = [&](amrex::Real rhs_time, N_Vector y_data, N_Vector y_rhs,
621 udata.
fe = [&](amrex::Real rhs_time, N_Vector y_data, N_Vector y_rhs,
635 udata.
ff = [&](amrex::Real rhs_time, N_Vector y_data, N_Vector y_rhs,
709 if (
type ==
"EX-MRI" ||
type ==
"IM-MRI" ||
type ==
"IMEX-MRI") {
712 MRIStepPrintAllStats(
arkode_mem, stdout, SUN_OUTPUTFORMAT_TABLE);
721 ARKStepPrintAllStats(
arkode_mem, stdout, SUN_OUTPUTFORMAT_TABLE);
729 SUNNonlinSolFree(
NLS);
736 amrex::Real
advance (T& S_old, T& S_new, amrex::Real time,
const amrex::Real dt)
override
738 amrex::Real tout = time + dt;
747 int flag = ARKStepEvolve(
arkode_mem, tout, y_new, &tret, ARK_ONE_STEP);
753 int flag = MRIStepEvolve(
arkode_mem, tout, y_new, &tret, ARK_ONE_STEP);
756 Error(
"SUNDIALS integrator type not specified.");
765 void evolve (T& S_out,
const amrex::Real time_out)
override
768 amrex::Real time_ret;
776 flag = ARKStepEvolve(
arkode_mem, time_out, y_out, &time_ret, ARK_NORMAL);
786 flag = MRIStepEvolve(
arkode_mem, time_out, y_out, &time_ret, ARK_NORMAL);
789 Error(
"SUNDIALS integrator type not specified.");
797 void map_data (std::function<
void(T&)> )
override {}
#define AMREX_ALWAYS_ASSERT(EX)
Definition AMReX_BLassert.H:50
amrex::ParmParse pp
Input file parser instance for the given namespace.
Definition AMReX_HypreIJIface.cpp:15
Long numPts() const noexcept
Returns the total number of cells contained in all boxes in the BoxArray.
Definition AMReX_BoxArray.cpp:387
int nGrow(int direction=0) const noexcept
Return the grow factor that defines the region of definition.
Definition AMReX_FabArrayBase.H:78
const DistributionMapping & DistributionMap() const noexcept
Return constant reference to associated DistributionMapping.
Definition AMReX_FabArrayBase.H:131
int nComp() const noexcept
Return number of variables (aka components) associated with each point.
Definition AMReX_FabArrayBase.H:83
const BoxArray & boxArray() const noexcept
Return a constant reference to the BoxArray that defines the valid region associated with this FabArr...
Definition AMReX_FabArrayBase.H:95
Definition AMReX_IntegratorBase.H:164
bool use_adaptive_fast_time_step
Flag to enable/disable adaptive time stepping at the fast time scale in multirate methods (bool)
Definition AMReX_IntegratorBase.H:246
amrex::Real fast_rel_tol
Relative tolerance for adaptive time stepping at the fast time scale (Real)
Definition AMReX_IntegratorBase.H:278
amrex::Real rel_tol
Relative tolerance for adaptive time stepping (Real)
Definition AMReX_IntegratorBase.H:267
std::function< void(T &, amrex::Real)> post_fast_stage_action
The post_stage_action function is called by the integrator on the computed stage just after it is com...
Definition AMReX_IntegratorBase.H:218
amrex::Real fast_abs_tol
Absolute tolerance for adaptive time stepping at the fast time scale (Real)
Definition AMReX_IntegratorBase.H:284
amrex::Real fast_time_step
Current integrator fast time scale time step size with multirate methods (Real)
Definition AMReX_IntegratorBase.H:252
std::function< void(T &rhs, T &state, const amrex::Real time)> RhsEx
RhsEx is the explicit right-hand-side function an ImEx integrator will use.
Definition AMReX_IntegratorBase.H:194
std::function< void(T &, amrex::Real)> post_step_action
The post_step_action function is called by the integrator on the computed state just after it is comp...
Definition AMReX_IntegratorBase.H:212
std::function< void(T &, amrex::Real)> post_fast_step_action
The post_step_action function is called by the integrator on the computed state just after it is comp...
Definition AMReX_IntegratorBase.H:224
std::function< void(T &rhs, T &state, const amrex::Real time)> RhsIm
RhsIm is the implicit right-hand-side function an ImEx integrator will use.
Definition AMReX_IntegratorBase.H:188
std::function< void(T &rhs, T &state, const amrex::Real time)> Rhs
Rhs is the right-hand-side function the integrator will use.
Definition AMReX_IntegratorBase.H:182
bool use_adaptive_time_step
Flag to enable/disable adaptive time stepping in single rate methods or at the slow time scale in mul...
Definition AMReX_IntegratorBase.H:230
std::function< void(T &, amrex::Real)> post_stage_action
The post_stage_action function is called by the integrator on the computed stage just after it is com...
Definition AMReX_IntegratorBase.H:206
amrex::Real time_step
Current integrator time step size (Real)
Definition AMReX_IntegratorBase.H:235
std::function< void(T &rhs, T &state, const amrex::Real time)> RhsFast
RhsFast is the fast timescale right-hand-side function a multirate integrator will use.
Definition AMReX_IntegratorBase.H:200
amrex::Real abs_tol
Absolute tolerance for adaptive time stepping (Real)
Definition AMReX_IntegratorBase.H:272
A collection (stored as an array) of FArrayBox objects.
Definition AMReX_MultiFab.H:38
static void Copy(MultiFab &dst, const MultiFab &src, int srccomp, int dstcomp, int numcomp, int nghost)
Copy from src to dst including nghost ghost cells. The two MultiFabs MUST have the same underlying Bo...
Definition AMReX_MultiFab.cpp:193
Parse Parameters From Command Line and Input Files.
Definition AMReX_ParmParse.H:343
int query(std::string_view name, bool &ref, int ival=FIRST) const
Same as querykth() but searches for the last occurrence of name.
Definition AMReX_ParmParse.cpp:1393
This class provides the user with a few print options.
Definition AMReX_Print.H:35
Definition AMReX_SundialsIntegrator.H:88
SUNLinearSolver fast_LS
Definition AMReX_SundialsIntegrator.H:139
bool use_mri
Definition AMReX_SundialsIntegrator.H:120
void time_interpolate(const T &, const T &, amrex::Real, T &) override
Definition AMReX_SundialsIntegrator.H:795
SundialsIntegrator()
Definition AMReX_SundialsIntegrator.H:567
void unpack_vector(N_Vector y_data, amrex::MultiFab &S_data)
Definition AMReX_SundialsIntegrator.H:531
N_Vector copy_data(const amrex::MultiFab &S_data)
Definition AMReX_SundialsIntegrator.H:547
std::string linear_solver
Definition AMReX_SundialsIntegrator.H:112
std::string type
Definition AMReX_SundialsIntegrator.H:93
void initialize(const T &S_data, const amrex::Real time=0.0)
Definition AMReX_SundialsIntegrator.H:574
SUNNonlinearSolver NLS
Definition AMReX_SundialsIntegrator.H:134
::sundials::Context sunctx
Definition AMReX_SundialsIntegrator.H:129
int fast_max_nonlinear_iters
Definition AMReX_SundialsIntegrator.H:109
std::string fast_nonlinear_solver
Definition AMReX_SundialsIntegrator.H:108
amrex::Long max_num_steps
Definition AMReX_SundialsIntegrator.H:147
int max_nonlinear_iters
Definition AMReX_SundialsIntegrator.H:106
int max_linear_iters
Definition AMReX_SundialsIntegrator.H:113
std::string fast_linear_solver
Definition AMReX_SundialsIntegrator.H:115
MRIStepInnerStepper fast_stepper
Definition AMReX_SundialsIntegrator.H:138
void unpack_vector(N_Vector y_data, amrex::Vector< amrex::MultiFab > &S_data)
Definition AMReX_SundialsIntegrator.H:455
N_Vector copy_data(const amrex::Vector< amrex::MultiFab > &S_data)
Definition AMReX_SundialsIntegrator.H:493
void SetupRK(amrex::Real time, N_Vector y_data)
Definition AMReX_SundialsIntegrator.H:190
std::string method_e
Definition AMReX_SundialsIntegrator.H:97
void * arkode_mem
Definition AMReX_SundialsIntegrator.H:132
std::string method_i
Definition AMReX_SundialsIntegrator.H:98
std::string fast_method
Definition AMReX_SundialsIntegrator.H:102
amrex::Real advance(T &S_old, T &S_new, amrex::Real time, const amrex::Real dt) override
Take a single time step from (time, S_old) to (time + dt, S_new) with the given step size.
Definition AMReX_SundialsIntegrator.H:736
std::string nonlinear_solver
Definition AMReX_SundialsIntegrator.H:105
void * arkode_fast_mem
Definition AMReX_SundialsIntegrator.H:137
void evolve(T &S_out, const amrex::Real time_out) override
Evolve the current (internal) integrator state to time_out.
Definition AMReX_SundialsIntegrator.H:765
amrex::Real stop_time
Definition AMReX_SundialsIntegrator.H:144
N_Vector wrap_data(amrex::MultiFab &S_data)
Definition AMReX_SundialsIntegrator.H:540
bool use_ark
Definition AMReX_SundialsIntegrator.H:119
virtual ~SundialsIntegrator()
Definition AMReX_SundialsIntegrator.H:706
int fast_max_linear_iters
Definition AMReX_SundialsIntegrator.H:116
SundialsUserData udata
Definition AMReX_SundialsIntegrator.H:123
void initialize_parameters()
Definition AMReX_SundialsIntegrator.H:149
N_Vector wrap_data(amrex::Vector< amrex::MultiFab > &S_data)
Definition AMReX_SundialsIntegrator.H:470
SUNNonlinearSolver fast_NLS
Definition AMReX_SundialsIntegrator.H:140
std::string method
Definition AMReX_SundialsIntegrator.H:96
SUNLinearSolver LS
Definition AMReX_SundialsIntegrator.H:133
SundialsIntegrator(const T &S_data, const amrex::Real time=0.0)
Definition AMReX_SundialsIntegrator.H:569
void SetupMRI(amrex::Real time, N_Vector y_data)
Definition AMReX_SundialsIntegrator.H:285
std::string fast_type
Definition AMReX_SundialsIntegrator.H:101
void map_data(std::function< void(T &)>) override
Definition AMReX_SundialsIntegrator.H:797
bool set_stop_time
Definition AMReX_SundialsIntegrator.H:143
This class is a thin wrapper around std::vector. Unlike vector, Vector::operator[] provides bound che...
Definition AMReX_Vector.H:28
Long size() const noexcept
Definition AMReX_Vector.H:53
MPI_Comm CommunicatorSub() noexcept
sub-communicator for current frame
Definition AMReX_ParallelContext.H:70
bool IOProcessor() noexcept
Is this CPU the I/O Processor? To get the rank number, call IOProcessorNumber()
Definition AMReX_ParallelDescriptor.H:275
static int fi(amrex::Real t, N_Vector y_data, N_Vector y_rhs, void *user_data)
Definition AMReX_SundialsIntegrator.H:50
static int fe(amrex::Real t, N_Vector y_data, N_Vector y_rhs, void *user_data)
Definition AMReX_SundialsIntegrator.H:55
static int post_fast_step(amrex::Real t, N_Vector y_data, void *user_data)
Definition AMReX_SundialsIntegrator.H:80
static int f(amrex::Real t, N_Vector y_data, N_Vector y_rhs, void *user_data)
Definition AMReX_SundialsIntegrator.H:45
static int post_step(amrex::Real t, N_Vector y_data, void *user_data)
Definition AMReX_SundialsIntegrator.H:70
static int ff(amrex::Real t, N_Vector y_data, N_Vector y_rhs, void *user_data)
Definition AMReX_SundialsIntegrator.H:60
static int post_stage(amrex::Real t, N_Vector y_data, void *user_data)
Definition AMReX_SundialsIntegrator.H:65
static int post_fast_stage(amrex::Real t, N_Vector y_data, void *user_data)
Definition AMReX_SundialsIntegrator.H:75
int MPI_Comm
Definition AMReX_ccse-mpi.H:51
N_Vector N_VMake_MultiFab(sunindextype length, amrex::MultiFab *v_mf, ::sundials::Context *sunctx)
Definition AMReX_NVector_MultiFab.cpp:103
N_Vector N_VNew_MultiFab(sunindextype length, const amrex::BoxArray &ba, const amrex::DistributionMapping &dm, sunindextype nComp, sunindextype nGhost, ::sundials::Context *sunctx)
Definition AMReX_NVector_MultiFab.cpp:78
amrex::MultiFab *& getMFptr(N_Vector v)
Definition AMReX_NVector_MultiFab.cpp:228
Definition AMReX_Amr.cpp:49
@ make_alias
Definition AMReX_MakeType.H:7
int nComp(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2852
DistributionMapping const & DistributionMap(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2867
void Error(const std::string &msg)
Print out message to cerr and exit via amrex::Abort().
Definition AMReX.cpp:224
int Verbose() noexcept
Definition AMReX.cpp:169
BoxArray const & boxArray(FabArrayBase const &fa)
Definition AMReX_FabArrayBase.cpp:2862
Definition AMReX_SundialsIntegrator.H:22
std::function< int(amrex::Real, N_Vector, N_Vector, void *)> fi
Definition AMReX_SundialsIntegrator.H:29
std::function< int(amrex::Real, N_Vector, void *)> post_fast_stage
Definition AMReX_SundialsIntegrator.H:40
std::function< int(amrex::Real, N_Vector, void *)> post_step
Definition AMReX_SundialsIntegrator.H:37
std::function< int(amrex::Real, N_Vector, N_Vector, void *)> f
Definition AMReX_SundialsIntegrator.H:25
std::function< int(amrex::Real, N_Vector, void *)> post_stage
Definition AMReX_SundialsIntegrator.H:36
std::function< int(amrex::Real, N_Vector, N_Vector, void *)> ff
Definition AMReX_SundialsIntegrator.H:33
std::function< int(amrex::Real, N_Vector, void *)> post_fast_step
Definition AMReX_SundialsIntegrator.H:41
std::function< int(amrex::Real, N_Vector, N_Vector, void *)> fe
Definition AMReX_SundialsIntegrator.H:30