3.4.2.1. 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 the ARKStep module. Most of the steps are independent of the NVECTOR, SUNMATRIX, SUNLINSOL and SUNNONLINSOL implementations used. For the steps that are not, refer to §9, §10, §11, and §12 for the specific name of the function to be called or macro to be referenced.

  1. Initialize parallel or multi-threaded environment, if appropriate.

    For example, call MPI_Init to initialize MPI if used, or set num_threads, the number of threads to use within the threaded vector functions, if used.

  2. Create the SUNDIALS simulation context object.

    Call SUNContext_Create() to allocate the SUNContext object.

  3. Set problem dimensions, etc.

    This generally includes the problem size, N, and may include the local vector length Nlocal.

    Note

    The variables N and Nlocal should be of type sunindextype.

  4. 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 realtype array ydata containing the initial values of \(y\) already exists. Otherwise, create a new vector by making a call of the form

    y0 = 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 §9 for details.

  5. Create ARKStep object

    Call arkode_mem = ARKStepCreate(...) to create the ARKStep memory block. ARKStepCreate() returns a void* pointer to this memory structure. See §3.4.2.2.1 for details.

  6. Specify integration tolerances

    Call ARKStepSStolerances() or ARKStepSVtolerances() to specify either a scalar relative tolerance and scalar absolute tolerance, or a scalar relative tolerance and a vector of absolute tolerances, respectively. Alternatively, call ARKStepWFtolerances() to specify a function which sets directly the weights used in evaluating WRMS vector norms. See §3.4.2.2.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 ARKStepResStolerance(), ARKStepResVtolerance(), or ARKStepResFtolerance().

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

  8. 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 §11.

  9. 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 §11 for details.

  10. 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 §3.4.2.2.3):

    ier = ARKStepSetLinearSolver(...);
    

    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 §3.4.2.2.3):

    ier = ARKStepSetMassLinearSolver(...);
    
  11. 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 §3.4.2.2.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 §12 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 §12.

  12. Attach nonlinear solver module

    If a nonlinear solver object was created above, then it must be attached to ARKStep using the call (for details see §3.4.2.2.5):

    ier = ARKStepSetNonlinearSolver(...);
    
  13. 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 ARKStep, otherwise the optional inputs will be overridden by ARKStep defaults. See §12 for more information on optional inputs.

  14. Set optional inputs

    Call ARKStepSet* functions to change any optional inputs that control the behavior of ARKStep from their default values. See §3.4.2.2.8 for details.

  15. Specify rootfinding problem

    Optionally, call ARKStepRootInit() to initialize a rootfinding problem to be solved during the integration of the ODE system. See §3.4.2.2.6 for general details, and §3.4.2.2.8 for relevant optional input calls.

  16. Advance solution in time

    For each point at which output is desired, call

    ier = ARKStepEvolve(arkode_mem, tout, yout, &tret, itask);
    

    Here, itask specifies the return mode. The vector yout (which can be the same as the vector y0 above) will contain \(y(t_\text{out})\). See §3.4.2.2.7 for details.

  17. Get optional outputs

    Call ARKStepGet* functions to obtain optional output. See §3.4.2.2.10 for details.

  18. Deallocate memory for solution vector

    Upon completion of the integration, deallocate memory for the vector y (or yout) by calling the destructor function:

    N_VDestroy(y);
    
  19. Free solver memory

    Call ARKStepFree() to free the memory allocated for the ARKStep module (and any nonlinear solver module).

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

  21. Free nonlinear solver memory

    If a user-supplied SUNNonlinearSolver was provided to ARKStep, then call SUNNonlinSolFree() to free any memory allocated for the nonlinear solver object created above.

  22. Finalize MPI, if used

    Call MPI_Finalize to terminate MPI.