2.4.2. A skeleton of the user’s main program
The following is a skeleton of the user’s main program (or calling program) for the integration of an ODE IVP using ARKODE. Most of the steps are independent of the NVECTOR, SUNMATRIX, SUNLINSOL and SUNNONLINSOL implementations used. For the steps that are not, refer to §8, §9, §10, and §11 for the specific name of the function to be called or macro to be referenced.
Initialize parallel or multi-threaded environment, if appropriate.
For example, call
MPI_Init
to initialize MPI if used, or setnum_threads
, the number of threads to use within the threaded vector functions, if used.Create the SUNDIALS simulation context object.
Call
SUNContext_Create()
to allocate theSUNContext
object.Set problem dimensions, etc.
This generally includes the problem size,
N
, and may include the local vector lengthNlocal
.Note
The variables
N
andNlocal
should be of typesunindextype
.Set vector of initial values
To set the vector
y0
of initial values, use the appropriate functions defined by the particular NVECTOR implementation.For native SUNDIALS vector implementations (except the CUDA and RAJA based ones), use a call of the form
y0 = N_VMake_***(..., ydata);
if the
sunrealtype
arrayydata
containing the initial values of \(y\) already exists. Otherwise, create a new vector by making a call of the formy0 = N_VNew_***(...);
and then set its elements by accessing the underlying data where it is located with a call of the form
ydata = N_VGetArrayPointer_***(y0);
For details on each of SUNDIALS’ provided vector implementations, see the corresponding sections in §8 for details.
Create ARKODE object
Call a stepper-specific constructor,
arkode_mem = *StepCreate(...)
, to create the ARKODE memory block. These routines return avoid*
pointer to this memory structure. See §2.4.7.1.1, §2.4.8.1.1, §2.4.11.2.1, or §2.4.13.1.1 for details.Specify integration tolerances
Call
ARKodeSStolerances()
orARKodeSVtolerances()
to specify either a scalar relative tolerance and scalar absolute tolerance, or a scalar relative tolerance and a vector of absolute tolerances, respectively. Alternatively, callARKodeWFtolerances()
to specify a function which sets directly the weights used in evaluating WRMS vector norms. See §2.4.3.2 for details.If a problem with non-identity mass matrix is used, and the solution units differ considerably from the equation units, absolute tolerances for the equation residuals (nonlinear and linear) may be specified separately through calls to
ARKodeResStolerance()
,ARKodeResVtolerance()
, orARKodeResFtolerance()
.Create matrix object
If a nonlinear solver requiring a linear solver will be used (e.g., a Newton iteration) and the linear solver will be a matrix-based linear solver, then a template Jacobian matrix must be created by using the appropriate functions defined by the particular SUNMATRIX implementation.
For the SUNDIALS-supplied SUNMATRIX implementations, the matrix object may be created using a call of the form
SUNMatrix A = SUNBandMatrix(..., sunctx);
or similar for the other matrix modules (see §9 for further information).
Similarly, if the problem involves a non-identity mass matrix, and the mass-matrix linear systems will be solved using a direct linear solver, then a template mass matrix must be created by using the appropriate functions defined by the particular SUNMATRIX implementation.
Create linear solver object
If a nonlinear solver requiring a linear solver will be used (e.g., a Newton iteration), or if the problem involves a non-identity mass matrix, then the desired linear solver object(s) must be created by using the appropriate functions defined by the particular SUNLINSOL implementation.
For any of the SUNDIALS-supplied SUNLINSOL implementations, the linear solver object may be created using a call of the form
SUNLinearSolver LS = SUNLinSol_*(...);
where
*
can be replaced with “Dense”, “SPGMR”, or other options, as discussed in §10.Set linear solver optional inputs
Call
*Set*
functions from the selected linear solver module to change optional inputs specific to that linear solver. See the documentation for each SUNLINSOL module in §10 for details.Attach linear solver module
If a linear solver was created above for implicit stage solves, initialize the ARKLS linear solver interface by attaching the linear solver object (and Jacobian matrix object, if applicable) with the call (for details see §2.4.3.3):
ier = ARKodeSetLinearSolver(...);
Similarly, if the problem involves a non-identity mass matrix, initialize the ARKLS mass matrix linear solver interface by attaching the mass linear solver object (and mass matrix object, if applicable) with the call (for details see §2.4.3.3):
ier = ARKodeSetMassLinearSolver(...);
Create nonlinear solver object
If the problem involves an implicit component, and if a non-default nonlinear solver object will be used for implicit stage solves (see §2.4.3.5), then the desired nonlinear solver object must be created by using the appropriate functions defined by the particular SUNNONLINSOL implementation (e.g.,
NLS = SUNNonlinSol_***(...);
where***
is the name of the nonlinear solver (see §11 for details).For the SUNDIALS-supplied SUNNONLINSOL implementations, the nonlinear solver object may be created using a call of the form
SUNNonlinearSolver NLS = SUNNonlinSol_*(...);
where
*
can be replaced with “Newton”, “FixedPoint”, or other options, as discussed in §11.Attach nonlinear solver module
If a nonlinear solver object was created above, then it must be attached to ARKODE using the call (for details see §2.4.3.5):
ier = ARKodeSetNonlinearSolver(...);
Set nonlinear solver optional inputs
Call the appropriate set functions for the selected nonlinear solver module to change optional inputs specific to that nonlinear solver. These must be called after attaching the nonlinear solver to ARKODE, otherwise the optional inputs will be overridden by ARKODE defaults. See §11 for more information on optional inputs.
Set optional inputs
Call
ARKodeSet*
functions to change any optional inputs that control the behavior of ARKODE from their default values. See §2.4.3.8 for details.Additionally, call
*StepSet*
routines to change any stepper-specific optional inputs from their default values. See §2.4.7.1.8, §2.4.8.1.5, §2.4.11.2.7, or §2.4.13.1.4 for details.Specify rootfinding problem
Optionally, call
ARKodeRootInit()
to initialize a rootfinding problem to be solved during the integration of the ODE system. See §2.4.3.6 for general details, and §2.4.3.8 for relevant optional input calls.Advance solution in time
For each point at which output is desired, call
ier = ARKodeEvolve(arkode_mem, tout, yout, &tret, itask);
Here,
itask
specifies the return mode. The vectoryout
(which can be the same as the vectory0
above) will contain \(y(t_\text{out})\). See §2.4.3.7 for details.Get optional outputs
Call
ARKodeGet*
functions to obtain optional output. See §2.4.3.10 for details.Additionally, call
*StepGet*
routines to retrieve any stepper-specific optional outputs. See §2.4.7.1.10, §2.4.8.1.7, §2.4.11.2.9, or §2.4.13.1.6 for details.Deallocate memory for solution vector
Upon completion of the integration, deallocate memory for the vector
y
(oryout
) by calling the destructor function:N_VDestroy(y);
Free solver memory
Call
ARKodeFree()
to free the memory allocated for the ARKODE module (and any nonlinear solver module).Free linear solver and matrix memory
Call
SUNLinSolFree()
and (possibly)SUNMatDestroy()
to free any memory allocated for the linear solver and matrix objects created above.Free nonlinear solver memory
If a user-supplied
SUNNonlinearSolver
was provided to ARKODE, then callSUNNonlinSolFree()
to free any memory allocated for the nonlinear solver object created above.Free the SUNContext object
Call
SUNContext_Free()
to free the memory allocated for theSUNContext
object.Finalize MPI, if used
Call
MPI_Finalize
to terminate MPI.