2.1. Introduction

The ARKODE infrastructure provides adaptive-step time integration modules for stiff, nonstiff and mixed stiff/nonstiff systems of ordinary differential equations (ODEs). ARKODE itself is structured to support a wide range of one-step (but multi-stage) methods, allowing for rapid development of parallel implementations of state-of-the-art time integration methods. At present, ARKODE is packaged with four time-stepping modules, ARKStep, ERKStep, SPRKStep, and MRIStep.

ARKStep supports ODE systems posed in split, linearly-implicit form,

(2.1)\[M(t)\, \dot{y} = f^E(t,y) + f^I(t,y), \qquad y(t_0) = y_0,\]

where \(t\) is the independent variable, \(y\) is the set of dependent variables (in \(\mathbb{R}^N\)), \(M\) is a user-specified, nonsingular operator from \(\mathbb{R}^N\) to \(\mathbb{R}^N\), and the right-hand side function is partitioned into up to two components:

  • \(f^E(t,y)\) contains the “nonstiff” time scale components to be integrated explicitly, and

  • \(f^I(t,y)\) contains the “stiff” time scale components to be integrated implicitly.

Either of these operators may be disabled, allowing for fully explicit, fully implicit, or combination implicit-explicit (ImEx) time integration.

The algorithms used in ARKStep are adaptive- and fixed-step additive Runge–Kutta methods. Such methods are defined through combining two complementary Runge–Kutta methods: one explicit (ERK) and the other diagonally implicit (DIRK). Through appropriately partitioning the ODE right-hand side into explicit and implicit components (2.1), such methods have the potential to enable accurate and efficient time integration of stiff, nonstiff, and mixed stiff/nonstiff systems of ordinary differential equations. A key feature allowing for high efficiency of these methods is that only the components in \(f^I(t,y)\) must be solved implicitly, allowing for splittings tuned for use with optimal implicit solver algorithms.

This framework allows for significant freedom over the constitutive methods used for each component, and ARKODE is packaged with a wide array of built-in methods for use. These built-in Butcher tables include adaptive explicit methods of orders 2-9, adaptive implicit methods of orders 2-5, and adaptive ImEx methods of orders 2-5.

ERKStep focuses specifically on problems posed in explicit form,

(2.2)\[\dot{y} = f(t,y), \qquad y(t_0) = y_0.\]

allowing for increased computational efficiency and memory savings. The algorithms used in ERKStep are adaptive- and fixed-step explicit Runge–Kutta methods. As with ARKStep, the ERKStep module is packaged with adaptive explicit methods of orders 2-9.

SPRKStep focuses on Hamiltonian systems posed in the form,

\[H(t, p, q) = T(t, p) + V(t, q)\]
(2.3)\[\dot{p} = f_1(t,q) = \frac{\partial V(t,q)}{\partial q}, \quad \dot{q} = f_2(t,p) = \frac{\partial T(t,p)}{\partial p},\]

allowing for conservation of quadratic invariants.

MRIStep focuses specifically on problems posed in additive form,

(2.4)\[\dot{y} = f^E(t,y) + f^I(t,y) + f^F(t,y), \qquad y(t_0) = y_0.\]

where here the right-hand side function is additively split into three components:

  • \(f^E(t,y)\) contains the “slow-nonstiff” components of the system (this will be integrated using an explicit method and a large time step \(h^S\)),

  • \(f^I(t,y)\) contains the “slow-stiff” components of the system (this will be integrated using an implicit method and a large time step \(h^S\)), and

  • \(f^F(t,y)\) contains the “fast” components of the system (this will be integrated using a possibly different method than the slow time scale and a small time step \(h^F \ll h^S\)).

For such problems, MRIStep provides fixed-step slow step multirate infinitesimal step (MIS), multirate infinitesimal GARK (MRI-GARK), and implicit-explicit MRI-GARK (IMEX-MRI-GARK) methods, allowing for evolution of the problem (2.4) using multirate methods having orders of accuracy 2-4.

For ARKStep or MRIStep problems that include nonzero implicit term \(f^I(t,y)\), the resulting implicit system (assumed nonlinear, unless specified otherwise) is solved approximately at each integration step, using a SUNNonlinearSolver module, supplied either by the user or from the underlying SUNDIALS infrastructure. For nonlinear solver algorithms that internally require a linear solver, ARKODE may use a variety of SUNLinearSolver modules provided with SUNDIALS, or again may utilize a user-supplied module.

2.1.1. Changes to SUNDIALS in release 6.7.0

New Features and Enhancements

The default number of stages for the SSP Runge-Kutta methods ARKODE_LSRK_SSP_S_2 and ARKODE_LSRK_SSP_S_3 in LSRKStep were changed from 10 and 9, respectively, to their minimum allowable values of 2 and 4. Users may revert to the previous values by calling LSRKStepSetNumSSPStages().

Added the optional function ARKodeInit() to ARKODE to enable data allocation before the first call to ARKodeEvolve() (but after all other optional input routines have been called), to support users who measure memory usage before beginning a simulation.

Added the function ARKodeGetStageIndex() that returns the index of the stage currently being processed, and the total number of stages in the method, for users who wish to compute auxiliary quantities in their IVP right-hand side functions during some stages and not others (e.g., in all but the first or last stage).

Added the functions ARKodeGetLastTime() and ARKodeGetLastState() to return the last successful time and state achieved by ARKODE, respectively.

ARKODE now allows users to supply functions that will be called before each internal time step attempt (ARKodeSetPreStepFn()), after each successful time step (ARKodeSetPostStepFn()), before right-hand side routines are called on an updated state (ARKodeSetPreRhsFn()), and/or once each internal step/stage is computed (ARKodeSetPostprocessStepFn()/ ARKodeSetPostprocessStageFn()). These are considered advanced functions, as they should treat the state vector as read-only, otherwise all theoretical guarantees of solution accuracy and stability will be lost. As a result of these new functions, the values of multiple ARKODE return codes (e.g., ARK_INTERP_FAIL) have been updated; users who key off of the named constants will not be affected, but users who rely on the values themselves should update their codes accordingly.

Note to users utilizing the previously undocumented ARKodeSetPostprocessStepFn() function, the supplied function is now called on the newly computed state vector for all step attempts not just successful steps. To obtain the previous behavior of only calling a function on successful steps, switch to using ARKodeSetPostStepFn().

Added SUNLogger_Set{Error,Warning,Info,Debug}File functions to allow setting logger output streams with a FILE*.

Updated the Kokkos N_Vector to support Kokkos 5.x versions.

Bug Fixes

Fixed a CMake bug where the SuperLU_MT interface would not be built and installed without setting the SUPERLUMT_WORKS option to TRUE.

Fixed the embedded coefficients for the ARKODE_TSITOURAS_7_4_5 Butcher table.

Fixed a bug in LSRKStep where an incorrect state vector could be passed to a user-supplied dominant eigenvalue function on the first step unless the output vector passed to ARKodeEvolve() contained the initial condition and when an eigenvalue estimate is requested on the first step in a subsequent call to ARKodeEvolve() unless the output vector passed contained the most recently returned solution.

Fixed a potential bug in LSRKStep’s ARKODE_LSRK_SSP_S_3 method, where a real number was used instead of an integer, potentially resulting in a rounding error.

Fixed a bug in MRIStep for estimating the first “slow” time step in an adaptive multirate calculation.

Fixed a bug in MRIStep when using a custom inner integrator that relies on the input state being the initial condition for the fast integration rather than retaining the result from the last inner integration or most recent reset call and the output vector passed to ARKodeEvolve() does not contain the initial condition on the first call or the last returned solution on subsequent calls.

Added a missing call to SUNNonlinSolSetup() in MRIStep when using an IMEX-MRI-SR method.

Fixed a bug in the ARKODE discrete adjoint checkpointing where an incorrect state would be stored on the first step if the output vector passed to ARKodeEvolve() did not contain the initial condition on the first call.

Removed extraneous copy of output vector when using ARKODE in ARK_ONE_STEP mode.

Removed an extraneous copy of the output vector in each step with SplittingStep.

Fixed a bug in logging output from ARKODE, where for some time stepping modules, the current “time” output in the logger was incorrect.

Fixed a bug where passing an empty string to SUNLogger_Set{Error,Warning,Info,Debug}Filename did not disable the corresponding logging stream Issue #844.

Deprecation Notices

The CVodeSetMonitorFn and CVodeSetMonitorFrequency functions have been deprecated and will be removed in the next major release.

Several CMake options have been deprecated in favor of namespaced versions prefixed with SUNDIALS_ to avoid naming collisions in applications that include SUNDIALS directly within their CMake builds. Additionally, a consistent naming convention (SUNDIALS_ENABLE) is now used for all boolean options. The table below lists the old CMake option names and the new replacements.

Old Option

New Option

BUILD_ARKODE

SUNDIALS_ENABLE_ARKODE

BUILD_CVODE

SUNDIALS_ENABLE_CVODE

BUILD_CVODES

SUNDIALS_ENABLE_CVODES

BUILD_IDA

SUNDIALS_ENABLE_IDA

BUILD_IDAS

SUNDIALS_ENABLE_IDAS

BUILD_KINSOL

SUNDIALS_ENABLE_KINSOL

ENABLE_MPI

SUNDIALS_ENABLE_MPI

ENABLE_OPENMP

SUNDIALS_ENABLE_OPENMP

ENABLE_OPENMP_DEVICE

SUNDIALS_ENABLE_OPENMP_DEVICE

OPENMP_DEVICE_WORKS

SUNDIALS_ENABLE_OPENMP_DEVICE_CHECKS

ENABLE_PTHREAD

SUNDIALS_ENABLE_PTHREAD

ENABLE_CUDA

SUNDIALS_ENABLE_CUDA

ENABLE_HIP

SUNDIALS_ENABLE_HIP

ENABLE_SYCL

SUNDIALS_ENABLE_SYCL

ENABLE_LAPACK

SUNDIALS_ENABLE_LAPACK

LAPACK_WORKS

SUNDIALS_ENABLE_LAPACK_CHECKS

ENABLE_GINKGO

SUNDIALS_ENABLE_GINKGO

GINKGO_WORKS

SUNDIALS_ENABLE_GINKGO_CHECKS

ENABLE_MAGMA

SUNDIALS_ENABLE_MAGMA

MAGMA_WORKS

SUNDIALS_ENABLE_MAGMA_CHECKS

ENABLE_SUPERLUDIST

SUNDIALS_ENABLE_SUPERLUDIST

SUPERLUDIST_WORKS

SUNDIALS_ENABLE_SUPERLUDIST_CHECKS

ENABLE_SUPERLUMT

SUNDIALS_ENABLE_SUPERLUMT

SUPERLUMT_WORKS

SUNDIALS_ENABLE_SUPERLUMT_CHECKS

ENABLE_KLU

SUNDIALS_ENABLE_KLU

KLU_WORKS

SUNDIALS_ENABLE_KLU_CHECKS

ENABLE_HYPRE

SUNDIALS_ENABLE_HYPRE

HYPRE_WORKS

SUNDIALS_ENABLE_HYPRE_CHECKS

ENABLE_PETSC

SUNDIALS_ENABLE_PETSC

PETSC_WORKS

SUNDIALS_ENABLE_PETSC_CHECKS

ENABLE_TRILINOS

SUNDIALS_ENABLE_TRILINOS

ENABLE_RAJA

SUNDIALS_ENABLE_RAJA

ENABLE_XBRAID

SUNDIALS_ENABLE_XBRAID

XBRAID_WORKS

SUNDIALS_ENABLE_XBRAID_CHECKS

ENABLE_ONEMKL

SUNDIALS_ENABLE_ONEMKL

ONEMKL_WORKS

SUNDIALS_ENABLE_ONEMKL_CHECKS

ENABLE_CALIPER

SUNDIALS_ENABLE_CALIPER

ENABLE_ADIAK

SUNDIALS_ENABLE_ADIAK

ENABLE_KOKKOS

SUNDIALS_ENABLE_KOKKOS

KOKKOS_WORKS

SUNDIALS_ENABLE_KOKKOS_CHECKS

ENABLE_KOKKOS_KERNELS

SUNDIALS_ENABLE_KOKKOS_KERNELS

KOKKOS_KERNELS_WORKS

SUNDIALS_ENABLE_KOKKOS_KERNELS_CHECKS

BUILD_FORTRAN_MODULE_INTERFACE

SUNDIALS_ENABLE_FORTRAN

SUNDIALS_BUILD_WITH_PROFILING

SUNDIALS_ENABLE_PROFILING

SUNDIALS_BUILD_WITH_MONITORING

SUNDIALS_ENABLE_MONITORING

SUNDIALS_BUILD_PACKAGE_FUSED_KERNELS

SUNDIALS_ENABLE_PACKAGE_FUSED_KERNELS

EXAMPLES_ENABLE_C

SUNDIALS_ENABLE_C_EXAMPLES

EXAMPLES_ENABLE_CXX

SUNDIALS_ENABLE_CXX_EXAMPLES

EXAMPLES_ENABLE_F2003

SUNDIALS_ENABLE_FORTRAN_EXAMPLES

EXAMPLES_ENABLE_CUDA

SUNDIALS_ENABLE_CUDA_EXAMPLES

EXAMPLES_INSTALL

SUNDIALS_ENABLE_EXAMPLES_INSTALL

EXAMPLES_INSTALL_PATH

SUNDIALS_EXAMPLES_INSTALL_PATH

BUILD_BENCHMARKS

SUNDIALS_ENABLE_BENCHMARKS

BENCHMARKS_INSTALL_PATH

SUNDIALS_BENCHMARKS_INSTALL_PATH

SUNDIALS_BENCHMARK_OUTPUT_DIR

SUNDIALS_BENCHMARKS_OUTPUT_DIR

SUNDIALS_BENCHMARK_CALIPER_OUTPUT_DIR

SUNDIALS_BENCHMARKS_CALIPER_OUTPUT_DIR

SUNDIALS_BENCHMARK_NUM_CPUS

SUNDIALS_BENCHMARKS_NUM_CPUS

SUNDIALS_BENCHMARK_NUM_GPUS

SUNDIALS_BENCHMARKS_NUM_GPUS

ENABLE_ALL_WARNINGS

SUNDIALS_ENABLE_ALL_WARNINGS

ENABLE_WARNINGS_AS_ERRORS

CMAKE_COMPILE_WARNING_AS_ERROR

ENABLE_ADDRESS_SANITIZER

SUNDIALS_ENABLE_ADDRESS_SANITIZER

ENABLE_MEMORY_SANITIZER

SUNDIALS_ENABLE_MEMORY_SANITIZER

ENABLE_LEAK_SANITIZER

SUNDIALS_ENABLE_LEAK_SANITIZER

Following the updated CMake options, the macros listed below have been deprecated and replaced with versions that align with the new CMake options.

Old Macro

New Macro

SUNDIALS_BUILD_WITH_PROFILING

SUNDIALS_ENABLE_PROFILING

SUNDIALS_BUILD_WITH_MONITORING

SUNDIALS_ENABLE_MONITORING

SUNDIALS_BUILD_PACKAGE_FUSED_KERNELS

SUNDIALS_ENABLE_PACKAGE_FUSED_KERNELS

For changes in prior versions of SUNDIALS see §18.

2.1.2. Reading this User Guide

This user guide is a combination of general usage instructions and specific example programs. We expect that some readers will want to concentrate on the general instructions, while others will refer mostly to the examples, and the organization is intended to accommodate both styles.

The structure of this document is as follows:

2.1.3. SUNDIALS License and Notices

All SUNDIALS packages are released open source, under the BSD 3-Clause license for more details see the LICENSE and NOTICE files provided with all SUNDIALS packages.