2.4.11.4.1. The MRIStepInnerStepper Class
As with other SUNDIALS classes, the MRIStepInnerStepper
abstract base
class is implemented using a C structure containing a content
pointer to the
derived class member data and a structure of function pointers (vtable) to the
derived class implementations of the base class virtual methods.
-
type MRIStepInnerStepper
An object for solving the fast (inner) ODE in an MRI method.
The actual definitions of the structure and the corresponding operations structure are kept private to allow for the object internals to change without impacting user code. The following sections describe the base (§2.4.11.4.1.1) and virtual methods (§2.4.11.4.1.2) that a must be provided by a derived class.
2.4.11.4.1.1. Base Class Methods
This section describes methods provided by the MRIStepInnerStepper
abstract base class that aid the user in implementing derived classes. This
includes functions for creating and destroying a generic base class object,
attaching and retrieving the derived class content
pointer, setting function
pointers to derived class method implementations, and accessing base class data
e.g., for computing the forcing term (2.77).
2.4.11.4.1.1.1. Creating and Destroying an Object
-
int MRIStepInnerStepper_Create(SUNContext sunctx, MRIStepInnerStepper *stepper)
This function creates an
MRIStepInnerStepper
object to which a user should attach the member data (content) pointer and method function pointers.- Parameters:
sunctx – the SUNDIALS simulation context.
stepper – a pointer to an inner stepper object.
- Return values:
ARK_SUCCESS – if successful
ARK_MEM_FAIL – if a memory allocation error occurs
Example usage:
/* create an instance of the base class */ MRIStepInnerStepper inner_stepper = NULL; flag = MRIStepInnerStepper_Create(&inner_stepper);
- Example codes:
examples/arkode/CXX_parallel/ark_diffusion_reaction_p.cpp
Note
See §2.4.11.4.1.1.2 and §2.4.11.4.1.1.3 for details on how to attach member data and method function pointers.
-
int MRIStepInnerStepper_CreateFromSUNStepper(SUNStepper sunstepper, MRIStepInnerStepper *stepper)
This utility function wraps a
SUNStepper
as anMRIStepInnerStepper
.- Parameters:
sunctx – the SUNDIALS simulation context.
sunstepper – the c:type:SUNStepper to wrap.
stepper – a pointer to an MRI inner stepper object.
- Return values:
ARK_SUCCESS – if successful
ARK_MEM_FAIL – if a memory allocation error occurs
Example usage:
SUNStepper sunstepper = NULL; SUNStepper_Create(ctx, &sunstepper); /* Attach content and functions to the SUNStepper... */ MRIStepInnerStepper inner_stepper = NULL; flag = MRIStepInnerStepper_CreateFromSUNStepper(sunstepper, &inner_stepper);
Added in version 6.2.0.
-
int MRIStepInnerStepper_Free(MRIStepInnerStepper *stepper)
This function destroys an
MRIStepInnerStepper
object.- Parameters:
stepper – a pointer to an inner stepper object.
- Return values:
ARK_SUCCESS – if successful
Example usage:
/* destroy an instance of the base class */ flag = MRIStepInnerStepper_Free(&inner_stepper);
- Example codes:
examples/arkode/CXX_parallel/ark_diffusion_reaction_p.cpp
Note
This function only frees memory allocated within the base class and the base class structure itself. The user is responsible for freeing any memory allocated for the member data (content).
2.4.11.4.1.1.2. Attaching and Accessing the Content Pointer
-
int MRIStepInnerStepper_SetContent(MRIStepInnerStepper stepper, void *content)
This function attaches a member data (content) pointer to an
MRIStepInnerStepper
object.- Parameters:
stepper – an inner stepper object.
content – a pointer to the stepper member data.
- Return values:
ARK_SUCCESS – if successful
ARK_ILL_INPUT – if the stepper is
NULL
Example usage:
/* set the inner stepper content pointer */ MyStepperContent my_object_data; flag = MRIStepInnerStepper_SetContent(inner_stepper, &my_object_data);
- Example codes:
examples/arkode/CXX_parallel/ark_diffusion_reaction_p.cpp
-
int MRIStepInnerStepper_GetContent(MRIStepInnerStepper stepper, void **content)
This function retrieves the member data (content) pointer from an
MRIStepInnerStepper
object.- Parameters:
stepper – an inner stepper object.
content – a pointer to set to the stepper member data pointer.
- Return values:
ARK_SUCCESS – if successful
ARK_ILL_INPUT – if the stepper is
NULL
Example usage:
/* get the inner stepper content pointer */ void *content; MyStepperContent *my_object_data; flag = MRIStepInnerStepper_GetContent(inner_stepper, &content); my_object_data = (MyStepperContent*) content;
- Example codes:
examples/arkode/CXX_parallel/ark_diffusion_reaction_p.cpp
2.4.11.4.1.1.3. Setting Member Functions
-
int MRIStepInnerStepper_SetEvolveFn(MRIStepInnerStepper stepper, MRIStepInnerEvolveFn fn)
This function attaches an
MRIStepInnerEvolveFn
function to anMRIStepInnerStepper
object.- Parameters:
stepper – an inner stepper object.
fn – the
MRIStepInnerStepper
function to attach.
- Return values:
ARK_SUCCESS – if successful
ARK_ILL_INPUT – if the stepper is
NULL
Example usage:
/* set the inner stepper evolve function */ flag = MRIStepInnerStepper_SetEvolveFn(inner_stepper, MyEvolve);
- Example codes:
examples/arkode/CXX_parallel/ark_diffusion_reaction_p.cpp
-
int MRIStepInnerStepper_SetFullRhsFn(MRIStepInnerStepper stepper, MRIStepInnerFullRhsFn fn)
This function attaches an
MRIStepInnerFullRhsFn
function to anMRIStepInnerStepper
object.- Parameters:
stepper – an inner stepper object.
fn – the
MRIStepInnerFullRhsFn
function to attach.
- Return values:
ARK_SUCCESS – if successful
ARK_ILL_INPUT – if the stepper is
NULL
Example usage:
/* set the inner stepper full right-hand side function */ flag = MRIStepInnerStepper_SetFullRhsFn(inner_stepper, MyFullRHS);
- Example codes:
examples/arkode/CXX_parallel/ark_diffusion_reaction_p.cpp
-
int MRIStepInnerStepper_SetResetFn(MRIStepInnerStepper stepper, MRIStepInnerResetFn fn)
This function attaches an
MRIStepInnerResetFn
function to anMRIStepInnerStepper
object.- Parameters:
stepper – an inner stepper object.
fn – the
MRIStepInnerResetFn
function to attach.
- Return values:
ARK_SUCCESS – if successful
ARK_ILL_INPUT – if the stepper is
NULL
Example usage:
/* set the inner stepper reset function */ flag = MRIStepInnerStepper_SetResetFn(inner_stepper, MyReset);
- Example codes:
examples/arkode/CXX_parallel/ark_diffusion_reaction_p.cpp
-
int MRIStepInnerStepper_SetAccumulatedErrorGetFn(MRIStepInnerStepper stepper, MRIStepInnerGetAccumulatedError fn)
This function attaches an
MRIStepInnerGetAccumulatedError
function to anMRIStepInnerStepper
object.- Parameters:
stepper – an inner stepper object.
fn – the
MRIStepInnerGetAccumulatedError
function to attach.
- Return values:
ARK_SUCCESS – if successful
ARK_ILL_INPUT – if the stepper is
NULL
-
int MRIStepInnerStepper_SetAccumulatedErrorResetFn(MRIStepInnerStepper stepper, MRIStepInnerResetAccumulatedError fn)
This function attaches an
MRIStepInnerResetAccumulatedError
function to anMRIStepInnerStepper
object.- Parameters:
stepper – an inner stepper object.
fn – the
MRIStepInnerResetAccumulatedError
function to attach.
- Return values:
ARK_SUCCESS – if successful
ARK_ILL_INPUT – if the stepper is
NULL
-
int MRIStepInnerStepper_SetRTolFn(MRIStepInnerStepper stepper, MRIStepInnerSetRTol fn)
This function attaches an
MRIStepInnerSetRTol
function to anMRIStepInnerStepper
object.- Parameters:
stepper – an inner stepper object.
fn – the
MRIStepInnerSetRTol
function to attach.
- Return values:
ARK_SUCCESS – if successful
ARK_ILL_INPUT – if the stepper is
NULL
2.4.11.4.1.1.4. Applying and Accessing Forcing Data
When integrating the ODE (2.76) the MRIStepInnerStepper
is
responsible for evaluating ODE right-hand side function \(f^F(t,v)\) as well
as computing and applying the forcing term (2.77) to obtain the
full right-hand side of the inner (fast) ODE (2.76). The functions in
this section can be used to either apply the inner (fast) forcing or access the
data necessary to construct the inner (fast) forcing polynomial. While the first of
these is less intrusive and may be used to package an existing black-box IVP solver
as an MRIStepInnerStepper, the latter may be more computationally efficient since it
does not traverse the data directly.
-
int MRIStepInnerStepper_AddForcing(MRIStepInnerStepper stepper, sunrealtype t, N_Vector ff)
This function computes the forcing term (2.77) at the input time t and adds it to input vector ff, i.e., the inner (fast) right-hand side vector.
- Parameters:
stepper – an inner stepper object.
t – the time at which the forcing should be evaluated.
f – the vector to which the forcing should be applied.
- Return values:
ARK_SUCCESS – if successful
ARK_ILL_INPUT – if the stepper is
NULL
- Example codes:
examples/arkode/CXX_parallel/ark_diffusion_reaction_p.cpp
-
int MRIStepInnerStepper_GetForcingData(MRIStepInnerStepper stepper, sunrealtype *tshift, sunrealtype *tscale, N_Vector **forcing, int *nforcing)
This function provides access to data necessary to compute the forcing term (2.77). This includes the shift and scaling factors for the normalized time \(\tau = (t - t_{n,i-1}^S)/(h^S \Delta c_i^S)\) and the array of polynomial coefficient vectors \(\hat{\gamma}^{i,k}\).
- Parameters:
stepper – an inner stepper object.
tshift – the time shift to apply to the current time when computing the forcing, \(t_{n,i-1}^S\).
tscale – the time scaling to apply to the current time when computing the forcing, \(h^S \Delta c_i^S\).
forcing – a pointer to an array of forcing vectors, \(\hat{\gamma}_{i,k}\).
nforcing – the number of forcing vectors.
- Return values:
ARK_SUCCESS – if successful
ARK_ILL_INPUT – if the stepper is
NULL
Example usage:
int k, flag; int nforcing_vecs; /* number of forcing vectors */ double tshift, tscale; /* time normalization values */ double tau; /* normalized time */ double tau_k; /* tau raised to the power k */ N_Vector *forcing_vecs; /* array of forcing vectors */ /* get the forcing data from the inner (fast) stepper */ flag = MRIStepInnerStepper_GetForcingData(inner_stepper, &tshift, &tscale, &forcing_vecs, &nforcing_vecs); /* compute the normalized time, initialize tau^k */ tau = (t - tshift) / tscale; tau_k = 1.0; /* compute the polynomial forcing terms and add them to fast RHS vector */ for (k = 0; k < nforcing_vecs; k++) { N_VLinearSum(1.0, f_fast, tau_k, forcing_vecs[k], f_fast); tau_k *= tau; }
- Example codes:
examples/arkode/CXX_parallel/ark_diffusion_reaction_p.cpp
2.4.11.4.1.2. Implementation Specific Methods
This section describes the required and optional virtual methods defined by the
MRIStepInnerStepper
abstract base class.
2.4.11.4.1.2.1. Required Member Functions
An MRIStepInnerStepper
must provide implementations of the following
member functions:
-
typedef int (*MRIStepInnerEvolveFn)(MRIStepInnerStepper stepper, sunrealtype t0, sunrealtype tout, N_Vector v)
This function advances the state vector v for the inner (fast) ODE system from time t0 to time tout.
- Arguments:
stepper – the inner stepper object.
t0 – the initial time for the inner (fast) integration.
tout – the final time for the inner (fast) integration.
v – on input the state at time t0 and, on output, the state at time tout.
- Return value:
An
MRIStepInnerEvolveFn
should return 0 if successful, a positive value if a recoverable error occurred (i.e., the function could be successful if called over a smaller time interval \([t0,tout]\)), or a negative value if it failed unrecoverably.- Example codes:
examples/arkode/CXX_parallel/ark_diffusion_reaction_p.cpp
2.4.11.4.1.2.2. Optional Member Functions
An MRIStepInnerStepper
may provide implementations of any of the
following member functions:
-
typedef int (*MRIStepInnerFullRhsFn)(MRIStepInnerStepper stepper, sunrealtype t, N_Vector v, N_Vector f, int mode)
This function computes the full right-hand side function of the inner (fast) ODE, \(f^F(t,v)\) in (2.76) for a given value of the independent variable t and state vector y. We note that this routine should not include contributions from the forcing term (2.77).
- Arguments:
stepper – the inner stepper object.
t – the current value of the independent variable.
y – the current value of the dependent variable vector.
f – the output vector that forms a portion the ODE right-hand side, \(f^F(t,y)\) in (2.15).
mode – a flag indicating the purpose for which the right-hand side function evaluation is called.
ARK_FULLRHS_START
– called at the beginning of the simulationARK_FULLRHS_END
– called at the end of a successful stepARK_FULLRHS_OTHER
– called elsewhere e.g., for dense output
- Return value:
An
MRIStepInnerFullRhsFn
should return 0 if successful, or a nonzero value upon failure.- Example codes:
examples/arkode/CXX_parallel/ark_diffusion_reaction_p.cpp
Changed in version v5.7.0: Supplying a full right-hand side function was made optional.
-
typedef int (*MRIStepInnerResetFn)(MRIStepInnerStepper stepper, sunrealtype tR, N_Vector vR)
This function resets the inner (fast) stepper state to the provided independent variable value and dependent variable vector.
If provided, the
MRIStepInnerResetFn
function will be called before a call toMRIStepInnerEvolveFn
when the state was updated at the slow timescale.- Arguments:
stepper – the inner stepper object.
tR – the value of the independent variable \(t_R\).
vR – the value of the dependent variable vector \(v(t_R)\).
- Return value:
An
MRIStepInnerResetFn
should return 0 if successful, or a nonzero value upon failure.- Example codes:
examples/arkode/CXX_parallel/ark_diffusion_reaction_p.cpp
-
typedef int (*MRIStepInnerGetAccumulatedError)(MRIStepInnerStepper stepper, sunrealtype *accum_error)
This function returns an estimate of the accumulated solution error arising from the inner stepper. Both the
MRIStepInnerGetAccumulatedError
andMRIStepInnerResetAccumulatedError
functions should be provided, or not; if only one is provided then MRIStep will disable multirate temporal adaptivity and call neither.- Arguments:
stepper – the inner stepper object.
accum_error – estimation of the accumulated solution error.
- Return value:
An
MRIStepInnerGetAccumulatedError
should return 0 if successful, a positive value if a recoverable error occurred (i.e., the function could be successful if called over a smaller time interval \([t0,tout]\)), or a negative value if it failed unrecoverably.
Note
This function is required when multirate temporal adaptivity has been enabled, using a
SUNAdaptController
module having typeSUN_ADAPTCONTROLLER_MRI_H_TOL
.If provided, the
MRIStepInnerGetAccumulatedError
function will always be called after a preceding call to theMRIStepInnerResetAccumulatedError
function.
-
typedef int (*MRIStepInnerResetAccumulatedError)(MRIStepInnerStepper stepper)
This function resets the inner stepper’s accumulated solution error to zero. This function performs a different role within MRIStep than the
MRIStepInnerResetFn
, and thus an implementation should make no assumptions about the frequency/ordering of calls to either.- Arguments:
stepper – the inner stepper object.
- Return value:
An
MRIStepInnerResetAccumulatedError
should return 0 if successful, or a nonzero value upon failure.
Note
This function is required when multirate temporal adaptivity has been enabled, using a
SUNAdaptController
module having typeSUN_ADAPTCONTROLLER_MRI_H_TOL
.The
MRIStepInnerResetAccumulatedError
function will always be called before any calls to theMRIStepInnerGetAccumulatedError
function.Both the
MRIStepInnerGetAccumulatedError
andMRIStepInnerResetAccumulatedError
functions should be provided, or not; if only one is provided then MRIStep will disable multirate temporal adaptivity and call neither.
-
typedef int (*MRIStepInnerSetRTol)(MRIStepInnerStepper stepper, sunrealtype rtol)
This function accepts a relative tolerance for the inner stepper to use in its upcoming adaptive solve. It is assumed that if the inner stepper supports absolute tolerances as well, then these have been set up directly by the user to indicate the “noise” level for solution components.
- Arguments:
stepper – the inner stepper object.
rtol – relative tolerance to use on the upcoming solve.
- Return value:
An
MRIStepInnerSetRTol
should return 0 if successful, or a nonzero value upon failure.
Note
This function is required when multirate temporal adaptivity has been enabled using a
SUNAdaptController
module having typeSUN_ADAPTCONTROLLER_MRI_H_TOL
.