2.4.11.2. MRIStep User-callable functions

This section describes the MRIStep-specific functions that may be called by the user to setup and then solve an IVP using the MRIStep time-stepping module. The large majority of these routines merely wrap underlying ARKODE functions, and are now deprecated – each of these are clearly marked. However, some of these user-callable functions are specific to MRIStep, as explained below.

As discussed in the main ARKODE user-callable function introduction, each of ARKODE’s time-stepping modules clarifies the categories of user-callable functions that it supports. MRIStep supports the following categories:

  • temporal adaptivity

  • implicit nonlinear and/or linear solvers

MRIStep also has forcing function support when converted to a SUNStepper or MRIStepInnerStepper. See ARKodeCreateSUNStepper() and ARKStepCreateMRIStepInnerStepper() for additional details.

2.4.11.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.

Parameters:
  • 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 §1.3)

Returns:

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) object */
MRIStepInnerStepper stepper = NULL;

/* create an ARKODE object, setting fast (inner) right-hand side
   functions and the initial condition */
inner_arkode_mem = *StepCreate(...);

/* configure the inner integrator */
retval = ARKodeSet*(inner_arkode_mem, ...);

/* create MRIStepInnerStepper wrapper for the ARKODE integrator */
flag = ARKodeCreateMRIStepInnerStepper(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

  • examples/arkode/CXX_serial/ark_test_kpr_nestedmri.cpp (uses MRIStep within itself)

void MRIStepFree(void **arkode_mem)

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

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

Deprecated since version 6.1.0: Use ARKodeFree() instead.

2.4.11.2.2. MRIStep tolerance specification functions

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

This function specifies scalar relative and absolute tolerances.

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

  • reltol – scalar relative tolerance.

  • abstol – scalar absolute tolerance.

Return values:
  • 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 had an illegal value (e.g. a negative tolerance).

Deprecated since version 6.1.0: Use ARKodeSStolerances() instead.

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).

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

  • reltol – scalar relative tolerance.

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

Return values:
  • 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 had an illegal value (e.g. a negative tolerance).

Deprecated since version 6.1.0: Use ARKodeSVtolerances() instead.

int MRIStepWFtolerances(void *arkode_mem, ARKEwtFn efun)

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

Parameters:
  • 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 values:
  • 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

Deprecated since version 6.1.0: Use ARKodeWFtolerances() instead.

2.4.11.2.3. Linear solver interface functions

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).

Parameters:
  • 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 values:
  • 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.

Note

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 §9 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.

Deprecated since version 6.1.0: Use ARKodeSetLinearSolver() instead.

2.4.11.2.4. Nonlinear solver interface functions

int MRIStepSetNonlinearSolver(void *arkode_mem, SUNNonlinearSolver NLS)

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

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

  • NLS – the SUNNonlinearSolver object to use.

Return values:
  • 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.

Note

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

Deprecated since version 6.1.0: Use ARKodeSetNonlinearSolver() instead.

2.4.11.2.5. Rootfinding initialization function

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().

Parameters:
  • 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 values:
  • 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.

Note

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 activated for the fast (inner) integrator.

Deprecated since version 6.1.0: Use ARKodeRootInit() instead.

2.4.11.2.6. MRIStep solver function

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

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

Parameters:
  • 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 §2.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 values:
  • 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 ARKodeGetRootInfo() 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.

Note

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.

Deprecated since version 6.1.0: Use ARKodeEvolve() instead.

2.4.11.2.7. Optional input functions

2.4.11.2.7.1. Optional inputs for MRIStep

int MRIStepSetDefaults(void *arkode_mem)

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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()).

Deprecated since version 6.1.0: Use ARKodeSetDefaults() instead.

int MRIStepSetInterpolantType(void *arkode_mem, int itype)

Deprecated since version 6.1.0: This function is now a wrapper to ARKodeSetInterpolantType(), see the documentation for that function instead.

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).

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

  • degree – requested polynomial degree.

Return values:
  • 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

Note

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.

Deprecated since version 6.1.0: Use ARKodeSetInterpolantDegree() instead.

int MRIStepSetDenseOrder(void *arkode_mem, int dord)

Deprecated since version 5.2.0: Use ARKodeSetInterpolantDegree() 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.

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

  • diagfp – pointer to the diagnostics output file.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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).

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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

Deprecated since version 6.1.0: Use ARKodeSetFixedStep() instead.

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.

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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

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

Deprecated since version 6.1.0: Use ARKodeSetMaxHnilWarns() instead.

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.

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

  • mxsteps – maximum allowed number of internal steps.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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

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

Deprecated since version 6.1.0: Use ARKodeSetMaxNumSteps() instead.

int MRIStepSetStopTime(void *arkode_mem, sunrealtype tstop)

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

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

  • tstop – stopping time for the integrator.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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 re-enabled 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().

Deprecated since version 6.1.0: Use ARKodeSetStopTime() instead.

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\)).

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

Added in version 5.6.0.

Deprecated since version 6.1.0: Use ARKodeSetInterpolateStopTime() instead.

int MRIStepClearStopTime(void *arkode_mem)

Disables the stop time set with MRIStepSetStopTime().

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

Note

The stop time can be re-enabled though a new call to MRIStepSetStopTime().

Added in version 5.5.1.

Deprecated since version 6.1.0: Use ARKodeClearStopTime() instead.

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.

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

  • user_data – pointer to the user data.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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.

Deprecated since version 6.1.0: Use ARKodeSetUserData() instead.

int MRIStepSetPreInnerFn(void *arkode_mem, MRIStepPreInnerFn prefn)

Specifies the function called before each inner integration.

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

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

Return values:
  • 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.

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

2.4.11.2.7.2. Optional inputs for IVP method selection

Table 2.4 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 1 will result in using the default.

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

  • ord – the method order.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

Deprecated since version 6.1.0: Use ARKodeSetOrder() instead.

int MRIStepSetCoupling(void *arkode_mem, MRIStepCoupling C)

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

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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

Warning

This should not be used with ARKodeSetOrder().

2.4.11.2.7.3. Optional inputs for implicit stage solves

int MRIStepSetLinear(void *arkode_mem, int timedepend)

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

Parameters:
  • 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 values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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.

Deprecated since version 6.1.0: Use ARKodeSetLinear() instead.

int MRIStepSetNonlinear(void *arkode_mem)

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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.

Deprecated since version 6.1.0: Use ARKodeSetNonlinear() instead.

int MRIStepSetPredictorMethod(void *arkode_mem, int method)

Specifies the method to use for predicting implicit solutions.

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

  • method

    the predictor method

    • 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 values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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

Warning

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

Deprecated since version 6.1.0: Use ARKodeSetPredictorMethod() instead.

int MRIStepSetMaxNonlinIters(void *arkode_mem, int maxcor)

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

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

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

Return values:
  • 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

Note

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

Deprecated since version 6.1.0: Use ARKodeSetMaxNonlinIters() instead.

int MRIStepSetNonlinConvCoef(void *arkode_mem, sunrealtype nlscoef)

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

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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

Deprecated since version 6.1.0: Use ARKodeSetNonlinConvCoef() instead.

int MRIStepSetNonlinCRDown(void *arkode_mem, sunrealtype crdown)

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

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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

Deprecated since version 6.1.0: Use ARKodeSetNonlinCRDown() instead.

int MRIStepSetNonlinRDiv(void *arkode_mem, sunrealtype rdiv)

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

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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

Deprecated since version 6.1.0: Use ARKodeSetNonlinRDiv() instead.

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.

Parameters:
  • 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 values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

Note

See §2.4.4.6 for more information on this user-supplied routine.

Deprecated since version 6.1.0: Use ARKodeSetStagePredictFn() instead.

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.

Parameters:
  • 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 values:
  • ARK_SUCCESS – if successful.

  • ARK_MEM_NULL – if the MRIStep memory was NULL.

Note

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().

Deprecated since version 6.1.0: Use ARKodeSetNlsRhsFn() instead.

int MRIStepSetDeduceImplicitRhs(void *arkode_mem, sunbooleantype deduce)

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

Parameters:
  • 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 values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

Added in version 5.2.0.

Deprecated since version 6.1.0: Use ARKodeSetDeduceImplicitRhs() instead.

2.4.11.2.7.4. Linear solver interface optional input functions

2.4.11.2.7.4.1. Optional inputs for the ARKLS linear solver interface

int MRIStepSetDeltaGammaMax(void *arkode_mem, sunrealtype dgmax)

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

Parameters:
  • 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 values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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

Deprecated since version 6.1.0: Use ARKodeSetDeltaGammaMax() instead.

int MRIStepSetLSetupFrequency(void *arkode_mem, int msbp)

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

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

  • msbp – the linear solver setup frequency.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

Note

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.

Deprecated since version 6.1.0: Use ARKodeSetLSetupFrequency() instead.

int MRIStepSetJacEvalFrequency(void *arkode_mem, long int msbj)

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

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

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

Return values:
  • ARKLS_SUCCESS – if successful.

  • ARKLS_MEM_NULL – if the MRIStep memory was NULL.

  • ARKLS_LMEM_NULL – if the linear solver memory was NULL.

Note

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().

Deprecated since version 6.1.0: Use ARKodeSetJacEvalFrequency() instead.

2.4.11.2.7.4.2. Optional inputs for matrix-based SUNLinearSolver modules

int MRIStepSetJacFn(void *arkode_mem, ARKLsJacFn jac)

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

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

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

Return values:
  • ARKLS_SUCCESS – if successful

  • ARKLS_MEM_NULL – if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL – if the linear solver memory was NULL

Note

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 §2.4.4.

Deprecated since version 6.1.0: Use ARKodeSetJacFn() instead.

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.

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

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

Return values:
  • ARKLS_SUCCESS – if successful

  • ARKLS_MEM_NULL – if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL – if the linear solver memory was NULL

Note

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 §2.4.4.

Deprecated since version 6.1.0: Use ARKodeSetLinSysFn() instead.

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 §10.2.1.

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

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

Return values:
  • 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

Note

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

Deprecated since version 6.1.0: Use ARKodeSetLinearSolutionScaling() instead.

2.4.11.2.7.4.3. Optional inputs for matrix-free SUNLinearSolver modules

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

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

Parameters:
  • 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 values:
  • 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.

Note

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 §2.4.4.

Deprecated since version 6.1.0: Use ARKodeSetJacTimes() instead.

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.

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

  • jtimesRhsFn – the name of the C function defining the alternative right-hand side function.

Return values:
  • 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.

Note

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().

Deprecated since version 6.1.0: Use ARKodeSetJacTimesRhsFn() instead.

2.4.11.2.7.4.4. Optional inputs for iterative SUNLinearSolver modules

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

Specifies the user-supplied preconditioner setup and solve functions.

Parameters:
  • 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 values:
  • 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.

Note

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 §2.4.4.

Deprecated since version 6.1.0: Use ARKodeSetPreconditioner() instead.

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.

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

  • eplifac – linear convergence safety factor.

Return values:
  • 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.

Note

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().

Deprecated since version 6.1.0: Use ARKodeSetEpsLin() instead.

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.

Parameters:
  • 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 values:
  • ARK_SUCCESS – if successful.

  • ARK_MEM_NULL – if the MRIStep memory was NULL.

Note

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

Deprecated since version 6.1.0: Use ARKodeSetLSNormFactor() instead.

2.4.11.2.7.5. Rootfinding optional input functions

int MRIStepSetRootDirection(void *arkode_mem, int *rootdir)

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

Parameters:
  • 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 values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

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

Deprecated since version 6.1.0: Use ARKodeSetRootDirection() instead.

int MRIStepSetNoInactiveRootWarn(void *arkode_mem)

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory is NULL

Note

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.

Deprecated since version 6.1.0: Use ARKodeSetNoInactiveRootWarn() instead.

2.4.11.2.8. Interpolated output function

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.

Parameters:
  • 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 values:
  • 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

Note

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.

Deprecated since version 6.1.0: Use ARKodeGetDky() instead.

2.4.11.2.9. Optional output functions

2.4.11.2.9.1. Main solver optional output functions

int MRIStepGetNumInnerStepperFails(void *arkode_mem, long int *inner_fails)

Returns the number of recoverable failures reported by the inner stepper (so far).

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

  • inner_fails – number of failed fast (inner) integrations.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Added in version 6.2.0.

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

Returns the MRIStep real and integer workspace sizes.

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

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetWorkSpace() instead.

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).

Parameters:
  • 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 values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetNumSteps() instead.

int MRIStepGetLastStep(void *arkode_mem, sunrealtype *hlast)

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

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetLastStep() instead.

int MRIStepGetCurrentTime(void *arkode_mem, sunrealtype *tcur)

Returns the current internal time reached by the solver.

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

  • tcur – current internal time reached.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetCurrentTime() instead.

int MRIStepGetCurrentState(void *arkode_mem, N_Vector *ycur)

Returns the current internal solution reached by the solver.

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

  • ycur – current internal solution.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Note

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.

Deprecated since version 6.1.0: Use ARKodeGetCurrentState() instead.

int MRIStepGetCurrentGamma(void *arkode_mem, sunrealtype *gamma)

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

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetCurrentGamma() instead.

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.

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetTolScaleFactor() instead.

int MRIStepGetErrWeights(void *arkode_mem, N_Vector eweight)

Returns the current error weight vector.

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

  • eweight – solution error weights at the current time.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Note

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

Deprecated since version 6.1.0: Use ARKodeGetErrWeights() instead.

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

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

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

  • outfile – pointer to output file.

  • fmt

    the output format:

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

  • ARK_MEM_NULL – if the MRIStep memory was NULL.

  • ARK_ILL_INPUT – if an invalid formatting option was provided.

Note

The Python module tools/suntools provides utilities to read and output the data from a SUNDIALS CSV output file using the key and value pair format.

Added in version 5.2.0.

Deprecated since version 6.1.0: Use ARKodePrintAllStats() instead.

char *MRIStepGetReturnFlagName(long int flag)

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

Parameters:
  • flag – a return flag from an MRIStep function.

Returns:

A string containing the name of the corresponding constant.

Deprecated since version 6.1.0: Use ARKodeGetReturnFlagName() instead.

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.

Parameters:
  • 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 values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Deprecated since version 6.2.0: Use ARKodeGetNumRhsEvals() instead.

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

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

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

  • ncnf – number of step failures.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetNumStepSolveFails() instead.

int MRIStepGetCurrentCoupling(void *arkode_mem, MRIStepCoupling *C)

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

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Note

The MRIStepCoupling data structure is defined in the header file arkode/arkode_mristep.h. For more details see §2.4.11.3.

int MRIStepGetLastInnerStepFlag(void *arkode_mem, int *flag)

Returns the last return value from the inner stepper.

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

  • flag – inner stepper return value.

Return values:
  • 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().

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

  • user_data – memory reference to a user data pointer

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the ARKStep memory was NULL

Added in version 5.3.0.

Deprecated since version 6.1.0: Use ARKodeGetUserData() instead.

2.4.11.2.9.2. Implicit solver optional output functions

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

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

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

  • nlinsetups – number of linear solver setup calls made.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Note

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.

Deprecated since version 6.1.0: Use ARKodeGetNumLinSolvSetups() instead.

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

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

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

  • nniters – number of nonlinear iterations performed.

Return values:
  • 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

Note

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.

Deprecated since version 6.1.0: Use ARKodeGetNumNonlinSolvIters() instead.

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

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

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

  • nncfails – number of nonlinear convergence failures.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Note

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.

Deprecated since version 6.1.0: Use ARKodeGetNumNonlinSolvConvFails() instead.

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

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

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

  • nniters – number of nonlinear iterations performed.

  • nncfails – number of nonlinear convergence failures.

Return values:
  • 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

Note

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.

Deprecated since version 6.1.0: Use ARKodeGetNonlinSolvStats() instead.

2.4.11.2.9.3. Rootfinding optional output functions

int MRIStepGetRootInfo(void *arkode_mem, int *rootsfound)

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

Parameters:
  • 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 values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Note

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\).

Deprecated since version 6.1.0: Use ARKodeGetRootInfo() instead.

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

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

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetNumGEvals() instead.

2.4.11.2.9.4. Linear solver interface optional output functions

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.

Deprecated since version 6.1.0: Use ARKodeGetJac() instead.

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

Deprecated since version 6.1.0: Use ARKodeGetJacTime() instead.

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

Deprecated since version 6.1.0: Use ARKodeGetJacNumSteps() instead.

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

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

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

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

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

Return values:
  • ARKLS_SUCCESS – if successful

  • ARKLS_MEM_NULL – if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL – if the linear solver memory was NULL

Note

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).

Deprecated since version 6.1.0: Use ARKodeGetLinWorkSpace() instead.

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

Returns the number of Jacobian evaluations.

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

  • njevals – number of Jacobian evaluations.

Return values:
  • ARKLS_SUCCESS – if successful

  • ARKLS_MEM_NULL – if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL – if the linear solver memory was NULL

Note

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.

Deprecated since version 6.1.0: Use ARKodeGetNumJacEvals() instead.

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.

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

  • npevals – the current number of calls to psetup.

Return values:
  • ARKLS_SUCCESS – if successful

  • ARKLS_MEM_NULL – if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL – if the linear solver memory was NULL

Note

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.

Deprecated since version 6.1.0: Use ARKodeGetNumPrecEvals() instead.

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

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

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

  • npsolves – the number of calls to psolve.

Return values:
  • ARKLS_SUCCESS – if successful

  • ARKLS_MEM_NULL – if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL – if the linear solver memory was NULL

Note

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.

Deprecated since version 6.1.0: Use ARKodeGetNumPrecSolves() instead.

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

Returns the cumulative number of linear iterations.

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

  • nliters – the current number of linear iterations.

Return values:
  • ARKLS_SUCCESS – if successful

  • ARKLS_MEM_NULL – if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL – if the linear solver memory was NULL

Note

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.

Deprecated since version 6.1.0: Use ARKodeGetNumLinIters() instead.

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

Returns the cumulative number of linear convergence failures.

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

  • nlcfails – the current number of linear convergence failures.

Return values:
  • ARKLS_SUCCESS – if successful

  • ARKLS_MEM_NULL – if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL – if the linear solver memory was NULL

Note

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.

Deprecated since version 6.1.0: Use ARKodeGetNumLinConvFails() instead.

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

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

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

  • njtsetup – the current number of calls to jtsetup.

Return values:
  • ARKLS_SUCCESS – if successful

  • ARKLS_MEM_NULL – if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL – if the linear solver memory was NULL

Note

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.

Deprecated since version 6.1.0: Use ARKodeGetNumJTSetupEvals() instead.

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

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

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

  • njvevals – the current number of calls to jtimes.

Return values:
  • ARKLS_SUCCESS – if successful

  • ARKLS_MEM_NULL – if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL – if the linear solver memory was NULL

Note

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.

Deprecated since version 6.1.0: Use ARKodeGetNumJtimesEvals() instead.

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.

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

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

Return values:
  • ARKLS_SUCCESS – if successful

  • ARKLS_MEM_NULL – if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL – if the linear solver memory was NULL

Note

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.

Deprecated since version 6.1.0: Use ARKodeGetNumLinRhsEvals() instead.

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

Returns the last return value from an ARKLS routine.

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

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

Return values:
  • ARKLS_SUCCESS – if successful

  • ARKLS_MEM_NULL – if the MRIStep memory was NULL

  • ARKLS_LMEM_NULL – if the linear solver memory was NULL

Note

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 SUNLS_PACKAGE_FAIL_UNREC.

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:

  • SUNLS_MEM_NULL, 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

  • SUNLS_PACKAGE_FAIL_UNREC, indicating an unrecoverable failure in an external iterative linear solver package.

Deprecated since version 6.1.0: Use ARKodeGetLastLinFlag() instead.

char *MRIStepGetLinReturnFlagName(long int lsflag)

Returns the name of the ARKLS constant corresponding to lsflag.

Parameters:
  • lsflag – a return flag from an ARKLS function.

Returns:

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”.

Deprecated since version 6.1.0: Use ARKodeGetLinReturnFlagName() instead.

2.4.11.2.9.5. General usability functions

int MRIStepWriteParameters(void *arkode_mem, FILE *fp)

Outputs all MRIStep solver parameters to the provided file pointer.

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

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

Return values:
  • ARKS_SUCCESS – if successful

  • ARKS_MEM_NULL – if the MRIStep memory was NULL

Note

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.

Deprecated since version 6.1.0: Use ARKodeWriteParameters() instead.

int MRIStepWriteCoupling(void *arkode_mem, FILE *fp)

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

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the MRIStep memory was NULL

Note

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.

Deprecated since version 6.1.0: Use MRIStepGetCurrentCoupling() and MRIStepCoupling_Write() instead.

2.4.11.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 MRIStep 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.

Parameters:
  • 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 values:
  • 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.

Note

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.

2.4.11.2.11. MRIStep reset function

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.

Parameters:
  • 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 values:
  • 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.

Note

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.

Deprecated since version 6.1.0: Use ARKodeReset() instead.

2.4.11.2.12. MRIStep system resize function

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

Re-initializes MRIStep with a different state vector.

Parameters:
  • 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 values:
  • 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.

Note

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().

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.

Example codes:

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

Deprecated since version 6.1.0: Use ARKodeResize() instead.