3.4.5.2. MRIStep User-callable functions

This section describes the functions that are called by the user to setup and then solve an IVP using the MRIStep time-stepping module. Some of these are required; however, starting with §3.4.5.2.7, the functions listed involve optional inputs/outputs or restarting, and those paragraphs may be skipped for a casual use of ARKODE’s MRIStep module. In any case, refer to the preceding section, §3.4.5.1, for the correct order of these calls.

On an error, each user-callable function returns a negative value (or NULL if the function returns a pointer) and sends an error message to the error handler routine, which prints the message to stderr by default. However, the user can set a file as error output or can provide their own error handler function (see §3.4.5.2.7 for details).

3.4.5.2.1. MRIStep initialization and deallocation functions

void *MRIStepCreate(ARKRhsFn fse, ARKRhsFn fsi, sunrealtype t0, N_Vector y0, MRIStepInnerStepper stepper, SUNContext sunctx)

This function allocates and initializes memory for a problem to be solved using the MRIStep time-stepping module in ARKODE.

Arguments:
  • fse – the name of the function (of type ARKRhsFn()) defining the explicit slow portion of the right-hand side function in \(\dot{y} = f^E(t,y) + f^I(t,y) + f^F(t,y)\).

  • fsi – the name of the function (of type ARKRhsFn()) defining the implicit slow portion of the right-hand side function in \(\dot{y} = f^E(t,y) + f^I(t,y) + f^F(t,y)\).

  • t0 – the initial value of \(t\).

  • y0 – the initial condition vector \(y(t_0)\).

  • stepper – an MRIStepInnerStepper for integrating the fast time scale.

  • sunctx – the SUNContext object (see §2.4)

Return value:

If successful, a pointer to initialized problem memory of type void*, to be passed to all user-facing MRIStep routines listed below. If unsuccessful, a NULL pointer will be returned, and an error message will be printed to stderr.

Example usage:

/* fast (inner) and slow (outer) ARKODE objects */
void *inner_arkode_mem = NULL;
void *outer_arkode_mem = NULL;

/* MRIStepInnerStepper to wrap the inner (fast) ARKStep object */
MRIStepInnerStepper stepper = NULL;

/* create an ARKStep object, setting fast (inner) right-hand side
   functions and the initial condition */
inner_arkode_mem = ARKStepCreate(ffe, ffi, t0, y0, sunctx);

/* setup ARKStep */
. . .

/* create MRIStepInnerStepper wrapper for the ARKStep memory block */
flag = ARKStepCreateMRIStepInnerStepper(inner_arkode_mem, &stepper);

/* create an MRIStep object, setting the slow (outer) right-hand side
   functions and the initial condition */
outer_arkode_mem = MRIStepCreate(fse, fsi, t0, y0, stepper, sunctx)
Example codes:
  • examples/arkode/C_serial/ark_brusselator_mri.c

  • examples/arkode/C_serial/ark_twowaycouple_mri.c

  • examples/arkode/C_serial/ark_brusselator_1D_mri.c

  • examples/arkode/C_serial/ark_onewaycouple_mri.c

  • examples/arkode/C_serial/ark_reaction_diffusion_mri.c

  • examples/arkode/C_serial/ark_kpr_mri.c

  • examples/arkode/CXX_parallel/ark_diffusion_reaction_p.cpp

void MRIStepFree(void **arkode_mem)

This function frees the problem memory arkode_mem created by MRIStepCreate().

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

Return value: None

3.4.5.2.2. MRIStep tolerance specification functions

These functions specify the integration tolerances. One of them should be called before the first call to MRIStepEvolve(); otherwise default values of reltol = 1e-4 and abstol = 1e-9 will be used, which may be entirely incorrect for a specific problem.

The integration tolerances reltol and abstol define a vector of error weights, ewt. In the case of MRIStepSStolerances(), this vector has components

ewt[i] = 1.0/(reltol*abs(y[i]) + abstol);

whereas in the case of MRIStepSVtolerances() the vector components are given by

ewt[i] = 1.0/(reltol*abs(y[i]) + abstol[i]);

This vector is used in all error tests, which use a weighted RMS norm on all error-like vectors \(v\):

\[\|v\|_{WRMS} = \left( \frac{1}{N} \sum_{i=1}^N (v_i\; ewt_i)^2 \right)^{1/2},\]

where \(N\) is the problem dimension.

Alternatively, the user may supply a custom function to supply the ewt vector, through a call to MRIStepWFtolerances().

int MRIStepSStolerances(void *arkode_mem, sunrealtype reltol, sunrealtype abstol)

This function specifies scalar relative and absolute tolerances.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • reltol – scalar relative tolerance.

  • abstol – scalar absolute tolerance.

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

  • ARK_NO_MALLOC if the MRIStep memory was not allocated by the time-stepping module

  • ARK_ILL_INPUT if an argument has an illegal value (e.g. a negative tolerance).

int MRIStepSVtolerances(void *arkode_mem, sunrealtype reltol, N_Vector abstol)

This function specifies a scalar relative tolerance and a vector absolute tolerance (a potentially different absolute tolerance for each vector component).

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • reltol – scalar relative tolerance.

  • abstol – vector containing the absolute tolerances for each solution component.

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

  • ARK_NO_MALLOC if the MRIStep memory was not allocated by the time-stepping module

  • ARK_ILL_INPUT if an argument has an illegal value (e.g. a negative tolerance).

int MRIStepWFtolerances(void *arkode_mem, ARKEwtFn efun)

This function specifies a user-supplied function efun to compute the error weight vector ewt.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • efun – the name of the function (of type ARKEwtFn()) that implements the error weight vector computation.

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

  • ARK_NO_MALLOC if the MRIStep memory was not allocated by the time-stepping module

3.4.5.2.2.1. General advice on the choice of tolerances

For many users, the appropriate choices for tolerance values in reltol and abstol are a concern. The following pieces of advice are relevant.

  1. The scalar relative tolerance reltol is to be set to control relative errors. So a value of \(10^{-4}\) means that errors are controlled to .01%. We do not recommend using reltol larger than \(10^{-3}\). On the other hand, reltol should not be so small that it is comparable to the unit roundoff of the machine arithmetic (generally around \(10^{-15}\) for double-precision).

  2. The absolute tolerances abstol (whether scalar or vector) need to be set to control absolute errors when any components of the solution vector \(y\) may be so small that pure relative error control is meaningless. For example, if \(y_i\) starts at some nonzero value, but in time decays to zero, then pure relative error control on \(y_i\) makes no sense (and is overly costly) after \(y_i\) is below some noise level. Then abstol (if scalar) or abstol[i] (if a vector) needs to be set to that noise level. If the different components have different noise levels, then abstol should be a vector. For example, see the example problem ark_robertson.c, and the discussion of it in the ARKODE Examples Documentation [99]. In that problem, the three components vary between 0 and 1, and have different noise levels; hence the atols vector therein. It is impossible to give any general advice on abstol values, because the appropriate noise levels are completely problem-dependent. The user or modeler hopefully has some idea as to what those noise levels are.

  3. Finally, it is important to pick all the tolerance values conservatively, because they control the error committed on each individual step. The final (global) errors are an accumulation of those per-step errors, where that accumulation factor is problem-dependent. A general rule of thumb is to reduce the tolerances by a factor of 10 from the actual desired limits on errors. So if you want .01% relative accuracy (globally), a good choice for reltol is \(10^{-5}\). In any case, it is a good idea to do a few experiments with the tolerances to see how the computed solution values vary as tolerances are reduced.

3.4.5.2.2.2. Advice on controlling nonphysical negative values

In many applications, some components in the true solution are always positive or non-negative, though at times very small. In the numerical solution, however, small negative (nonphysical) values can then occur. In most cases, these values are harmless, and simply need to be controlled, not eliminated, but in other cases any value that violates a constraint may cause a simulation to halt. For both of these scenarios the following pieces of advice are relevant.

  1. The best way to control the size of unwanted negative computed values is with tighter absolute tolerances. Again this requires some knowledge of the noise level of these components, which may or may not be different for different components. Some experimentation may be needed.

  2. If output plots or tables are being generated, and it is important to avoid having negative numbers appear there (for the sake of avoiding a long explanation of them, if nothing else), then eliminate them, but only in the context of the output medium. Then the internal values carried by the solver are unaffected. Remember that a small negative value in \(y\) returned by MRIStep, with magnitude comparable to abstol or less, is equivalent to zero as far as the computation is concerned.

  3. The user’s right-hand side routine \(f^I\) should never change a negative value in the solution vector \(y\) to a non-negative value in attempt to “fix” this problem, since this can lead to numerical instability. If the \(f^I\) routine cannot tolerate a zero or negative value (e.g. because there is a square root or log), then the offending value should be changed to zero or a tiny positive number in a temporary variable (not in the input \(y\) vector) for the purposes of computing \(f^I(t, y)\).

3.4.5.2.3. Linear solver interface functions

As previously explained, the Newton iterations used in solving implicit systems within MRIStep require the solution of linear systems of the form

\[\mathcal{A}\left(z_i^{(m)}\right) \delta^{(m+1)} = -G\left(z_i^{(m)}\right)\]

where

\[\mathcal{A} \approx I - \gamma J, \qquad J = \frac{\partial f^I}{\partial y}.\]

ARKODE’s ARKLS linear solver interface supports all valid SUNLinearSolver modules for this task.

Matrix-based SUNLinearSolver modules utilize SUNMatrix objects to store the approximate Jacobian matrix \(J\), the Newton matrix \(\mathcal{A}\), and, when using direct solvers, the factorizations used throughout the solution process.

Matrix-free SUNLinearSolver modules instead use iterative methods to solve the Newton systems of equations, and only require the action of the matrix on a vector, \(\mathcal{A}v\). With most of these methods, preconditioning can be done on the left only, on the right only, on both the left and the right, or not at all. The exceptions to this rule are SPFGMR that supports right preconditioning only and PCG that performs symmetric preconditioning. For the specification of a preconditioner, see the iterative linear solver portions of §3.4.5.2.7 and §3.4.6.

If preconditioning is done, user-supplied functions should be used to define left and right preconditioner matrices \(P_1\) and \(P_2\) (either of which could be the identity matrix), such that the product \(P_{1}P_{2}\) approximates the Newton matrix \(\mathcal{A} = I - \gamma J\).

To specify a generic linear solver for MRIStep to use for the Newton systems, after the call to MRIStepCreate() but before any calls to MRIStepEvolve(), the user’s program must create the appropriate SUNLinearSolver object and call the function MRIStepSetLinearSolver(), as documented below. To create the SUNLinearSolver object, the user may call one of the SUNDIALS-packaged SUNLinSol module constructor routines via a call of the form

SUNLinearSolver LS = SUNLinSol_*(...);

The current list of SUNDIALS-packaged SUNLinSol modules, and their constructor routines, may be found in chapter §11. Alternately, a user-supplied SUNLinearSolver module may be created and used. Specific information on how to create such user-provided modules may be found in §11.1.8.

Once this solver object has been constructed, the user should attach it to MRIStep via a call to MRIStepSetLinearSolver(). The first argument passed to this function is the MRIStep memory pointer returned by MRIStepCreate(); the second argument is the SUNLinearSolver object created above. The third argument is an optional SUNMatrix object to accompany matrix-based SUNLinearSolver inputs (for matrix-free linear solvers, the third argument should be NULL). A call to this function initializes the ARKLS linear solver interface, linking it to the MRIStep integrator, and allows the user to specify additional parameters and routines pertinent to their choice of linear solver.

int MRIStepSetLinearSolver(void *arkode_mem, SUNLinearSolver LS, SUNMatrix J)

This function specifies the SUNLinearSolver object that MRIStep should use, as well as a template Jacobian SUNMatrix object (if applicable).

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • LS – the SUNLinearSolver object to use.

  • J – the template Jacobian SUNMatrix object to use (or NULL if not applicable).

Return value:
  • ARKLS_SUCCESS if successful

  • ARKLS_MEM_NULL if the MRIStep memory was NULL

  • ARKLS_MEM_FAIL if there was a memory allocation failure

  • ARKLS_ILL_INPUT if ARKLS is incompatible with the provided LS or J input objects, or the current N_Vector module.

Notes: If LS is a matrix-free linear solver, then the J argument should be NULL.

If LS is a matrix-based linear solver, then the template Jacobian matrix J will be used in the solve process, so if additional storage is required within the SUNMatrix object (e.g. for factorization of a banded matrix), ensure that the input object is allocated with sufficient size (see the documentation of the particular SUNMATRIX type in §10 for further information).

When using sparse linear solvers, it is typically much more efficient to supply J so that it includes the full sparsity pattern of the Newton system matrices \(\mathcal{A} = I-\gamma J\), even if J itself has zeros in nonzero locations of \(I\). The reasoning for this is that \(\mathcal{A}\) is constructed in-place, on top of the user-specified values of J, so if the sparsity pattern in J is insufficient to store \(\mathcal{A}\) then it will need to be resized internally by MRIStep.

3.4.5.2.4. Nonlinear solver interface functions

When changing the nonlinear solver in MRIStep, after the call to MRIStepCreate() but before any calls to MRIStepEvolve(), the user’s program must create the appropriate SUNNonlinSol object and call MRIStepSetNonlinearSolver(), as documented below. If any calls to MRIStepEvolve() have been made, then MRIStep will need to be reinitialized by calling MRIStepReInit() to ensure that the nonlinear solver is initialized correctly before any subsequent calls to MRIStepEvolve().

The first argument passed to the routine MRIStepSetNonlinearSolver() is the MRIStep memory pointer returned by MRIStepCreate(); the second argument passed to this function is the desired SUNNonlinearSolver object to use for solving the nonlinear system for each implicit stage. A call to this function attaches the nonlinear solver to the main MRIStep integrator.

int MRIStepSetNonlinearSolver(void *arkode_mem, SUNNonlinearSolver NLS)

This function specifies the SUNNonlinearSolver object that MRIStep should use for implicit stage solves.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • NLS – the SUNNonlinearSolver object to use.

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

  • ARK_MEM_FAIL if there was a memory allocation failure

  • ARK_ILL_INPUT if MRIStep is incompatible with the provided NLS input object.

Notes: MRIStep will use the Newton SUNNonlinearSolver module by default; a call to this routine replaces that module with the supplied NLS object.

3.4.5.2.5. Rootfinding initialization function

As described in the section §3.2.12, while solving the IVP, ARKODE’s time-stepping modules have the capability to find the roots of a set of user-defined functions. In the MRIStep module root finding is performed between slow solution time steps only (i.e., it is not performed within the sub-stepping a fast time scales). To activate the root-finding algorithm, call the following function. This is normally called only once, prior to the first call to MRIStepEvolve(), but if the rootfinding problem is to be changed during the solution, MRIStepRootInit() can also be called prior to a continuation call to MRIStepEvolve().

int MRIStepRootInit(void *arkode_mem, int nrtfn, ARKRootFn g)

Initializes a rootfinding problem to be solved during the integration of the ODE system. It must be called after MRIStepCreate(), and before MRIStepEvolve().

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • nrtfn – number of functions \(g_i\), an integer \(\ge\) 0.

  • g – name of user-supplied function, of type ARKRootFn(), defining the functions \(g_i\) whose roots are sought.

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

  • ARK_MEM_FAIL if there was a memory allocation failure

  • ARK_ILL_INPUT if nrtfn is greater than zero but g = NULL.

Notes: To disable the rootfinding feature after it has already been initialized, or to free memory associated with MRIStep’s rootfinding module, call MRIStepRootInit with nrtfn = 0.

Similarly, if a new IVP is to be solved with a call to MRIStepReInit(), where the new IVP has no rootfinding problem but the prior one did, then call MRIStepRootInit with nrtfn = 0.

Rootfinding is only supported for the slow (outer) integrator and should not be actived for the fast (inner) integrator.

3.4.5.2.6. MRIStep solver function

This is the central step in the solution process – the call to perform the integration of the IVP. The input argument itask specifies one of two modes as to where MRIStep is to return a solution. These modes are modified if the user has set a stop time (with a call to the optional input function MRIStepSetStopTime()) or has requested rootfinding.

int MRIStepEvolve(void *arkode_mem, sunrealtype tout, N_Vector yout, sunrealtype *tret, int itask)

Integrates the ODE over an interval in \(t\).

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • tout – the next time at which a computed solution is desired.

  • yout – the computed solution vector.

  • tret – the time corresponding to yout (output).

  • itask – a flag indicating the job of the solver for the next user step.

    The ARK_NORMAL option causes the solver to take internal steps until it has just overtaken a user-specified output time, tout, in the direction of integration, i.e. \(t_{n-1} <\) tout \(\le t_{n}\) for forward integration, or \(t_{n} \le\) tout \(< t_{n-1}\) for backward integration. It will then compute an approximation to the solution \(y(tout)\) by interpolation (as described in §3.2.2).

    The ARK_ONE_STEP option tells the solver to only take a single internal step, \(y_{n-1} \to y_{n}\), and return the solution at that point, \(y_{n}\), in the vector yout.

Return value:
  • ARK_SUCCESS if successful.

  • ARK_ROOT_RETURN if MRIStepEvolve() succeeded, and found one or more roots. If the number of root functions, nrtfn, is greater than 1, call MRIStepGetRootInfo() to see which \(g_i\) were found to have a root at (*tret).

  • ARK_TSTOP_RETURN if MRIStepEvolve() succeeded and returned at tstop.

  • ARK_MEM_NULL if the arkode_mem argument was NULL.

  • ARK_NO_MALLOC if arkode_mem was not allocated.

  • ARK_ILL_INPUT if one of the inputs to MRIStepEvolve() is illegal, or some other input to the solver was either illegal or missing. Details will be provided in the error message. Typical causes of this failure:

    1. A component of the error weight vector became zero during internal time-stepping.

    2. The linear solver initialization function (called by the user after calling ARKStepCreate()) failed to set the linear solver-specific lsolve field in arkode_mem.

    3. A root of one of the root functions was found both at a point \(t\) and also very near \(t\).

  • ARK_TOO_MUCH_WORK if the solver took mxstep internal steps but could not reach tout. The default value for mxstep is MXSTEP_DEFAULT = 500.

  • ARK_CONV_FAILURE if convergence test failures occurred too many times (ark_maxncf) during one internal time step.

  • ARK_LINIT_FAIL if the linear solver’s initialization function failed.

  • ARK_LSETUP_FAIL if the linear solver’s setup routine failed in an unrecoverable manner.

  • ARK_LSOLVE_FAIL if the linear solver’s solve routine failed in an unrecoverable manner.

  • ARK_VECTOROP_ERR a vector operation error occurred.

  • ARK_INNERSTEP_FAILED if the inner stepper returned with an unrecoverable error. The value returned from the inner stepper can be obtained with MRIStepGetLastInnerStepFlag().

  • ARK_INVALID_TABLE if an invalid coupling table was provided.

Notes:

The input vector yout can use the same memory as the vector y0 of initial conditions that was passed to MRIStepCreate().

In ARK_ONE_STEP mode, tout is used only on the first call, and only to get the direction and a rough scale of the independent variable.

All failure return values are negative and so testing the return argument for negative values will trap all MRIStepEvolve() failures.

Since interpolation may reduce the accuracy in the reported solution, if full method accuracy is desired the user should issue a call to MRIStepSetStopTime() before the call to MRIStepEvolve() to specify a fixed stop time to end the time step and return to the user. Upon return from MRIStepEvolve(), a copy of the internal solution \(y_{n}\) will be returned in the vector yout. Once the integrator returns at a tstop time, any future testing for tstop is disabled (and can be re-enabled only though a new call to MRIStepSetStopTime()).

On any error return in which one or more internal steps were taken by MRIStepEvolve(), the returned values of tret and yout correspond to the farthest point reached in the integration. On all other error returns, tret and yout are left unchanged from those provided to the routine.

3.4.5.2.7. Optional input functions

There are numerous optional input parameters that control the behavior of MRIStep, each of which may be modified from its default value through calling an appropriate input function. The following tables list all optional input functions, grouped by which aspect of MRIStep they control. Detailed information on the calling syntax and arguments for each function are then provided following each table.

The optional inputs are grouped into the following categories:

For the most casual use of MRIStep, relying on the default set of solver parameters, the reader can skip to the section on user-supplied functions, §3.4.6.

We note that, on an error return, all of the optional input functions send an error message to the error handler function. All error return values are negative, so a test on the return arguments for negative values will catch all errors. Finally, a call to an MRIStepSet*** function can generally be made from the user’s calling program at any time and, if successful, takes effect immediately. MRIStepSet*** functions that cannot be called at any time note this in the “Notes:” section of the function documentation.

3.4.5.2.7.1. Optional inputs for MRIStep

Table 3.16 Optional inputs for MRIStep

Optional input

Function name

Default

Return MRIStep solver parameters to their defaults

MRIStepSetDefaults()

internal

Set dense output interpolation type

MRIStepSetInterpolantType()

ARK_INTERP_HERMITE

Set dense output polynomial degree

MRIStepSetInterpolantDegree()

5

Supply a pointer to a diagnostics output file

MRIStepSetDiagnostics()

NULL

Run with fixed-step sizes

MRIStepSetFixedStep()

required

Maximum no. of warnings for \(t_n+h = t_n\)

MRIStepSetMaxHnilWarns()

10

Maximum no. of internal steps before tout

MRIStepSetMaxNumSteps()

500

Set a value for \(t_{stop}\)

MRIStepSetStopTime()

undefined

Interpolate at \(t_{stop}\)

MRIStepSetInterpolateStopTime()

SUNFALSE

Disable the stop time

MRIStepClearStopTime()

N/A

Supply a pointer for user data

MRIStepSetUserData()

NULL

Supply a function to be called prior to the inner integration

MRIStepSetPreInnerFn()

NULL

Supply a function to be called after the inner integration

MRIStepSetPostInnerFn()

NULL

int MRIStepSetDefaults(void *arkode_mem)

Resets all optional input parameters to MRIStep’s original default values.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes: This function does not change problem-defining function pointers fs and ff or the user_data pointer. It also does not affect any data structures or options related to root-finding (those can be reset using MRIStepRootInit()).

int MRIStepSetInterpolantType(void *arkode_mem, int itype)

Specifies use of the Lagrange or Hermite interpolation modules (used for dense output – interpolation of solution output values and implicit method predictors).

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • itype – requested interpolant type (ARK_INTERP_HERMITE or ARK_INTERP_LAGRANGE)

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_MEM_FAIL if the interpolation module cannot be allocated

  • ARK_ILL_INPUT if the itype argument is not recognized or the interpolation module has already been initialized

Notes: The Hermite interpolation module is described in §3.2.2.1, and the Lagrange interpolation module is described in §3.2.2.2.

This routine frees any previously-allocated interpolation module, and re-creates one according to the specified argument. Thus any previous calls to MRIStepSetInterpolantDegree() will be nullified.

This routine must be called after the call to MRIStepCreate(). After the first call to MRIStepEvolve() the interpolation type may not be changed without first calling MRIStepReInit().

If this routine is not called, the Hermite interpolation module will be used.

int MRIStepSetInterpolantDegree(void *arkode_mem, int degree)

Specifies the degree of the polynomial interpolant used for dense output (i.e. interpolation of solution output values and implicit method predictors).

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • degree – requested polynomial degree.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory or interpolation module are NULL

  • ARK_INTERP_FAIL if this is called after MRIStepEvolve()

  • ARK_ILL_INPUT if an argument has an illegal value or the interpolation module has already been initialized

Notes: Allowed values are between 0 and 5.

This routine should be called after MRIStepCreate() and before MRIStepEvolve(). After the first call to MRIStepEvolve() the interpolation degree may not be changed without first calling MRIStepReInit().

If a user calls both this routine and MRIStepSetInterpolantType(), then MRIStepSetInterpolantType() must be called first.

Since the accuracy of any polynomial interpolant is limited by the accuracy of the time-step solutions on which it is based, the actual polynomial degree that is used by MRIStep will be the minimum of \(q-1\) and the input degree, for \(q > 1\) where \(q\) is the order of accuracy for the time integration method.

Changed in version 5.5.1: When \(q=1\), a linear interpolant is the default to ensure values obtained by the integrator are returned at the ends of the time interval.

int MRIStepSetDenseOrder(void *arkode_mem, int dord)

This function is deprecated, and will be removed in a future release. Users should transition to calling MRIStepSetInterpolantDegree() instead.

int MRIStepSetDiagnostics(void *arkode_mem, FILE *diagfp)

Specifies the file pointer for a diagnostics file where all MRIStep step adaptivity and solver information is written.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • diagfp – pointer to the diagnostics output file.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes: This parameter can be stdout or stderr, although the suggested approach is to specify a pointer to a unique file opened by the user and returned by fopen. If not called, or if called with a NULL file pointer, all diagnostics output is disabled.

When run in parallel, only one process should set a non-NULL value for this pointer, since statistics from all processes would be identical.

Deprecated since version 5.2.0: Use SUNLogger_SetInfoFilename() instead.

int MRIStepSetFixedStep(void *arkode_mem, sunrealtype hs)

Set the slow step size used within MRIStep for the following internal step(s).

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • hs – value of the outer (slow) step size.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes:

The step sizes used by the inner (fast) stepper may be controlled through calling the appropriate “Set” routines on the inner integrator.

int MRIStepSetMaxHnilWarns(void *arkode_mem, int mxhnil)

Specifies the maximum number of messages issued by the solver to warn that \(t+h=t\) on the next internal step, before MRIStep will instead return with an error.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • mxhnil – maximum allowed number of warning messages \((>0)\).

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes: The default value is 10; set mxhnil to zero to specify this default.

A negative value indicates that no warning messages should be issued.

int MRIStepSetMaxNumSteps(void *arkode_mem, long int mxsteps)

Specifies the maximum number of steps to be taken by the solver in its attempt to reach the next output time, before MRIStep will return with an error.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • mxsteps – maximum allowed number of internal steps.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes: Passing mxsteps = 0 results in MRIStep using the default value (500).

Passing mxsteps < 0 disables the test (not recommended).

int MRIStepSetStopTime(void *arkode_mem, sunrealtype tstop)

Specifies the value of the independent variable \(t\) past which the solution is not to proceed.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • tstop – stopping time for the integrator.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes:

The default is that no stop time is imposed.

Once the integrator returns at a stop time, any future testing for tstop is disabled (and can be reenabled only though a new call to MRIStepSetStopTime()).

A stop time not reached before a call to MRIStepReInit() or MRIStepReset() will remain active but can be disabled by calling MRIStepClearStopTime().

int MRIStepSetInterpolateStopTime(void *arkode_mem, sunbooleantype interp)

Specifies that the output solution should be interpolated when the current \(t\) equals the specified tstop (instead of merely copying the internal solution \(y_n\)).

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • interp – flag indicating to use interpolation (1) or copy (0).

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the ARKStep memory is NULL

New in version 5.6.0.

int MRIStepClearStopTime(void *arkode_mem)

Disables the stop time set with MRIStepSetStopTime().

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

Notes:

The stop time can be reenabled though a new call to MRIStepSetStopTime().

New in version 5.5.1.

int MRIStepSetUserData(void *arkode_mem, void *user_data)

Specifies the user data block user_data for the outer integrator and attaches it to the main MRIStep memory block.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • user_data – pointer to the user data.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes: If specified, the pointer to user_data is passed to all user-supplied functions called by the outer integrator for which it is an argument; otherwise NULL is passed.

To attach a user data block to the inner integrator call the appropriate SetUserData function for the inner integrator memory structure (e.g., ARKStepSetUserData() if the inner stepper is ARKStep). This pointer may be the same as or different from the pointer attached to the outer integrator depending on what is required by the user code.

int MRIStepSetPreInnerFn(void *arkode_mem, MRIStepPreInnerFn prefn)

Specifies the function called before each inner integration.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • prefn – the name of the C function (of type MRIStepPreInnerFn()) defining pre inner integration function.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

int MRIStepSetPostInnerFn(void *arkode_mem, MRIStepPostInnerFn postfn)

Specifies the function called after each inner integration.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • postfn – the name of the C function (of type MRIStepPostInnerFn()) defining post inner integration function.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

3.4.5.2.7.2. Optional inputs for IVP method selection

Table 3.17 Optional inputs for IVP method selection

Optional input

Function name

Default

Select the default MRI method of a given order

MRIStepSetOrder()

3

Set MRI coupling coefficients

MRIStepSetCoupling()

internal

int MRIStepSetOrder(void *arkode_mem, int ord)

Select the default MRI method of a given order.

The default order is 3. An order less than 3 or greater than 4 will result in using the default.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • ord – the method order.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

int MRIStepSetCoupling(void *arkode_mem, MRIStepCoupling C)

Specifies a customized set of slow-to-fast coupling coefficients for the MRI method.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • C – the table of coupling coefficients for the MRI method.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes:

For a description of the MRIStepCoupling type and related functions for creating Butcher tables see §3.4.5.3.

3.4.5.2.7.3. Optional inputs for implicit stage solves

The mathematical explanation for the nonlinear solver strategies used by MRIStep, including how each of the parameters below is used within the code, is provided in §3.2.11.1.

Optional input

Function name

Default

Specify linearly implicit \(f^I\)

MRIStepSetLinear()

SUNFALSE

Specify nonlinearly implicit \(f^I\)

MRIStepSetNonlinear()

SUNTRUE

Implicit predictor method

MRIStepSetPredictorMethod()

0

Maximum number of nonlinear iterations

MRIStepSetMaxNonlinIters()

3

Coefficient in the nonlinear convergence test

MRIStepSetNonlinConvCoef()

0.1

Nonlinear convergence rate constant

MRIStepSetNonlinCRDown()

0.3

Nonlinear residual divergence ratio

MRIStepSetNonlinRDiv()

2.3

User-provided implicit stage predictor

MRIStepSetStagePredictFn()

NULL

RHS function for nonlinear system evaluations

MRIStepSetNlsRhsFn()

NULL

Specify if \(f^I\) is deduced after a nonlinear solve

MRIStepSetDeduceImplicitRhs()

SUNFALSE

int MRIStepSetLinear(void *arkode_mem, int timedepend)

Specifies that the implicit slow right-hand side function, \(f^I(t,y)\) is linear in \(y\).

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • timedepend – flag denoting whether the Jacobian of \(f^I(t,y)\) is time-dependent (1) or not (0). Alternately, when using a matrix-free iterative linear solver this flag denotes time dependence of the preconditioner.

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes: Tightens the linear solver tolerances and takes only a single Newton iteration. Calls MRIStepSetDeltaGammaMax() to enforce Jacobian recomputation when the step size ratio changes by more than 100 times the unit roundoff (since nonlinear convergence is not tested). Only applicable when used in combination with the modified or inexact Newton iteration (not the fixed-point solver).

The only SUNDIALS-provided SUNNonlinearSolver module that is compatible with the MRIStepSetLinear() option is the Newton solver.

int MRIStepSetNonlinear(void *arkode_mem)

Specifies that the implicit slow right-hand side function, \(f^I(t,y)\) is nonlinear in \(y\).

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes: This is the default behavior of MRIStep, so the function is primarily useful to undo a previous call to MRIStepSetLinear(). Calls MRIStepSetDeltaGammaMax() to reset the step size ratio threshold to the default value.

int MRIStepSetPredictorMethod(void *arkode_mem, int method)

Specifies the method to use for predicting implicit solutions.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • method – method choice (0 \(\le\) method \(\le\) 4):

    • 0 is the trivial predictor,

    • 1 is the maximum order (dense output) predictor,

    • 2 is the variable order predictor, that decreases the polynomial degree for more distant RK stages,

    • 3 is the cutoff order predictor, that uses the maximum order for early RK stages, and a first-order predictor for distant RK stages,

    • 4 is the bootstrap predictor, that uses a second-order predictor based on only information within the current step. deprecated

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes: The default value is 0. If method is set to an undefined value, this default predictor will be used.

The “bootstrap” predictor (option 4 above) has been deprecated, and will be removed from a future release.

int MRIStepSetMaxNonlinIters(void *arkode_mem, int maxcor)

Specifies the maximum number of nonlinear solver iterations permitted per slow MRI stage within each time step.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • maxcor – maximum allowed solver iterations per stage \((>0)\).

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value or if the SUNNONLINSOL module is NULL

  • ARK_NLS_OP_ERR if the SUNNONLINSOL object returned a failure flag

Notes: The default value is 3; set maxcor \(\le 0\) to specify this default.

int MRIStepSetNonlinConvCoef(void *arkode_mem, sunrealtype nlscoef)

Specifies the safety factor used within the nonlinear solver convergence test.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • nlscoef – coefficient in nonlinear solver convergence test \((>0.0)\).

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes: The default value is 0.1; set nlscoef \(\le 0\) to specify this default.

int MRIStepSetNonlinCRDown(void *arkode_mem, sunrealtype crdown)

Specifies the constant used in estimating the nonlinear solver convergence rate.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • crdown – nonlinear convergence rate estimation constant (default is 0.3).

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes: Any non-positive parameter will imply a reset to the default value.

int MRIStepSetNonlinRDiv(void *arkode_mem, sunrealtype rdiv)

Specifies the nonlinear correction threshold beyond which the iteration will be declared divergent.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • rdiv – tolerance on nonlinear correction size ratio to declare divergence (default is 2.3).

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes: Any non-positive parameter will imply a reset to the default value.

int MRIStepSetStagePredictFn(void *arkode_mem, ARKStagePredictFn PredictStage)

Sets the user-supplied function to update the implicit stage predictor prior to execution of the nonlinear or linear solver algorithms that compute the implicit stage solution.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • PredictStage – name of user-supplied predictor function. If NULL, then any previously-provided stage prediction function will be disabled.

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

Notes: See §3.4.6.6 for more information on this user-supplied routine.

int MRIStepSetNlsRhsFn(void *arkode_mem, ARKRhsFn nls_fs)

Specifies an alternative implicit slow right-hand side function for evaluating \(f^I(t,y)\) within nonlinear system function evaluations.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • nls_fs – the alternative C function for computing the right-hand side function \(f^I(t,y)\) in the ODE.

Return value:
  • ARK_SUCCESS if successful.

  • ARK_MEM_NULL if the MRIStep memory was NULL.

Notes: The default is to use the implicit slow right-hand side function provided to MRIStepCreate() in nonlinear system functions. If the input implicit slow right-hand side function is NULL, the default is used.

When using a non-default nonlinear solver, this function must be called after MRIStepSetNonlinearSolver().

int MRIStepSetDeduceImplicitRhs(void *arkode_mem, sunbooleantype deduce)

Specifies if implicit stage derivatives are deduced without evaluating \(f^I\). See §3.2.11.1 for more details.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • deduce – If SUNFALSE (default), the stage derivative is obtained by evaluating \(f^I\) with the stage solution returned from the nonlinear solver. If SUNTRUE, the stage derivative is deduced without an additional evaluation of \(f^I\).

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

New in version 5.2.0.

3.4.5.2.7.4. Linear solver interface optional input functions

The mathematical explanation of the linear solver methods available to MRIStep is provided in §3.2.11.2. We group the user-callable routines into four categories: general routines concerning the update frequency for matrices and/or preconditioners, optional inputs for matrix-based linear solvers, optional inputs for matrix-free linear solvers, and optional inputs for iterative linear solvers. We note that the matrix-based and matrix-free groups are mutually exclusive, whereas the “iterative” tag can apply to either case.

3.4.5.2.7.4.1. Optional inputs for the ARKLS linear solver interface

As discussed in §3.2.11.2.3, ARKODE strives to reuse matrix and preconditioner data for as many solves as possible to amortize the high costs of matrix construction and factorization. To that end, MRIStep provides user-callable routines to modify this behavior. Recall that the Newton system matrices that arise within an implicit stage solve are \({\mathcal A}(t,z) \approx I - \gamma J(t,z)\), where the implicit right-hand side function has Jacobian matrix \(J(t,z) = \frac{\partial f^I(t,z)}{\partial z}\).

The matrix or preconditioner for \({\mathcal A}\) can only be updated within a call to the linear solver ‘setup’ routine. In general, the frequency with which the linear solver setup routine is called may be controlled with the msbp argument to MRIStepSetLSetupFrequency(). When this occurs, the validity of \({\mathcal A}\) for successive time steps intimately depends on whether the corresponding \(\gamma\) and \(J\) inputs remain valid.

At each call to the linear solver setup routine the decision to update \(\mathcal{A}\) with a new value of \(\gamma\), and to reuse or reevaluate Jacobian information, depends on several factors including:

  • the success or failure of previous solve attempts,

  • the success or failure of the previous time step attempts,

  • the change in \(\gamma\) from the value used when constructing \(\mathcal{A}\), and

  • the number of steps since Jacobian information was last evaluated.

The frequency with which to update Jacobian information can be controlled with the msbj argument to MRIStepSetJacEvalFrequency(). We note that this is only checked within calls to the linear solver setup routine, so values msbj \(<\) msbp do not make sense. For linear-solvers with user-supplied preconditioning the above factors are used to determine whether to recommend updating the Jacobian information in the preconditioner (i.e., whether to set jok to SUNFALSE in calling the user-supplied ARKLsPrecSetupFn()). For matrix-based linear solvers these factors determine whether the matrix \(J(t,y) = \frac{\partial f^I(t,y)}{\partial y}\) should be updated (either with an internal finite difference approximation or a call to the user-supplied ARKLsJacFn); if not then the previous value is reused and the system matrix \({\mathcal A}(t,y) \approx I - \gamma J(t,y)\) is recomputed using the current \(\gamma\) value.

Optional input

Function name

Default

Max change in step signaling new \(J\)

MRIStepSetDeltaGammaMax()

0.2

Linear solver setup frequency

MRIStepSetLSetupFrequency()

20

Jacobian / preconditioner update frequency

MRIStepSetJacEvalFrequency()

51

int MRIStepSetDeltaGammaMax(void *arkode_mem, sunrealtype dgmax)

Specifies a scaled step size ratio tolerance, beyond which the linear solver setup routine will be signaled.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • dgmax – tolerance on step size ratio change before calling linear solver setup routine (default is 0.2).

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes: Any non-positive parameter will imply a reset to the default value.

int MRIStepSetLSetupFrequency(void *arkode_mem, int msbp)

Specifies the frequency of calls to the linear solver setup routine.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • msbp – the linear solver setup frequency.

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

Notes: Positive values of msbp specify the linear solver setup frequency. For example, an input of 1 means the setup function will be called every time step while an input of 2 means it will be called called every other time step. If msbp is 0, the default value of 20 will be used. A negative value forces a linear solver step at each implicit stage.

int MRIStepSetJacEvalFrequency(void *arkode_mem, long int msbj)

Specifies the frequency for recomputing the Jacobian or recommending a preconditioner update.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • msbj – the Jacobian re-computation or preconditioner update frequency.

Return value:
  • ARKLS_SUCCESS if successful.

  • ARKLS_MEM_NULL if the MRIStep memory was NULL.

  • ARKLS_LMEM_NULL if the linear solver memory was NULL.

Notes: The Jacobian update frequency is only checked within calls to the linear solver setup routine, as such values of msbj \(<\) msbp will result in recomputing the Jacobian every msbp steps. See MRIStepSetLSetupFrequency() for setting the linear solver setup frequency msbp.

Passing a value msbj \(\le 0\) indicates to use the default value of 50.

This function must be called after the ARKLS system solver interface has been initialized through a call to MRIStepSetLinearSolver().

3.4.5.2.7.4.2. Optional inputs for matrix-based SUNLinearSolver modules

Optional input

Function name

Default

Jacobian function

MRIStepSetJacFn()

DQ

Linear system function

MRIStepSetLinSysFn()

internal

Enable or disable linear solution scaling

MRIStepSetLinearSolutionScaling()

on

When using matrix-based linear solver modules, the ARKLS solver interface needs a function to compute an approximation to the Jacobian matrix \(J(t,y)\) or the linear system \(I - \gamma J\). The function to evaluate the Jacobian must be of type ARKLsJacFn(). The user can supply a custom Jacobian function, or if using a dense or banded \(J\) can use the default internal difference quotient approximation that comes with the ARKLS interface. At present, we do not supply a corresponding routine to approximate Jacobian entries in sparse matrices \(J\). To specify a user-supplied Jacobian function jac, MRIStep provides the function MRIStepSetJacFn(). Alternatively, a function of type ARKLsLinSysFn() can be provided to evaluate the matrix \(I - \gamma J\). By default, ARKLS uses an internal linear system function leveraging the SUNMATRIX API to form the matrix \(I - \gamma J\). To specify a user-supplied linear system function linsys, MRIStep provides the function MRIStepSetLinSysFn(). In either case the matrix information will be updated infrequently to reduce matrix construction and, with direct solvers, factorization costs. As a result the value of \(\gamma\) may not be current and a scaling factor is applied to the solution of the linear system to account for lagged value of \(\gamma\). See §11.2.1 for more details. The function MRIStepSetLinearSolutionScaling() can be used to disable this scaling when necessary, e.g., when providing a custom linear solver that updates the matrix using the current \(\gamma\) as part of the solve.

The ARKLS interface passes the user data pointer to the Jacobian and linear system functions. This allows the user to create an arbitrary structure with relevant problem data and access it during the execution of the user-supplied Jacobian or linear system functions, without using global data in the program. The user data pointer may be specified through MRIStepSetUserData().

int MRIStepSetJacFn(void *arkode_mem, ARKLsJacFn jac)

Specifies the Jacobian approximation routine to be used for the matrix-based solver with the ARKLS interface.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • jac – name of user-supplied Jacobian approximation function.

Return value:
  • ARKLS_SUCCESS if successful

  • ARKLS_MEM_NULL if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL if the linear solver memory was NULL

Notes: This routine must be called after the ARKLS linear solver interface has been initialized through a call to MRIStepSetLinearSolver().

By default, ARKLS uses an internal difference quotient function for dense and band matrices. If NULL is passed in for jac, this default is used. An error will occur if no jac is supplied when using other matrix types.

The function type ARKLsJacFn() is described in §3.4.6.

int MRIStepSetLinSysFn(void *arkode_mem, ARKLsLinSysFn linsys)

Specifies the linear system approximation routine to be used for the matrix-based solver with the ARKLS interface.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • linsys – name of user-supplied linear system approximation function.

Return value:
  • ARKLS_SUCCESS if successful

  • ARKLS_MEM_NULL if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL if the linear solver memory was NULL

Notes: This routine must be called after the ARKLS linear solver interface has been initialized through a call to MRIStepSetLinearSolver().

By default, ARKLS uses an internal linear system function that leverages the SUNMATRIX API to form the system \(I - \gamma J\). If NULL is passed in for linsys, this default is used.

The function type ARKLsLinSysFn() is described in §3.4.6.

int MRIStepSetLinearSolutionScaling(void *arkode_mem, sunbooleantype onoff)

Enables or disables scaling the linear system solution to account for a change in \(\gamma\) in the linear system. For more details see §11.2.1.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • onoff – flag to enable (SUNTRUE) or disable (SUNFALSE) scaling

Return value:
  • ARKLS_SUCCESS if successful

  • ARKLS_MEM_NULL if the MRIStep memory was NULL

  • ARKLS_ILL_INPUT if the attached linear solver is not matrix-based

Notes: Linear solution scaling is enabled by default when a matrix-based linear solver is attached.

3.4.5.2.7.4.3. Optional inputs for matrix-free SUNLinearSolver modules

Optional input

Function name

Default

\(Jv\) functions (jtimes and jtsetup)

MRIStepSetJacTimes()

DQ, none

\(Jv\) DQ rhs function (jtimesRhsFn)

MRIStepSetJacTimesRhsFn()

fs

As described in §3.2.11.2, when solving the Newton linear systems with matrix-free methods, the ARKLS interface requires a jtimes function to compute an approximation to the product between the Jacobian matrix \(J(t,y)\) and a vector \(v\). The user can supply a custom Jacobian-times-vector approximation function, or use the default internal difference quotient function that comes with the ARKLS interface.

A user-defined Jacobian-vector function must be of type ARKLsJacTimesVecFn and can be specified through a call to MRIStepSetJacTimes() (see §3.4.6 for specification details). As with the user-supplied preconditioner functions, the evaluation and processing of any Jacobian-related data needed by the user’s Jacobian-times-vector function is done in the optional user-supplied function of type ARKLsJacTimesSetupFn (see §3.4.6 for specification details). As with the preconditioner functions, a pointer to the user-defined data structure, user_data, specified through MRIStepSetUserData() (or a NULL pointer otherwise) is passed to the Jacobian-times-vector setup and product functions each time they are called.

int MRIStepSetJacTimes(void *arkode_mem, ARKLsJacTimesSetupFn jtsetup, ARKLsJacTimesVecFn jtimes)

Specifies the Jacobian-times-vector setup and product functions.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • jtsetup – user-defined Jacobian-vector setup function. Pass NULL if no setup is necessary.

  • jtimes – user-defined Jacobian-vector product function.

Return value:
  • ARKLS_SUCCESS if successful.

  • ARKLS_MEM_NULL if the MRIStep memory was NULL.

  • ARKLS_LMEM_NULL if the linear solver memory was NULL.

  • ARKLS_ILL_INPUT if an input has an illegal value.

  • ARKLS_SUNLS_FAIL if an error occurred when setting up the Jacobian-vector product in the SUNLinearSolver object used by the ARKLS interface.

Notes: The default is to use an internal finite difference quotient for jtimes and to leave out jtsetup. If NULL is passed to jtimes, these defaults are used. A user may specify non-NULL jtimes and NULL jtsetup inputs.

This function must be called after the ARKLS system solver interface has been initialized through a call to MRIStepSetLinearSolver().

The function types ARKLsJacTimesSetupFn and ARKLsJacTimesVecFn are described in §3.4.6.

When using the internal difference quotient the user may optionally supply an alternative implicit right-hand side function for use in the Jacobian-vector product approximation by calling MRIStepSetJacTimesRhsFn(). The alternative implicit right-hand side function should compute a suitable (and differentiable) approximation to the \(f^I\) function provided to MRIStepCreate(). For example, as done in [46], the alternative function may use lagged values when evaluating a nonlinearity in \(f^I\) to avoid differencing a potentially non-differentiable factor.

int MRIStepSetJacTimesRhsFn(void *arkode_mem, ARKRhsFn jtimesRhsFn)

Specifies an alternative implicit right-hand side function for use in the internal Jacobian-vector product difference quotient approximation.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • jtimesRhsFn – the name of the C function (of type ARKRhsFn()) defining the alternative right-hand side function.

Return value:
  • ARKLS_SUCCESS if successful.

  • ARKLS_MEM_NULL if the MRIStep memory was NULL.

  • ARKLS_LMEM_NULL if the linear solver memory was NULL.

  • ARKLS_ILL_INPUT if an input has an illegal value.

Notes: The default is to use the implicit right-hand side function provided to MRIStepCreate() in the internal difference quotient. If the input implicit right-hand side function is NULL, the default is used.

This function must be called after the ARKLS system solver interface has been initialized through a call to MRIStepSetLinearSolver().

3.4.5.2.7.4.4. Optional inputs for iterative SUNLinearSolver modules

Optional input

Function name

Default

Newton preconditioning functions

MRIStepSetPreconditioner()

NULL, NULL

Newton linear and nonlinear tolerance ratio

MRIStepSetEpsLin()

0.05

Newton linear solve tolerance conversion factor

MRIStepSetLSNormFactor()

vector length

As described in §3.2.11.2, when using an iterative linear solver the user may supply a preconditioning operator to aid in solution of the system. This operator consists of two user-supplied functions, psetup and psolve, that are supplied to MRIStep using the function MRIStepSetPreconditioner(). The psetup function supplied to these routines should handle evaluation and preprocessing of any Jacobian data needed by the user’s preconditioner solve function, psolve. The user data pointer received through MRIStepSetUserData() (or a pointer to NULL if user data was not specified) is passed to the psetup and psolve functions. This allows the user to create an arbitrary structure with relevant problem data and access it during the execution of the user-supplied preconditioner functions without using global data in the program.

Also, as described in §3.2.11.3.2, the ARKLS interface requires that iterative linear solvers stop when the norm of the preconditioned residual satisfies

\[\|r\| \le \frac{\epsilon_L \epsilon}{10}\]

where the default \(\epsilon_L = 0.05\), which may be modified by the user through the MRIStepSetEpsLin() function.

int MRIStepSetPreconditioner(void *arkode_mem, ARKLsPrecSetupFn psetup, ARKLsPrecSolveFn psolve)

Specifies the user-supplied preconditioner setup and solve functions.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • psetup – user defined preconditioner setup function. Pass NULL if no setup is needed.

  • psolve – user-defined preconditioner solve function.

Return value:
  • ARKLS_SUCCESS if successful.

  • ARKLS_MEM_NULL if the MRIStep memory was NULL.

  • ARKLS_LMEM_NULL if the linear solver memory was NULL.

  • ARKLS_ILL_INPUT if an input has an illegal value.

  • ARKLS_SUNLS_FAIL if an error occurred when setting up preconditioning in the SUNLinearSolver object used by the ARKLS interface.

Notes: The default is NULL for both arguments (i.e., no preconditioning).

This function must be called after the ARKLS system solver interface has been initialized through a call to MRIStepSetLinearSolver().

Both of the function types ARKLsPrecSetupFn() and ARKLsPrecSolveFn() are described in §3.4.6.

int MRIStepSetEpsLin(void *arkode_mem, sunrealtype eplifac)

Specifies the factor by which the tolerance on the nonlinear iteration is multiplied to get a tolerance on the linear iteration.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • eplifac – linear convergence safety factor.

Return value:
  • ARKLS_SUCCESS if successful.

  • ARKLS_MEM_NULL if the MRIStep memory was NULL.

  • ARKLS_LMEM_NULL if the linear solver memory was NULL.

  • ARKLS_ILL_INPUT if an input has an illegal value.

Notes: Passing a value eplifac \(\le 0\) indicates to use the default value of 0.05.

This function must be called after the ARKLS system solver interface has been initialized through a call to MRIStepSetLinearSolver().

int MRIStepSetLSNormFactor(void *arkode_mem, sunrealtype nrmfac)

Specifies the factor to use when converting from the integrator tolerance (WRMS norm) to the linear solver tolerance (L2 norm) for Newton linear system solves e.g., tol_L2 = fac * tol_WRMS.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • nrmfac – the norm conversion factor. If nrmfac is:

    \(> 0\) then the provided value is used.

    \(= 0\) then the conversion factor is computed using the vector length i.e., nrmfac = sqrt(N_VGetLength(y)) (default).

    \(< 0\) then the conversion factor is computed using the vector dot product i.e., nrmfac = sqrt(N_VDotProd(v,v)) where all the entries of v are one.

Return value:
  • ARK_SUCCESS if successful.

  • ARK_MEM_NULL if the MRIStep memory was NULL.

Notes: This function must be called after the ARKLS system solver interface has been initialized through a call to MRIStepSetLinearSolver().

3.4.5.2.7.5. Rootfinding optional input functions

The following functions can be called to set optional inputs to control the rootfinding algorithm, the mathematics of which are described in the section §3.2.12.

Optional input

Function name

Default

Direction of zero-crossings to monitor

MRIStepSetRootDirection()

both

Disable inactive root warnings

MRIStepSetNoInactiveRootWarn()

enabled

int MRIStepSetRootDirection(void *arkode_mem, int *rootdir)

Specifies the direction of zero-crossings to be located and returned.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • rootdir – state array of length nrtfn, the number of root functions \(g_i\) (the value of nrtfn was supplied in the call to MRIStepRootInit()). If rootdir[i] == 0 then crossing in either direction for \(g_i\) should be reported. A value of +1 or -1 indicates that the solver should report only zero-crossings where \(g_i\) is increasing or decreasing, respectively.

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

  • ARK_ILL_INPUT if an argument has an illegal value

Notes: The default behavior is to monitor for both zero-crossing directions.

int MRIStepSetNoInactiveRootWarn(void *arkode_mem)

Disables issuing a warning if some root function appears to be identically zero at the beginning of the integration.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory is NULL

Notes: MRIStep will not report the initial conditions as a possible zero-crossing (assuming that one or more components \(g_i\) are zero at the initial time). However, if it appears that some \(g_i\) is identically zero at the initial time (i.e., \(g_i\) is zero at the initial time and after the first step), MRIStep will issue a warning which can be disabled with this optional input function.

3.4.5.2.8. Interpolated output function

An optional function MRIStepGetDky() is available to obtain additional values of solution-related quantities. This function should only be called after a successful return from MRIStepEvolve(), as it provides interpolated values either of \(y\) or of its derivatives (up to the 3rd derivative) interpolated to any value of \(t\) in the last internal step taken by MRIStepEvolve(). Internally, this “dense output” or “continuous extension” algorithm is identical to the algorithm used for the maximum order implicit predictors, described in §3.2.11.5.2, except that derivatives of the polynomial model may be evaluated upon request.

int MRIStepGetDky(void *arkode_mem, sunrealtype t, int k, N_Vector dky)

Computes the k-th derivative of the function \(y\) at the time t, i.e. \(y^{(k)}(t)\), for values of the independent variable satisfying \(t_n-h_n \le t \le t_n\), with \(t_n\) as current internal time reached, and \(h_n\) is the last internal step size successfully used by the solver. This routine uses an interpolating polynomial of degree min(degree, 5), where degree is the argument provided to MRIStepSetInterpolantDegree(). The user may request k in the range {0,…, min(degree, kmax)} where kmax depends on the choice of interpolation module. For Hermite interpolants kmax = 5 and for Lagrange interpolants kmax = 3.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • t – the value of the independent variable at which the derivative is to be evaluated.

  • k – the derivative order requested.

  • dky – output vector (must be allocated by the user).

Return value:

  • ARK_SUCCESS if successful

  • ARK_BAD_K if k is not in the range {0,…, min(degree, kmax)}.

  • ARK_BAD_T if t is not in the interval \([t_n-h_n, t_n]\)

  • ARK_BAD_DKY if the dky vector was NULL

  • ARK_MEM_NULL if the MRIStep memory is NULL

Notes: It is only legal to call this function after a successful return from MRIStepEvolve().

A user may access the values \(t_n\) and \(h_n\) via the functions MRIStepGetCurrentTime() and MRIStepGetLastStep(), respectively.

3.4.5.2.9. Optional output functions

MRIStep provides an extensive set of functions that can be used to obtain solver performance information. We organize these into groups:

  1. General MRIStep output routines are in §3.4.5.2.9.1,

  2. MRIStep implicit solver output routines are in §3.4.5.2.9.2,

  3. Linear solver output routines are in §3.4.5.2.9.4 and

  4. General usability routines (e.g. to print the current MRIStep parameters, or output the current coupling table) are in §3.4.5.2.9.5.

  5. Output routines regarding root-finding results are in §3.4.5.2.9.3,

Following each table, we elaborate on each function.

Some of the optional outputs, especially the various counters, can be very useful in determining the efficiency of various methods inside MRIStep. For example:

  • The number of steps and right-hand side evaluations at both the slow and fast time scales provide a rough measure of the overall cost of a given run, and can be compared between runs with different solver options to suggest which set of options is the most efficient.

  • The ratio nniters/nsteps measures the performance of the nonlinear iteration in solving the nonlinear systems at each implicit stage, providing a measure of the degree of nonlinearity in the problem. Typical values of this for a Newton solver on a general problem range from 1.1 to 1.8.

  • When using a Newton nonlinear solver, the ratio njevals/nniters (when using a direct linear solver), and the ratio nliters/nniters (when using an iterative linear solver) can indicate the quality of the approximate Jacobian or preconditioner being used. For example, if this ratio is larger for a user-supplied Jacobian or Jacobian-vector product routine than for the difference-quotient routine, it can indicate that the user-supplied Jacobian is inaccurate.

It is therefore recommended that users retrieve and output these statistics following each run, and take some time to investigate alternate solver options that will be more optimal for their particular problem of interest.

3.4.5.2.9.1. Main solver optional output functions

Table 3.18 Main solver optional output functions

Optional output

Function name

Size of MRIStep real and integer workspaces

MRIStepGetWorkSpace()

Cumulative number of internal steps

MRIStepGetNumSteps()

Step size used for the last successful step

MRIStepGetLastStep()

Current internal time reached by the solver

MRIStepGetCurrentTime()

Current internal solution reached by the solver

MRIStepGetCurrentState()

Current \(\gamma\) value used by the solver

MRIStepGetCurrentGamma()

Error weight vector for state variables

MRIStepGetErrWeights()

Suggested factor for tolerance scaling

MRIStepGetTolScaleFactor()

Print all statistics

MRIStepPrintAllStats()

Name of constant associated with a return flag

MRIStepGetReturnFlagName()

No. of calls to the \(f^E\) and \(f^I\)

MRIStepGetNumRhsEvals()

No. of failed steps due to a nonlinear solver failure

MRIStepGetNumStepSolveFails()

Current MRI coupling tables

MRIStepGetCurrentCoupling()

Last inner stepper return value

MRIStepGetLastInnerStepFlag()

Retrieve a pointer for user data

MRIStepGetUserData()

int MRIStepGetWorkSpace(void *arkode_mem, long int *lenrw, long int *leniw)

Returns the MRIStep real and integer workspace sizes.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • lenrw – the number of sunrealtype values in the MRIStep workspace.

  • leniw – the number of integer values in the MRIStep workspace.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

int MRIStepGetNumSteps(void *arkode_mem, long int *nssteps, long int *nfsteps)

Returns the cumulative number of slow and fast internal steps taken by the solver (so far).

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • nssteps – number of slow steps taken in the solver.

  • nfsteps – number of fast steps taken in the solver.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

int MRIStepGetLastStep(void *arkode_mem, sunrealtype *hlast)

Returns the integration step size taken on the last successful internal step.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • hlast – step size taken on the last internal step.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

int MRIStepGetCurrentTime(void *arkode_mem, sunrealtype *tcur)

Returns the current internal time reached by the solver.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • tcur – current internal time reached.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

int MRIStepGetCurrentState(void *arkode_mem, N_Vector *ycur)

Returns the current internal solution reached by the solver.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • ycur – current internal solution.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

Notes: Users should exercise extreme caution when using this function, as altering values of ycur may lead to undesirable behavior, depending on the particular use case and on when this routine is called.

int MRIStepGetCurrentGamma(void *arkode_mem, sunrealtype *gamma)

Returns the current internal value of \(\gamma\) used in the implicit solver Newton matrix (see equation (3.32)).

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • gamma – current step size scaling factor in the Newton system.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

int MRIStepGetTolScaleFactor(void *arkode_mem, sunrealtype *tolsfac)

Returns a suggested factor by which the user’s tolerances should be scaled when too much accuracy has been requested for some internal step.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • tolsfac – suggested scaling factor for user-supplied tolerances.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

int MRIStepGetErrWeights(void *arkode_mem, N_Vector eweight)

Returns the current error weight vector.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • eweight – solution error weights at the current time.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

Notes: The user must allocate space for eweight, that will be filled in by this function.

int MRIStepPrintAllStats(void *arkode_mem, FILE *outfile, SUNOutputFormat fmt)

Outputs all of the integrator, nonlinear solver, linear solver, and other statistics.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • outfile – pointer to output file.

  • fmt – the output format:

Return value:
  • ARK_SUCCESS – if the output was successfully.

  • CV_MEM_NULL – if the MRIStep memory was NULL.

  • CV_ILL_INPUT – if an invalid formatting option was provided.

Note

The file scripts/sundials_csv.py provides python utility functions to read and output the data from a SUNDIALS CSV output file using the key and value pair format.

New in version 5.2.0.

char *MRIStepGetReturnFlagName(long int flag)

Returns the name of the MRIStep constant corresponding to flag. See Appendix: ARKODE Constants.

Arguments:

  • flag – a return flag from an MRIStep function.

Return value: The return value is a string containing the name of the corresponding constant.

int MRIStepGetNumRhsEvals(void *arkode_mem, long int *nfse_evals, long int *nfsi_evals)

Returns the number of calls to the user’s outer (slow) right-hand side functions, \(f^E\) and \(f^I\), so far.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • nfse_evals – number of calls to the user’s \(f^E(t,y)\) function.

  • nfsi_evals – number of calls to the user’s \(f^I(t,y)\) function.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

int MRIStepGetNumStepSolveFails(void *arkode_mem, long int *ncnf)

Returns the number of failed steps due to a nonlinear solver failure (so far).

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • ncnf – number of step failures.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

int MRIStepGetCurrentCoupling(void *arkode_mem, MRIStepCoupling *C)

Returns the MRI coupling table currently in use by the solver.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • C – pointer to slow-to-fast MRI coupling structure.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

Notes: The MRIStepCoupling data structure is defined in the header file arkode/arkode_mristep.h. It is defined as a pointer to the following C structure:

struct MRIStepCouplingMem {

   int nmat;        /* number of MRI coupling matrices             */
   int stages;      /* size of coupling matrices (stages * stages) */
   int q;           /* method order of accuracy                    */
   int p;           /* embedding order of accuracy                 */
   sunrealtype ***G;   /* coupling matrices [nmat][stages][stages]    */
   sunrealtype *c;     /* abscissae                                   */

 };
 typedef MRIStepCouplingMem *MRIStepCoupling;

For more details see §3.4.5.3.

int MRIStepGetLastInnerStepFlag(void *arkode_mem, int *flag)

Returns the last return value from the inner stepper.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • flag – inner stepper return value.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

int MRIStepGetUserData(void *arkode_mem, void **user_data)

Returns the user data pointer previously set with MRIStepSetUserData().

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • user_data – memory reference to a user data pointer

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the ARKStep memory was NULL

New in version 5.3.0.

3.4.5.2.9.2. Implicit solver optional output functions

Table 3.19 Implicit solver optional output functions

Optional output

Function name

No. of calls to linear solver setup function

MRIStepGetNumLinSolvSetups()

No. of nonlinear solver iterations

MRIStepGetNumNonlinSolvIters()

No. of nonlinear solver convergence failures

MRIStepGetNumNonlinSolvConvFails()

Single accessor to all nonlinear solver statistics

MRIStepGetNonlinSolvStats()

int MRIStepGetNumLinSolvSetups(void *arkode_mem, long int *nlinsetups)

Returns the number of calls made to the linear solver’s setup routine (so far).

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • nlinsetups – number of linear solver setup calls made.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

Notes: This is only accumulated for the “life” of the nonlinear solver object; the counter is reset whenever a new nonlinear solver module is “attached” to MRIStep, or when MRIStep is resized.

int MRIStepGetNumNonlinSolvIters(void *arkode_mem, long int *nniters)

Returns the number of nonlinear solver iterations performed (so far).

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • nniters – number of nonlinear iterations performed.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

  • ARK_NLS_OP_ERR if the SUNNONLINSOL object returned a failure flag

Notes: This is only accumulated for the “life” of the nonlinear solver object; the counter is reset whenever a new nonlinear solver module is “attached” to MRIStep, or when MRIStep is resized.

int MRIStepGetNumNonlinSolvConvFails(void *arkode_mem, long int *nncfails)

Returns the number of nonlinear solver convergence failures that have occurred (so far).

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • nncfails – number of nonlinear convergence failures.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

Notes: This is only accumulated for the “life” of the nonlinear solver object; the counter is reset whenever a new nonlinear solver module is “attached” to MRIStep, or when MRIStep is resized.

int MRIStepGetNonlinSolvStats(void *arkode_mem, long int *nniters, long int *nncfails)

Returns all of the nonlinear solver statistics in a single call.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • nniters – number of nonlinear iterations performed.

  • nncfails – number of nonlinear convergence failures.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

  • ARK_NLS_OP_ERR if the SUNNONLINSOL object returned a failure flag

Notes: These are only accumulated for the “life” of the nonlinear solver object; the counters are reset whenever a new nonlinear solver module is “attached” to MRIStep, or when MRIStep is resized.

3.4.5.2.9.3. Rootfinding optional output functions

Optional output

Function name

Array showing roots found

MRIStepGetRootInfo()

No. of calls to user root function

MRIStepGetNumGEvals()

int MRIStepGetRootInfo(void *arkode_mem, int *rootsfound)

Returns an array showing which functions were found to have a root.

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • rootsfound – array of length nrtfn with the indices of the user functions \(g_i\) found to have a root (the value of nrtfn was supplied in the call to MRIStepRootInit()). For \(i = 0 \ldots\) nrtfn-1, rootsfound[i] is nonzero if \(g_i\) has a root, and 0 if not.

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

Notes: The user must allocate space for rootsfound prior to calling this function.

For the components of \(g_i\) for which a root was found, the sign of rootsfound[i] indicates the direction of zero-crossing. A value of +1 indicates that \(g_i\) is increasing, while a value of -1 indicates a decreasing \(g_i\).

int MRIStepGetNumGEvals(void *arkode_mem, long int *ngevals)

Returns the cumulative number of calls made to the user’s root function \(g\).

Arguments:
  • arkode_mem – pointer to the MRIStep memory block.

  • ngevals – number of calls made to \(g\) so far.

Return value:
  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

3.4.5.2.9.4. Linear solver interface optional output functions

A variety of optional outputs are available from the ARKLS interface, as listed in the following table and elaborated below. We note that where the name of an output would otherwise conflict with the name of an optional output from the main solver, a suffix LS (for Linear Solver) has been added here (e.g. lenrwLS).

Table 3.20 Linear solver interface optional output functions

Optional output

Function name

Stored Jacobian of the ODE RHS function

MRIStepGetJac()

Time at which the Jacobian was evaluated

MRIStepGetJacTime()

Step number at which the Jacobian was evaluated

MRIStepGetJacNumSteps()

Size of real and integer workspaces

MRIStepGetLinWorkSpace()

No. of Jacobian evaluations

MRIStepGetNumJacEvals()

No. of preconditioner evaluations

MRIStepGetNumPrecEvals()

No. of preconditioner solves

MRIStepGetNumPrecSolves()

No. of linear iterations

MRIStepGetNumLinIters()

No. of linear convergence failures

MRIStepGetNumLinConvFails()

No. of Jacobian-vector setup evaluations

MRIStepGetNumJTSetupEvals()

No. of Jacobian-vector product evaluations

MRIStepGetNumJtimesEvals()

No. of fs calls for finite diff. \(J\) or \(Jv\) evals.

MRIStepGetNumLinRhsEvals()

Last return from a linear solver function

MRIStepGetLastLinFlag()

Name of constant associated with a return flag

MRIStepGetLinReturnFlagName()

int MRIStepGetJac(void *arkode_mem, SUNMatrix *J)

Returns the internally stored copy of the Jacobian matrix of the ODE implicit slow right-hand side function.

Parameters:
  • arkode_mem – the MRIStep memory structure

  • J – the Jacobian matrix

Return values:
  • ARKLS_SUCCESS – the output value has been successfully set

  • ARKLS_MEM_NULLarkode_mem was NULL

  • ARKLS_LMEM_NULL – the linear solver interface has not been initialized

Warning

This function is provided for debugging purposes and the values in the returned matrix should not be altered.

int MRIStepGetJacTime(void *arkode_mem, sunrealtype *t_J)

Returns the time at which the internally stored copy of the Jacobian matrix of the ODE implicit slow right-hand side function was evaluated.

Parameters:
  • arkode_mem – the MRIStep memory structure

  • t_J – the time at which the Jacobian was evaluated

Return values:
  • ARKLS_SUCCESS – the output value has been successfully set

  • ARKLS_MEM_NULLarkode_mem was NULL

  • ARKLS_LMEM_NULL – the linear solver interface has not been initialized

int MRIStepGetJacNumSteps(void *arkode_mem, long int *nst_J)

Returns the value of the internal step counter at which the internally stored copy of the Jacobian matrix of the ODE implicit slow right-hand side function was evaluated.

Parameters:
  • arkode_mem – the MRIStep memory structure

  • nst_J – the value of the internal step counter at which the Jacobian was evaluated

Return values:
  • ARKLS_SUCCESS – the output value has been successfully set

  • ARKLS_MEM_NULLarkode_mem was NULL

  • ARKLS_LMEM_NULL – the linear solver interface has not been initialized

int MRIStepGetLinWorkSpace(void *arkode_mem, long int *lenrwLS, long int *leniwLS)

Returns the real and integer workspace used by the ARKLS linear solver interface.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • lenrwLS – the number of sunrealtype values in the ARKLS workspace.

  • leniwLS – the number of integer values in the ARKLS workspace.

Return value:

  • ARKLS_SUCCESS if successful

  • ARKLS_MEM_NULL if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL if the linear solver memory was NULL

Notes: The workspace requirements reported by this routine correspond only to memory allocated within this interface and to memory allocated by the SUNLinearSolver object attached to it. The template Jacobian matrix allocated by the user outside of ARKLS is not included in this report.

In a parallel setting, the above values are global (i.e., summed over all processors).

int MRIStepGetNumJacEvals(void *arkode_mem, long int *njevals)

Returns the number of Jacobian evaluations.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • njevals – number of Jacobian evaluations.

Return value:

  • ARKLS_SUCCESS if successful

  • ARKLS_MEM_NULL if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL if the linear solver memory was NULL

Notes: This is only accumulated for the “life” of the linear solver object; the counter is reset whenever a new linear solver module is “attached” to MRIStep, or when MRIStep is resized.

int MRIStepGetNumPrecEvals(void *arkode_mem, long int *npevals)

Returns the total number of preconditioner evaluations, i.e., the number of calls made to psetup with jok = SUNFALSE and that returned *jcurPtr = SUNTRUE.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • npevals – the current number of calls to psetup.

Return value:

  • ARKLS_SUCCESS if successful

  • ARKLS_MEM_NULL if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL if the linear solver memory was NULL

Notes: This is only accumulated for the “life” of the linear solver object; the counter is reset whenever a new linear solver module is “attached” to MRIStep, or when MRIStep is resized.

int MRIStepGetNumPrecSolves(void *arkode_mem, long int *npsolves)

Returns the number of calls made to the preconditioner solve function, psolve.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • npsolves – the number of calls to psolve.

Return value:

  • ARKLS_SUCCESS if successful

  • ARKLS_MEM_NULL if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL if the linear solver memory was NULL

Notes: This is only accumulated for the “life” of the linear solver object; the counter is reset whenever a new linear solver module is “attached” to MRIStep, or when MRIStep is resized.

int MRIStepGetNumLinIters(void *arkode_mem, long int *nliters)

Returns the cumulative number of linear iterations.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • nliters – the current number of linear iterations.

Return value:

  • ARKLS_SUCCESS if successful

  • ARKLS_MEM_NULL if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL if the linear solver memory was NULL

Notes: This is only accumulated for the “life” of the linear solver object; the counter is reset whenever a new linear solver module is “attached” to MRIStep, or when MRIStep is resized.

int MRIStepGetNumLinConvFails(void *arkode_mem, long int *nlcfails)

Returns the cumulative number of linear convergence failures.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • nlcfails – the current number of linear convergence failures.

Return value:

  • ARKLS_SUCCESS if successful

  • ARKLS_MEM_NULL if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL if the linear solver memory was NULL

Notes: This is only accumulated for the “life” of the linear solver object; the counter is reset whenever a new linear solver module is “attached” to MRIStep, or when MRIStep is resized.

int MRIStepGetNumJTSetupEvals(void *arkode_mem, long int *njtsetup)

Returns the cumulative number of calls made to the user-supplied Jacobian-vector setup function, jtsetup.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • njtsetup – the current number of calls to jtsetup.

Return value:

  • ARKLS_SUCCESS if successful

  • ARKLS_MEM_NULL if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL if the linear solver memory was NULL

Notes: This is only accumulated for the “life” of the linear solver object; the counter is reset whenever a new linear solver module is “attached” to MRIStep, or when MRIStep is resized.

int MRIStepGetNumJtimesEvals(void *arkode_mem, long int *njvevals)

Returns the cumulative number of calls made to the Jacobian-vector product function, jtimes.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • njvevals – the current number of calls to jtimes.

Return value:

  • ARKLS_SUCCESS if successful

  • ARKLS_MEM_NULL if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL if the linear solver memory was NULL

Notes: This is only accumulated for the “life” of the linear solver object; the counter is reset whenever a new linear solver module is “attached” to MRIStep, or when MRIStep is resized.

int MRIStepGetNumLinRhsEvals(void *arkode_mem, long int *nfevalsLS)

Returns the number of calls to the user-supplied implicit right-hand side function \(f^I\) for finite difference Jacobian or Jacobian-vector product approximation.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • nfevalsLS – the number of calls to the user implicit right-hand side function.

Return value:

  • ARKLS_SUCCESS if successful

  • ARKLS_MEM_NULL if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL if the linear solver memory was NULL

Notes: The value nfevalsLS is incremented only if the default internal difference quotient function is used.

This is only accumulated for the “life” of the linear solver object; the counter is reset whenever a new linear solver module is “attached” to MRIStep, or when MRIStep is resized.

int MRIStepGetLastLinFlag(void *arkode_mem, long int *lsflag)

Returns the last return value from an ARKLS routine.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • lsflag – the value of the last return flag from an ARKLS function.

Return value:

  • ARKLS_SUCCESS if successful

  • ARKLS_MEM_NULL if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL if the linear solver memory was NULL

Notes: If the ARKLS setup function failed when using the SUNLINSOL_DENSE or SUNLINSOL_BAND modules, then the value of lsflag is equal to the column index (numbered from one) at which a zero diagonal element was encountered during the LU factorization of the (dense or banded) Jacobian matrix. For all other failures, lsflag is negative.

Otherwise, if the ARKLS setup function failed (MRIStepEvolve() returned ARK_LSETUP_FAIL), then lsflag will be SUNLS_PSET_FAIL_UNREC, SUNLS_ASET_FAIL_UNREC or SUN_ERR_EXT_FAIL.

If the ARKLS solve function failed (MRIStepEvolve() returned ARK_LSOLVE_FAIL), then lsflag contains the error return flag from the SUNLinearSolver object, which will be one of: SUN_ERR_ARG_CORRUPTRRUPT, indicating that the SUNLinearSolver memory is NULL; SUNLS_ATIMES_NULL, indicating that a matrix-free iterative solver was provided, but is missing a routine for the matrix-vector product approximation, SUNLS_ATIMES_FAIL_UNREC, indicating an unrecoverable failure in the \(Jv\) function; SUNLS_PSOLVE_NULL, indicating that an iterative linear solver was configured to use preconditioning, but no preconditioner solve routine was provided, SUNLS_PSOLVE_FAIL_UNREC, indicating that the preconditioner solve function failed unrecoverably; SUNLS_GS_FAIL, indicating a failure in the Gram-Schmidt procedure (SPGMR and SPFGMR only); SUNLS_QRSOL_FAIL, indicating that the matrix \(R\) was found to be singular during the QR solve phase (SPGMR and SPFGMR only); or SUN_ERR_EXT_FAIL, indicating an unrecoverable failure in an external iterative linear solver package.

char *MRIStepGetLinReturnFlagName(long int lsflag)

Returns the name of the ARKLS constant corresponding to lsflag.

Arguments:

  • lsflag – a return flag from an ARKLS function.

Return value: The return value is a string containing the name of the corresponding constant. If using the SUNLINSOL_DENSE or SUNLINSOL_BAND modules, then if 1 \(\le\) lsflag \(\le n\) (LU factorization failed), this routine returns “NONE”.

3.4.5.2.9.5. General usability functions

The following optional routines may be called by a user to inquire about existing solver parameters or write the current MRI coupling table. While neither of these would typically be called during the course of solving an initial value problem, these may be useful for users wishing to better understand MRIStep.

Table 3.21 General usability functions

Optional routine

Function name

Output all MRIStep solver parameters

MRIStepWriteParameters()

Output the current MRI coupling table

MRIStepWriteCoupling()

int MRIStepWriteParameters(void *arkode_mem, FILE *fp)

Outputs all MRIStep solver parameters to the provided file pointer.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • fp – pointer to use for printing the solver parameters.

Return value:

  • ARKS_SUCCESS if successful

  • ARKS_MEM_NULL if the MRIStep memory was NULL

Notes: The fp argument can be stdout or stderr, or it may point to a specific file created using fopen.

When run in parallel, only one process should set a non-NULL value for this pointer, since parameters for all processes would be identical.

int MRIStepWriteCoupling(void *arkode_mem, FILE *fp)

Outputs the current MRI coupling table to the provided file pointer.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • fp – pointer to use for printing the Butcher tables.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

Notes: The fp argument can be stdout or stderr, or it may point to a specific file created using fopen.

When run in parallel, only one process should set a non-NULL value for this pointer, since tables for all processes would be identical.

3.4.5.2.10. MRIStep re-initialization function

To reinitialize the MRIStep module for the solution of a new problem, where a prior call to MRIStepCreate() has been made, the user must call the function MRIStepReInit(). The new problem must have the same size as the previous one. This routine retains the current settings for all ARKstep module options and performs the same input checking and initializations that are done in MRIStepCreate(), but it performs no memory allocation as is assumes that the existing internal memory is sufficient for the new problem. A call to this re-initialization routine deletes the solution history that was stored internally during the previous integration, and deletes any previously-set tstop value specified via a call to MRIStepSetStopTime(). Following a successful call to MRIStepReInit(), call MRIStepEvolve() again for the solution of the new problem.

The use of MRIStepReInit() requires that the number of Runge–Kutta stages for both the slow and fast methods be no larger for the new problem than for the previous problem.

One important use of the MRIStepReInit() function is in the treating of jump discontinuities in the RHS functions. Except in cases of fairly small jumps, it is usually more efficient to stop at each point of discontinuity and restart the integrator with a readjusted ODE model, using a call to this routine. To stop when the location of the discontinuity is known, simply make that location a value of tout. To stop when the location of the discontinuity is determined by the solution, use the rootfinding feature. In either case, it is critical that the RHS functions not incorporate the discontinuity, but rather have a smooth extension over the discontinuity, so that the step across it (and subsequent rootfinding, if used) can be done efficiently. Then use a switch within the RHS functions (communicated through user_data) that can be flipped between the stopping of the integration and the restart, so that the restarted problem uses the new values (which have jumped). Similar comments apply if there is to be a jump in the dependent variable vector.

int MRIStepReInit(void *arkode_mem, ARKRhsFn fse, ARKRhsFn fsi, sunrealtype t0, N_Vector y0)

Provides required problem specifications and re-initializes the MRIStep outer (slow) stepper.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • fse – the name of the function (of type ARKRhsFn()) defining the explicit slow portion of the right-hand side function in \(\dot{y} = f^E(t,y) + f^I(t,y) + f^F(t,y)\).

  • fsi – the name of the function (of type ARKRhsFn()) defining the implicit slow portion of the right-hand side function in \(\dot{y} = f^E(t,y) + f^I(t,y) + f^F(t,y)\).

  • t0 – the initial value of \(t\).

  • y0 – the initial condition vector \(y(t_0)\).

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

  • ARK_MEM_FAIL if a memory allocation failed

  • ARK_ILL_INPUT if an argument has an illegal value.

Notes: If the inner (fast) stepper also needs to be reinitialized, its reinitialization function should be called before calling MRIStepReInit() to reinitialize the outer stepper.

All previously set options are retained but may be updated by calling the appropriate “Set” functions.

If an error occurred, MRIStepReInit() also sends an error message to the error handler function.

3.4.5.2.11. MRIStep reset function

To reset the MRIStep module to a particular state \((t_R,y(t_R))\) for the continued solution of a problem, where a prior call to MRIStepCreate() has been made, the user must call the function MRIStepReset(). Like MRIStepReInit() this routine retains the current settings for all MRIStep module options and performs no memory allocations but, unlike MRIStepReInit(), this routine performs only a subset of the input checking and initializations that are done in MRIStepCreate(). In particular this routine retains all internal counter values and the step size/error history and does not reinitialize the linear and/or nonlinear solver but it does indicate that a linear solver setup is necessary in the next step. Like MRIStepReInit(), a call to MRIStepReset() will delete any previously-set tstop value specified via a call to MRIStepSetStopTime(). Following a successful call to MRIStepReset(), call MRIStepEvolve() again to continue solving the problem. By default the next call to MRIStepEvolve() will use the step size computed by MRIStep prior to calling MRIStepReset(). To set a different step size or have MRIStep estimate a new step size use MRIStepSetInitStep().

One important use of the MRIStepReset() function is in the treating of jump discontinuities in the RHS functions. Except in cases of fairly small jumps, it is usually more efficient to stop at each point of discontinuity and restart the integrator with a readjusted ODE model, using a call to MRIStepReset(). To stop when the location of the discontinuity is known, simply make that location a value of tout. To stop when the location of the discontinuity is determined by the solution, use the rootfinding feature. In either case, it is critical that the RHS functions not incorporate the discontinuity, but rather have a smooth extension over the discontinuity, so that the step across it (and subsequent rootfinding, if used) can be done efficiently. Then use a switch within the RHS functions (communicated through user_data) that can be flipped between the stopping of the integration and the restart, so that the restarted problem uses the new values (which have jumped). Similar comments apply if there is to be a jump in the dependent variable vector.

int MRIStepReset(void *arkode_mem, sunrealtype tR, N_Vector yR)

Resets the current MRIStep outer (slow) time-stepper module state to the provided independent variable value and dependent variable vector.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • tR – the value of the independent variable \(t\).

  • yR – the value of the dependent variable vector \(y(t_R)\).

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

  • ARK_MEM_FAIL if a memory allocation failed

  • ARK_ILL_INPUT if an argument has an illegal value.

Notes: If the inner (fast) stepper also needs to be reset, its reset function should be called before calling MRIStepReset() to reset the outer stepper.

All previously set options are retained but may be updated by calling the appropriate “Set” functions.

If an error occurred, MRIStepReset() also sends an error message to the error handler function.

Changed in version 5.3.0: This now calls the corresponding MRIStepInnerResetFn with the same (tR, yR) arguments for the MRIStepInnerStepper object that is used to evolve the MRI “fast” time scale subproblems.

3.4.5.2.12. MRIStep system resize function

For simulations involving changes to the number of equations and unknowns in the ODE system (e.g. when using spatially-adaptive PDE simulations under a method-of-lines approach), the MRIStep integrator may be “resized” between slow integration steps, through calls to the MRIStepResize() function. This function modifies MRIStep’s internal memory structures to use the new problem size.

To aid in the vector resize operation, the user can supply a vector resize function that will take as input a vector with the previous size, and transform it in-place to return a corresponding vector of the new size. If this function (of type ARKVecResizeFn()) is not supplied (i.e., is set to NULL), then all existing vectors internal to MRIStep will be destroyed and re-cloned from the new input vector.

int MRIStepResize(void *arkode_mem, N_Vector yR, sunrealtype tR, ARKVecResizeFn resize, void *resize_data)

Re-initializes MRIStep with a different state vector.

Arguments:

  • arkode_mem – pointer to the MRIStep memory block.

  • yR – the newly-sized solution vector, holding the current dependent variable values \(y(t_R)\).

  • tR – the current value of the independent variable \(t_R\) (this must be consistent with yR).

  • resize – the user-supplied vector resize function (of type ARKVecResizeFn().

  • resize_data – the user-supplied data structure to be passed to resize when modifying internal MRIStep vectors.

Return value:

  • ARK_SUCCESS if successful

  • ARK_MEM_NULL if the MRIStep memory was NULL

  • ARK_NO_MALLOC if arkode_mem was not allocated.

  • ARK_ILL_INPUT if an argument has an illegal value.

Notes: If an error occurred, MRIStepResize() also sends an error message to the error handler function.

Resizing the linear solver:

When using any of the SUNDIALS-provided linear solver modules, the linear solver memory structures must also be resized. At present, none of these include a solver-specific “resize” function, so the linear solver memory must be destroyed and re-allocated following each call to MRIStepResize(). Moreover, the existing ARKLS interface should then be deleted and recreated by attaching the updated SUNLinearSolver (and possibly SUNMatrix) object(s) through calls to MRIStepSetLinearSolver().

If any user-supplied routines are provided to aid the linear solver (e.g. Jacobian construction, Jacobian-vector product, mass-matrix-vector product, preconditioning), then the corresponding “set” routines must be called again following the solver re-specification.

Resizing the absolute tolerance array:

If using array-valued absolute tolerances, the absolute tolerance vector will be invalid after the call to MRIStepResize(), so the new absolute tolerance vector should be re-set following each call to MRIStepResize() through a new call to MRIStepSVtolerances() and possibly MRIStepResVtolerance() if applicable.

If scalar-valued tolerances or a tolerance function was specified through either MRIStepSStolerances() or MRIStepWFtolerances(), then these will remain valid and no further action is necessary.

Note

For an example showing usage of the similar ARKStepResize() routine, see the supplied serial C example problem, ark_heat1D_adapt.c.