2.4.9.1. SPRKStep User-callable functions

This section describes the SPRKStep-specific functions that may be called by the user to setup and then solve an IVP using the SPRKStep 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 ERKStep, 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. SPRKStep supports only the basic set of user-callable functions, and does not support any of the restricted groups (time adaptivity, implicit solvers, etc.).

2.4.9.1.1. SPRKStep initialization and deallocation functions

void *SPRKStepCreate(ARKRhsFn f1, ARKRhsFn f2, sunrealtype t0, N_Vector y0, SUNContext sunctx)

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

Parameters:
  • f1 – the name of the C function (of type ARKRhsFn()) defining \(f_1(t,q) = \frac{\partial V(t,q)}{\partial q}\)

  • f2 – the name of the C function (of type ARKRhsFn()) defining \(f_2(t,p) = \frac{\partial T(t,p)}{\partial p}\)

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

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

  • sunctx – the SUNContext object (see §1.4)

Returns:

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

Warning

SPRKStep requires a partitioned problem where f1 should only modify the q variables and f2 should only modify the p variables (or vice versa). However, the vector passed to these functions is the full vector with both p and q. The ordering of the variables is determined implicitly by the user when they set the initial conditions.

void SPRKStepFree(void **arkode_mem)

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

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

Deprecated since version 6.1.0: Use ARKodeFree() instead.

2.4.9.1.2. Rootfinding initialization function

int SPRKStepRootInit(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 SPRKStepCreate(), and before SPRKStepEvolve().

To disable the rootfinding feature after it has already been initialized, or to free memory associated with SPRKStep’s rootfinding module, call SPRKStepRootInit() with nrtfn = 0.

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

Parameters:
  • arkode_mem – pointer to the SPRKStep 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 SPRKStep 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.

Deprecated since version 6.1.0: Use ARKodeRootInit() instead.

2.4.9.1.3. SPRKStep solver function

int SPRKStepEvolve(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 SPRKStep 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 (using one of the dense output routines 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 SPRKStepEvolve() succeeded, and found one or more roots. If the number of root functions, nrtfn, is greater than 1, call SPRKStepGetRootInfo() to see which \(g_i\) were found to have a root at (*tret).

  • ARK_TSTOP_RETURN – if SPRKStepEvolve() 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 SPRKStepEvolve() 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 are 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_ERR_FAILURE – if error test failures occurred either too many times (ark_maxnef) during one internal time step or occurred with \(|h| = h_{min}\).

  • ARK_VECTOROP_ERR – a vector operation error occurred.

Note

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

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

Since interpolation may reduce the accuracy in the reported solution, if full method accuracy is desired the user should issue a call to SPRKStepSetStopTime() before the call to SPRKStepEvolve() to specify a fixed stop time to end the time step and return to the user. Upon return from SPRKStepEvolve(), 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 SPRKStepSetStopTime()). Interpolated outputs may or may not conserve the Hamiltonian. Our testing has shown that Lagrange interpolation typically performs well in this regard, while Hermite interpolation does not. As such, SPRKStep uses the Lagrange interpolation module by default.

On any error return in which one or more internal steps were taken by SPRKStepEvolve(), 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.9.1.4. Optional input functions

2.4.9.1.4.1. Optional inputs for SPRKStep

int SPRKStepSetDefaults(void *arkode_mem)

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

Does not change problem-defining function pointer f or the user_data pointer.

Also leaves alone any data structures or options related to root-finding (those can be reset using SPRKStepRootInit()).

Deprecated since version 6.1.0: Use ARKodeSetDefaults() instead.

int SPRKStepSetInterpolantType(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 SPRKStepSetInterpolantDegree(void *arkode_mem, int degree)

Specifies the degree of the polynomial interpolant used for dense output (i.e. interpolation of solution output values). Allowed values are between 0 and 5.

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

  • degree – requested polynomial degree.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory or interpolation module are NULL

  • ARK_INTERP_FAIL – if this is called after SPRKStepEvolve()

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

Note

This routine should be called after SPRKStepCreate() and before SPRKStepEvolve(). After the first call to SPRKStepEvolve() the interpolation degree may not be changed without first calling SPRKStepReInit().

If a user calls both this routine and SPRKStepSetInterpolantType(), then SPRKStepSetInterpolantType() 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 SPRKStep 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.

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 SPRKStepSetFixedStep(void *arkode_mem, sunrealtype hfixed)

Sets the time step size used within SPRKStep.

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

  • hfixed – value of the fixed step size to use.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Deprecated since version 6.1.0: Use ARKodeSetFixedStep() instead.

int SPRKStepSetMaxNumSteps(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 SPRKStep will return with an error.

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

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

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

  • mxsteps – maximum allowed number of internal steps.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Deprecated since version 6.1.0: Use ARKodeSetMaxNumSteps() instead.

int SPRKStepSetStopTime(void *arkode_mem, sunrealtype tstop)

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

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

A stop time not reached before a call to SPRKStepReInit() or SPRKStepReset() will remain active but can be disabled by calling SPRKStepClearStopTime().

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

  • tstop – stopping time for the integrator.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Deprecated since version 6.1.0: Use ARKodeSetStopTime() instead.

int SPRKStepClearStopTime(void *arkode_mem)

Disables the stop time set with SPRKStepSetStopTime().

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory is NULL

Deprecated since version 6.1.0: Use ARKodeClearStopTime() instead.

int SPRKStepSetUserData(void *arkode_mem, void *user_data)

Specifies the user data block user_data and attaches it to the main SPRKStep memory block.

If specified, the pointer to user_data is passed to all user-supplied functions for which it is an argument; otherwise NULL is passed.

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

  • user_data – pointer to the user data.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Deprecated since version 6.1.0: Use ARKodeSetUserData() instead.

2.4.9.1.4.2. Optional inputs for IVP method selection

Table 2.4 Optional inputs for IVP method selection

Optional input

Function name

Default

Set integrator method order

SPRKStepSetOrder()

4

Set SPRK method

SPRKStepSetMethod()

ARKODE_SPRK_MCLACHLAN_4_4

Set SPRK method by name

SPRKStepSetMethodName()

“ARKODE_SPRK_MCLACHLAN_4_4”

Use compensated summation

SPRKStepSetUseCompensatedSums()

false

int SPRKStepSetOrder(void *arkode_mem, int ord)

Specifies the order of accuracy for the SPRK integration method.

The allowed values are \(1,2,3,4,5,6,8,10\). Any illegal input will result in the default value of 4.

Since ord affects the memory requirements for the internal SPRKStep memory block, it cannot be changed after the first call to SPRKStepEvolve(), unless SPRKStepReInit() is called.

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

  • ord – requested order of accuracy.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Warning

This overrides any previously set method so it should not be used with SPRKStepSetMethod() or SPRKStepSetMethodName().

Deprecated since version 6.1.0: Use ARKodeSetOrder() instead.

int SPRKStepSetMethod(void *arkode_mem, ARKodeSPRKTable sprk_table)

Specifies the SPRK method.

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

  • sprk_table – the SPRK method table.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Note

No error checking is performed on the coefficients contained in the table to ensure its declared order of accuracy.

Warning

This should not be used with ARKodeSetOrder().

int SPRKStepSetMethodName(void *arkode_mem, const char *method)

Specifies the SPRK method by its name.

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

  • method – the SPRK method name.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Warning

This should not be used with ARKodeSetOrder().

int SPRKStepSetUseCompensatedSums(void *arkode_mem, sunbooleantype onoff)

Specifies if compensated summation (and the incremental form) should be used where applicable.

This increases the computational cost by 2 extra vector operations per stage and an additional 5 per time step. It also requires one extra vector to be stored. However, it is signficantly more robust to roundoff error accumulation.

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

  • onoff – should compensated summation be used (1) or not (0)

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

2.4.9.1.4.3. Rootfinding optional input functions

int SPRKStepSetRootDirection(void *arkode_mem, int *rootdir)

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

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

Parameters:
  • arkode_mem – pointer to the SPRKStep 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 SPRKStepRootInit()). 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 SPRKStep memory is NULL

  • ARK_ILL_INPUT – if an argument has an illegal value

Deprecated since version 6.1.0: Use ARKodeSetRootDirection() instead.

int SPRKStepSetNoInactiveRootWarn(void *arkode_mem)

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

SPRKStep 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), SPRKStep will issue a warning which can be disabled with this optional input function.

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory is NULL

Deprecated since version 6.1.0: Use ARKodeSetNoInactiveRootWarn() instead.

2.4.9.1.5. Interpolated output function

int SPRKStepGetDky(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. A user may access the values \(t_n\) and \(h_n\) via the functions SPRKStepGetCurrentTime() and SPRKStepGetLastStep(), respectively.

This routine uses an interpolating polynomial of degree min(degree, 5), where degree is the argument provided to SPRKStepSetInterpolantDegree(). 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 SPRKStep 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 SPRKStep memory is NULL

Note

Dense outputs may or may not conserve the Hamiltonian. Our testing has shown that Lagrange interpolation typically performs well in this regard, while Hermite interpolation does not.

Warning

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

Deprecated since version 6.1.0: Use ARKodeGetDky() instead.

2.4.9.1.6. Optional output functions

2.4.9.1.6.1. Main solver optional output functions

int SPRKStepGetNumSteps(void *arkode_mem, long int *nsteps)

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

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

  • nsteps – number of steps taken in the solver.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetNumSteps() instead.

int SPRKStepGetLastStep(void *arkode_mem, sunrealtype *hlast)

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

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetLastStep() instead.

int SPRKStepGetCurrentStep(void *arkode_mem, sunrealtype *hcur)

Returns the integration step size to be attempted on the next internal step.

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

  • hcur – step size to be attempted on the next internal step.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetCurrentStep() instead.

int SPRKStepGetCurrentTime(void *arkode_mem, sunrealtype *tcur)

Returns the current internal time reached by the solver.

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

  • tcur – current internal time reached.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetCurrentTime() instead.

int SPRKStepGetCurrentState(void *arkode_mem, N_Vector *ycur)

Returns the current internal solution reached by the solver.

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

  • ycur – current internal solution

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory was NULL

Warning

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 SPRKStepGetStepStats(void *arkode_mem, long int *nsteps, sunrealtype *hinused, sunrealtype *hlast, sunrealtype *hcur, sunrealtype *tcur)

Returns many of the most useful optional outputs in a single call.

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

  • nsteps – number of steps taken in the solver.

  • hinused – actual value of initial step size.

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

  • hcur – step size to be attempted on the next internal step.

  • tcur – current internal time reached.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetStepStats() instead.

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

Outputs all of the integrator statistics.

Parameters:
  • arkode_mem – pointer to the SPRKStep 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 SPRKStep memory was NULL.

  • ARK_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.

Deprecated since version 6.1.0: Use ARKodePrintAllStats() instead.

char *SPRKStepGetReturnFlagName(long int flag)

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

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

Returns:

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

Deprecated since version 6.1.0: Use ARKodeGetReturnFlagName() instead.

int SPRKStepGetNumStepAttempts(void *arkode_mem, long int *step_attempts)

Returns the cumulative number of steps attempted by the solver (so far).

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

  • step_attempts – number of steps attempted by solver.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetNumStepAttempts() instead.

int SPRKStepGetNumRhsEvals(void *arkode_mem, long int *nf1, long int *nf2)

Returns the number of calls to the user’s right-hand side functions, \(f_1\) and \(f_2\) (so far).

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

  • nf1 – number of calls to the user’s \(f_1(t,p)\) function.

  • nf2 – number of calls to the user’s \(f_2(t,q)\) function.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory was NULL

int SPRKStepGetCurrentMethod(void *arkode_mem, ARKodeSPRKTable *sprk_table)

Returns the SPRK method coefficient table currently in use by the solver.

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

  • sprk_table – pointer to the SPRK method table.

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory was NULL

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

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

Parameters:
  • arkode_mem – pointer to the SPRKStep 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

Deprecated since version 6.1.0: Use ARKodeGetUserData() instead.

2.4.9.1.6.2. Rootfinding optional output functions

int SPRKStepGetRootInfo(void *arkode_mem, int *rootsfound)

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

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

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

Parameters:
  • arkode_mem – pointer to the SPRKStep 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 SPRKStepRootInit()). 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 SPRKStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetRootInfo() instead.

int SPRKStepGetNumGEvals(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 SPRKStep memory block.

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeGetNumGEvals() instead.

2.4.9.1.6.3. General usability functions

int SPRKStepWriteParameters(void *arkode_mem, FILE *fp)

Outputs all SPRKStep solver parameters to the provided file pointer.

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.

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

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

Return values:
  • ARK_SUCCESS – if successful

  • ARK_MEM_NULL – if the SPRKStep memory was NULL

Deprecated since version 6.1.0: Use ARKodeWriteParameters() instead.

2.4.9.1.7. SPRKStep re-initialization function

To reinitialize the SPRKStep module for the solution of a new problem, where a prior call to SPRKStepCreate() has been made, the user must call the function SPRKStepReInit(). The new problem must have the same size as the previous one. This routine retains the current settings for all SPRKStep module options and performs the same input checking and initializations that are done in SPRKStepCreate(), but it performs no memory allocation as it 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 SPRKStepSetStopTime(). Following a successful call to SPRKStepReInit(), call SPRKStepEvolve() again for the solution of the new problem.

The use of SPRKStepReInit() requires that the number of Runge–Kutta stages, denoted by s, be no larger for the new problem than for the previous problem. This condition is automatically fulfilled if the method order q is left unchanged.

One potential use of the SPRKStepReInit() function is in the treating of jump discontinuities in the RHS function [125]. In lieu of including if statements within the RHS function to handle discontinuities, it may be more computationally efficient to stop at each point of discontinuity (e.g., through use of tstop or the rootfinding feature) and restart the integrator with a readjusted ODE model, using a call to this routine. We note that for the solution to retain temporal accuracy, the RHS function should not incorporate the discontinuity.

int SPRKStepReInit(void *arkode_mem, ARKRhsFn f1, ARKRhsFn f2, sunrealtype t0, N_Vector y0)

Provides required problem specifications and re-initializes the SPRKStep time-stepper module.

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

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

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

  • f1 – the name of the C function (of type ARKRhsFn()) defining \(f1(t,q) = \frac{\partial V(t,q)}{\partial q}\)

  • f2 – the name of the C function (of type ARKRhsFn()) defining \(f2(t,p) = \frac{\partial T(t,p)}{\partial p}\)

  • 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 SPRKStep memory was NULL

  • ARK_MEM_FAIL – if a memory allocation failed

  • ARK_ILL_INPUT – if an argument has an illegal value.

2.4.9.1.8. SPRKStep reset function

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

Resets the current SPRKStep time-stepper module state to the provided independent variable value and dependent variable vector.

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

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

Parameters:
  • arkode_mem – pointer to the SPRKStep 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 SPRKStep memory was NULL

  • ARK_MEM_FAIL – if a memory allocation failed

  • ARK_ILL_INPUTL – if an argument has an illegal value.

Note

By default the next call to SPRKStepEvolve() will use the step size computed by SPRKStep prior to calling SPRKStepReset().

Deprecated since version 6.1.0: Use ARKodeReset() instead.