3.4.2.5. Multigrid Reduction in Time with XBraid

The prior sections discuss using ARKStep in a traditional sequential time integration setting i.e., the solution is advanced from one time to the next where all parallelism resides within the evaluation of a step e.g., the computation of the right-hand side, (non)linear solves, vector operations etc. For example, when discretizing a partial differential equation using a method of lines approach the spatially-discretized equations comprise a large set of ordinary differential equations that can be evolved with ARKStep. In this case the parallelization lies in decomposing the spatial domain unknowns across distributed computational nodes. Considering the strong scaling case at a given spatial resolution, as the problem is spread across greater numbers of computational nodes scalability in the spatial dimension is exhausted and sequential time integration becomes a bottleneck. This bottleneck is largely driven by the hardware shift from faster clock speeds to greater concurrency to achieve performance gains. In this case, at the spatial scaling limit and with stagnant clock speeds, more time steps will lead to an increased runtime.

An alternative approach to sequential time integration is to solve for all time values simultaneously. One such approach is multigrid reduction in time [50] (MGRIT) which uses a highly parallel iterative method to expose parallelism in the time domain in addition to the spatial parallelization. Starting with an initial temporal grid the multilevel algorithm constructs successively coarser time grids and uses each coarse grid solution to improve the solution at the next finer scale. In the two level case the MGRIT algorithm is as follows:

  1. Relax the solution on the fine grid (parallel-in-time)

  2. Restrict the solution to the fine grid (time re-discretization).

  3. Solve the residual equation on the coarse grid (serial-in-time).

  4. Correct the fine grid solution (parallel-in-time).

Applying this algorithm recursively for the solve step above leads to the multilevel algorithm.

The XBraid library [1] implements the MGRIT algorithm in a non-intrusive manner, enabling the reuse of existing software for sequential time integration. The following sections describe the ARKStep + XBraid interface and the steps necessary to modify an existing code that already uses ARKStep to also use XBraid.

3.4.2.5.1. SUNBraid Interface

Interfacing ARKStep with XBraid requires defining two data structures. The first is the XBraid application data structure that contains the data necessary for carrying out a time step and is passed to every interface function (much like the user data pointer in SUNDIALS packages). For this structure the SUNBraid interface defines the generic SUNBraidApp structure described below that serves as the basis for creating integrator-specific or user-defined interfaces to XBraid. The second structure holds the problem state data at a certain time value. This structure is defined by the SUNBraidVector structure and simply contains an N_Vector. In addition to the two data structures several functions defined by the XBraid API are required. These functions include vector operations (e.g., computing vector sums or norms) as well as functions to initialize the problem state, access the current solution, and take a time step.

The ARKBraid interface, built on the SUNBraidApp and SUNBraidVector structures, provides all the functionaly needed combine ARKStep and XBraid for parallel-in-time integration. As such, only a minimal number of changes are necessary to update an exsting code that uses ARKStep to also use XBraid.

3.4.2.5.1.1. SUNBraidApp

As mentioned above the SUNBraid interface defines the SUNBraidApp structure to hold the data necessary to compute a time step. This structure, like other SUNDIALS generic objects, is defined as a structure consisting of an implementation specific content field and an operations structure comprised of a set of function pointers for implmentation-defined operations on the object. Specifically the SUNBraidApp type is defined as

/* Define XBraid App structure */
struct _braid_App_struct
{
  void        *content;
  SUNBraidOps ops;
};

/* Pointer to the interface object (same as braid_App) */
typedef struct _braid_App_struct *SUNBraidApp;

Here, the SUNBraidOps structure is defined as

/* Structure containing function pointers to operations */
struct _SUNBraidOps
{
  int (*getvectmpl)(braid_App app, N_Vector *tmpl);
};

/* Pointer to operations structure */
typedef struct _SUNBraidOps *SUNBraidOps;

The generic SUNBraidApp defines and implements the generic operations acting on a SUNBraidApp obejct. These generic functions are nothing but wrappers to access the specific implementation through the object’s operations structure. To illustrate this point we show below the implementation of the SUNBraidApp_GetVecTmpl() function:

/* Get a template vector from the integrator */
int SUNBraidApp_GetVecTmpl(braid_App app, N_Vector *y)
{
  if (app->ops->getvectmpl == NULL) return SUNBRAID_OPNULL;
  return app->ops->getvectmpl(app, y);
}

The SUNBraidApp operations are define below in §3.4.2.5.1.2.

3.4.2.5.1.2. SUNBraidOps

In this section we define the SUNBraidApp operations and, for each operation, we give the function signature, a description of the expected behavior, and an example usage of the function.

int SUNBraidApp_GetVecTmpl(braid_App app, N_Vector *y)

This function returns a vector to use as a template for creating new vectors with N_VClone().

Arguments:
  • app – input, a SUNBraidApp instance (XBraid app structure).

  • y – output, the template vector.

Return value:

If this function is not implemented by the SUNBraidApp implementation (i.e., the function pointer is NULL) then this function will return SUNBRAID_OPNULL. Otherwise the return value depends on the particular SUNBraidApp implementation. Users are encouraged to utilize the return codes defined in sundials/sundials_xbraid.h and listed in Table 3.2.

Usage:

/* Get template vector */
flag = SUNBraidApp_GetVecTmpl(app, y_ptr);
if (flag != SUNBRAID_SUCCESS) return flag;

3.4.2.5.1.3. SUNBraidApp Utility Functions

In addition to the generic SUNBraidApp operations the following utility functions are provided to assist in creating and destroying a SUNBraidApp instance.

int SUNBraidApp_NewEmpty(braid_App *app)

This function creates a new SUNBraidApp instance with the content and operations initialized to NULL.

Arguments:
  • app – output, an empty SUNBraidApp instance (XBraid app structure).

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ALLOCFAIL if a memory allocation failed.

Usage:

/* Create empty XBraid interface object */
flag = SUNBraidApp_NewEmpty(app_ptr);
if (flag != SUNBRAID_SUCCESS) return flag;
int SUNBraidApp_FreeEmpty(braid_App *app)

This function destroys an empty SUNBraidApp instance.

Arguments:
  • app – input, an empty SUNBraidApp instance (XBraid app structure).

Return value:
  • SUNBRAID_SUCCESS if successful.

Usage:

/* Free empty XBraid interface object */
flag = SUNBraidApp_FreeEmpty(app_ptr);

Warning

This function does not free the SUNBraidApp object’s content structure. An implementation should free its content before calling SUNBraidApp_FreeEmpty() to deallocate the base SUNBraidApp structure.

3.4.2.5.1.4. SUNBraidVector

As mentioned above the SUNBraid interface defines the SUNBraidVector structure to store a snapshot of solution data at a single point in time and this structure simply contains an N_Vector. Specifically, the structure is defined as follows:

struct _braid_Vector_struct
{
  N_Vector y;
};

/* Poiner to vector wrapper (same as braid_Vector) */
typedef struct _braid_Vector_struct *SUNBraidVector;

To assist in creating creating and destroying this structure the following utility functions are provided.

int SUNBraidVector_New(N_Vector y, SUNBraidVector *u)

This function creates a new SUNBraidVector wrapping the N_Vector y.

Arguments:
  • y – input, the N_Vector to wrap.

  • u – output, the SUNBraidVector wrapping y.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if y is NULL.

  • SUNBRAID_ALLOCFAIL if a memory allocation fails.

Usage:

/* Create new vector wrapper */
flag = SUNBraidVector_New(y, u_ptr);
if (flag != SUNBRAID_SUCCESS) return flag;

Warning

The SUNBraidVector takes ownership of the wrapped N_Vector and as such the wrapped N_Vector is destroyed when the SUNBraidVector is freed with SUNBraidVector_Free().

int SUNBraidVector_GetNVector(SUNBraidVector u, N_Vector *y)

This function retrieves the wrapped N_Vector from the SUNBraidVector.

Arguments:
  • u – input, the SUNBraidVector wrapping y.

  • y – output, the wrapped N_Vector.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if u is NULL.

  • SUNBRAID_MEMFAIL if y is NULL.

Usage:

/* Create new vector wrapper */
flag = SUNBraidVector_GetNVector(u, y_ptr);
if (flag != SUNBRAID_SUCCESS) return flag;

Finally, the SUNBraid interface defines the following vector operations acting on SUNBraidVectors, that consist of thin wrappers to compatible SUNDIALS N_Vector operations.

int SUNBraidVector_Clone(braid_App app, braid_Vector u, braid_Vector *v_ptr)

This function creates a clone of the input SUNBraidVector and copies the values of the input vector u into the output vector v_ptr using N_VClone() and N_VScale().

Arguments:
  • app – input, a SUNBraidApp instance (XBraid app structure).

  • u – input, the SUNBraidVector to clone.

  • v_ptr – output, the new SUNBraidVector.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if u is NULL.

  • SUNBRAID_MEMFAIL if the N_Vector y wrapped by u is NULL.

  • SUNBRAID_ALLOCFAIL if a memory allocation fails.

int SUNBraidVector_Free(braid_App app, braid_Vector u)

This function destroys the SUNBraidVector and the wrapped N_Vector using N_VDestroy().

Arguments:
  • app – input, a SUNBraidApp instance (XBraid app structure).

  • u – input, the SUNBraidVector to destroy.

Return value:
  • SUNBRAID_SUCCESS if successful.

int SUNBraidVector_Sum(braid_App app, braid_Real alpha, braid_Vector x, braid_Real beta, braid_Vector y)

This function computes the vector sum \(\alpha x + \beta y \rightarrow y\) using N_VLinearSum().

Arguments:
  • app – input, a SUNBraidApp instance (XBraid app structure).

  • alpha – input, the constant \(\alpha\).

  • x – input, the vector \(x\).

  • beta – input, the constant \(\beta\).

  • y – input/output, the vector \(y\).

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if x or y is NULL.

  • SUNBRAID_MEMFAIL if either of the wrapped N_Vectors are NULL.

int SUNBraidVector_SpatialNorm(braid_App app, braid_Vector u, braid_Real *norm_ptr)

This function computes the 2-norm of the vector u using N_VDotProd().

Arguments:
  • app – input, a SUNBraidApp instance (XBraid app structure).

  • u – input, the vector u.

  • norm_ptr – output, the L2 norm of u.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if u is NULL.

  • SUNBRAID_MEMFAIL if the wrapped N_Vector is NULL.

int SUNBraidVector_BufSize(braid_App app, braid_Int *size_ptr, braid_BufferStatus bstatus)

This function returns the buffer size for messages to exchange vector data using SUNBraidApp_GetVecTmpl() and N_VBufSize().

Arguments:
  • app – input, a SUNBraidApp instance (XBraid app structure).

  • size_ptr – output, the buffer size.

  • bstatus – input, a status object to query for information on the message type.

Return value:
int SUNBraidVector_BufPack(braid_App app, braid_Vector u, void *buffer, braid_BufferStatus bstatus)

This function packs the message buffer for exchanging vector data using N_VBufPack().

Arguments:
  • app – input, a SUNBraidApp instance (XBraid app structure).

  • u – input, the vector to pack into the exchange buffer.

  • buffer – output, the packed exchange buffer to pack.

  • bstatus – input, a status object to query for information on the message type.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if u is NULL.

  • An error flag from N_VBufPack().

int SUNBraidVector_BufUnpack(braid_App app, void *buffer, braid_Vector *u_ptr, braid_BufferStatus bstatus)

This function unpacks the message buffer and creates a new N_Vector and SUNBraidVector with the buffer data using N_VBufUnpack(), SUNBraidApp_GetVecTmpl(), and N_VClone().

Arguments:
  • app – input, a SUNBraidApp instance (XBraid app structure).

  • buffer – input, the exchange buffer to unpack.

  • u_ptr – output, a new SUNBraidVector containing the buffer data.

  • bstatus – input, a status object to query for information on the message type.

Return value:

3.4.2.5.1.5. SUNBraid Return Codes

The SUNBraid interface return values are given in Table 3.2.

Table 3.2 SUNBraid Return Codes

Return value name

Value

Meaning

SUNBRAID_SUCCESS

\(0\)

The call/operation was successful.

SUNBRAID_ALLOCFAIL

\(-1\)

A memory allocation failed.

SUNBRAID_MEMFAIL

\(-2\)

A memory access fail.

SUNBRAID_OPNULL

\(-3\)

The SUNBraid operation is NULL.

SUNBRAID_ILLINPUT

\(-4\)

An invalid input was provided.

SUNBRAID_BRAIDFAIL

\(-5\)

An XBraid function failed.

SUNBRAID_SUNFAIL

\(-6\)

A SUNDIALS function failed.

3.4.2.5.2. ARKBraid Interface

This section describes the ARKBraid implementation of a SUNBraidApp for using the ARKStep integration module with XBraid. The following section §3.4.2.5.2.1 describes routines for creating, initializing, and destroying the ARKStep + XBraid interface, routines for setting optional inputs, and routines for retrieving data from an ARKBraid instance. As noted above, interfacing with XBraid requires providing functions to initialize the problem state, access the current solution, and take a time step. The default ARKBraid functions for each of these actions are defined in §3.4.2.5.2.4 and may be overridden by user-defined if desired. A skeleton of the user’s main or calling program for using the ARKBraid interface is given in §3.4.2.5.3. Finally, for advanced users that wish to create their own SUNBraidApp implementation using ARKStep, §3.4.2.5.4 describes some helpful functions available to the user.

3.4.2.5.2.1. ARKBraid Initialization and Deallocation Functions

This section describes the functions that are called by the user to create, initialize, and destroy an ARKBraid instance. Each user-callable function returns SUNBRAID_SUCCESS (i.e., 0) on a successful call and a negative value if an error occurred. The possible return codes are given in Table 3.2.

int ARKBraid_Create(void *arkode_mem, braid_App *app)

This function creates a SUNBraidApp object, sets the content pointer to the private ARKBraid interface structure, and attaches the necessary SUNBraidOps implementations.

Arguments:
  • arkode_mem – input, a pointer to an ARKStep memory structure.

  • app – output, an ARKBraid instance (XBraid app structure).

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT arkode_mem is NULL.

  • SUNBRAID_ALLOCFAIL if a memory allocation failed.

Warning

The ARKBraid interface is ARKStep-specific. Although one could eventually construct an XBraid interface to either ERKStep or MRIStep, those are not supported by this implementation.

int ARKBraid_BraidInit(MPI_Comm comm_w, MPI_Comm comm_t, sunrealtype tstart, sunrealtype tstop, sunindextype ntime, braid_App app, braid_Core *core)

This function wraps the XBraid braid_Init() function to create the XBraid core memory structure and initializes XBraid with the ARKBraid and SUNBraidVector interface functions.

Arguments:
  • comm_w – input, the global MPI communicator for space and time.

  • comm_t – input, the MPI communicator for the time dimension.

  • tstart – input, the initial time value.

  • tstop – input, the final time value.

  • ntime – input, the initial number of grid points in time.

  • app – input, an ARKBraid instance.

  • core – output, the XBraid core memory structure.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if either MPI communicator is MPI_COMM_NULL, if ntime < 2, or if app or its content is NULL.

  • SUNBRAID_BRAIDFAIL if the braid_Init() call fails. The XBraid return value can be retrieved with ARKBraid_GetLastBraidFlag().

Note

If desired, the default functions for vector initialization, accessing the solution, taking a time step, and computing the spatial norm should be overridden before calling this function. See §3.4.2.5.2.2 for more details.

Warning

The user is responsible for deallocating the XBraid core memory structure with the XBraid function braid_Destroy().

int ARKBraid_Free(braid_App *app)

This function deallocates an ARKBraid instance.

Arguments:
  • app – input, a pointer to an ARKBraid instance.

Return value:
  • SUNBRAID_SUCCESS if successful.

3.4.2.5.2.2. ARKBraid Set Functions

This section describes the functions that are called by the user to set optional inputs to control the behavior of an ARKBraid instance or to provide alternative XBraid interface functions. Each user-callable function returns SUNBRAID_SUCCESS (i.e., 0) on a successful call and a negative value if an error occurred. The possible return codes are given in Table 3.2.

int ARKBraid_SetStepFn(braid_App app, braid_PtFcnStep step)

This function sets the step function provided to XBraid (default ARKBraid_Step()).

Arguments:
  • app – input, an ARKBraid instance.

  • step – input, an XBraid step function. If step is NULL, the default function will be used.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if app is NULL.

  • SUNBRAID_MEMFAIL if the app content is NULL.

Note

This function must be called prior to ARKBraid_BraidInit().

int ARKBraid_SetInitFn(braid_App app, braid_PtFcnInit init)

This function sets the vector initialization function provided to XBraid (default ARKBraid_Init()).

Arguments:
  • app – input, an ARKBraid instance.

  • init – input, an XBraid vector initialization function. If init is NULL, the default function will be used.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if app is NULL.

  • SUNBRAID_MEMFAIL if the app content is NULL.

Note

This function must be called prior to ARKBraid_BraidInit().

int ARKBraid_SetSpatialNormFn(braid_App app, braid_PtFcnSpatialNorm snorm)

This function sets the spatial norm function provided to XBraid (default SUNBraid_SpatialNorm()).

Arguments:
  • app – input, an ARKBraid instance.

  • snorm – input, an XBraid spatial norm function. If snorm is NULL, the default function will be used.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if app is NULL.

  • SUNBRAID_MEMFAIL if the app content is NULL.

Note

This function must be called prior to ARKBraid_BraidInit().

int ARKBraid_SetAccessFn(braid_App app, braid_PtFcnAccess access)

This function sets the user access function provided to XBraid (default ARKBraid_Access()).

Arguments:
  • app – input, an ARKBraid instance.

  • init – input, an XBraid user access function. If access is NULL, the default function will be used.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if app is NULL.

  • SUNBRAID_MEMFAIL if the app content is NULL.

Note

This function must be called prior to ARKBraid_BraidInit().

3.4.2.5.2.3. ARKBraid Get Functions

This section describes the functions that are called by the user to retrieve data from an ARKBraid instance. Each user-callable function returns SUNBRAID_SUCCESS (i.e., 0) on a successful call and a negative value if an error occurred. The possible return codes are given in Table 3.2.

int ARKBraid_GetVecTmpl(braid_App app, N_Vector *tmpl)

This function returns a vector from the ARKStep memory to use as a template for creating new vectors with N_VClone() i.e., this is the ARKBraid implementation of SUNBraidVector_GetVecTmpl().

Arguments:
  • app – input, an ARKBraid instance.

  • tmpl – output, a template vector.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if app is NULL.

  • SUNBRAID_MEMFAIL if the app content or ARKStep memory is NULL.

int ARKBraid_GetARKStepMem(braid_App app, void **arkode_mem)

This function returns the ARKStep memory structure pointer attached with ARKBraid_Create().

Arguments:
  • app – input, an ARKBraid instance.

  • arkode_mem – output, a pointer to the ARKStep memory structure.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if app is NULL.

  • SUNBRAID_MEMFAIL if the app content or ARKStep memory is NULL.

int ARKBraid_GetUserData(braid_App app, void **user_data)

This function returns the user data pointer attached with ARKStepSetUserData().

Arguments:
  • app – input, an ARKBraid instance.

  • user_data – output, a pointer to the user data structure.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if app is NULL.

  • SUNBRAID_MEMFAIL if the app content or ARKStep memory is NULL.

int ARKBraid_GetLastBraidFlag(braid_App app, int *last_flag)

This function returns the return value from the most recent XBraid function call.

Arguments:
  • app – input, an ARKBraid instance.

  • last_flag – output, the XBraid return value.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if app is NULL.

  • SUNBRAID_MEMFAIL if the app content is NULL.

int ARKBraid_GetLastARKStepFlag(braid_App app, int *last_flag)

This function returns the return value from the most recent ARKStep function call.

Arguments:
  • app – input, an ARKBraid instance.

  • last_flag – output, the ARKStep return value.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if app is NULL.

  • SUNBRAID_MEMFAIL if the app content is NULL.

int ARKBraid_GetSolution(braid_App app, sunrealtype *tout, N_Vector yout)

This function returns final time and state stored with the default access function ARKBraid_Access().

Arguments:
  • app – input, an ARKBraid instance.

  • last_flag – output, the ARKStep return value.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if app is NULL.

  • SUNBRAID_MEMFAIL if the app content or the stored vector is NULL.

Warning

If providing a non-default access function the final time and state are not stored within the ARKBraid structure and this function will return an error. In this case the user should allocate space to store any desired output within the user data pointer attached to ARKStep with ARKStepSetUserData(). This user data pointer can be retrieved from the ARKBraid structure with ARKBraid_GetUserData().

3.4.2.5.2.4. ARKBraid Interface Functions

This section describes the default XBraid interface functions provided by ARKBraid and called by XBraid to perform certain actions. Any or all of these functions may be overridden by supplying a user-defined function through the set functions defined in §3.4.2.5.2.2. Each default interface function returns SUNBRAID_SUCCESS (i.e., 0) on a successful call and a negative value if an error occurred. The possible return codes are given in Table 3.2.

int ARKBraid_Step(braid_App app, braid_Vector ustop, braid_Vector fstop, braid_Vector u, braid_StepStatus status)

This is the default step function provided to XBraid. The step function is called by XBraid to advance the vector u from one time to the next using the ARStep memory structure provided to ARKBraid_Create(). A user-defined step function may be set with ARKBraid_SetStepFn().

Arguments:
  • app – input, an ARKBraid instance.

  • ustop – input, u vector at the new time tstop.

  • fstop – input, the right-hand side vector at the new time tstop.

  • u - input/output, on input the vector at the start time and on return the vector at the new time.

  • status – input, a status object to query for information about u and to steer XBraid e.g., for temporal refinement.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if app is NULL.

  • SUNBRAID_MEMFAIL if the app content or ARKStep memory is NULL.

  • SUNBRAID_BRAIDFAIL if an XBraid function fails. The return value can be retrieved with ARKBraid_GetLastBraidFlag().

  • SUNBRAID_SUNFAIL if a SUNDIALS function fails. The return value can be retrieved with ARKBraid_GetLastARKStepFlag().

Note

If providing a non-default implementation of the step function the utility function ARKBraid_TakeStep() should be used to advance the input vector u to the new time.

int ARKBraid_Init(braid_App app, sunrealtype t, braid_Vector *u_ptr)

This is the default vector initialization function provided to XBraid. The initialization function is called by XBraid to create a new vector and set the initial guess for the solution at time \(t\). When using this default function the initial guess at all time values is the initial condition provided to ARKStepCreate(). A user-defined init function may be set with ARKBraid_SetInitFn().

Arguments:
  • app – input, an ARKBraid instance.

  • t – input, the initialization time for the output vector.

  • u_ptr – output, the new and initialized SUNBraidVector.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if app is NULL.

  • SUNBRAID_MEMFAIL if the app content or ARKStep memory is NULL.

  • SUNBRAID_ALLOCFAIL if a memory allocation failed.

Note

If providing a non-default implementation of the vector initialization function the utility functions SUNBraidApp_GetVecTmpl() and SUNBraidVector_New() can be helpful when creating the new vector returned by this function.

int ARKBraid_Access(braid_App app, braid_Vector u, braid_AccessStatus astatus)

This is the default access function provided to XBraid. The access function is called by XBraid to retrieve the current solution. When using this default function the final solution time and state are stored within the ARKBraid structure. This information can be retrieved with ARKBraid_GetSolution(). A user-defined access function may be set with ARKBraid_SetAccessFn().

Arguments:
  • app – input, an ARKBraid instance.

  • u – input, the vector to be accessed.

  • status – input, a status object to query for information about u.

Return value:
  • SUNBRAID_SUCCESS if successful.

  • SUNBRAID_ILLINPUT if any of the inputs are NULL.

  • SUNBRAID_MEMFAIL if the app content, the wrapped N_Vector, or the ARKStep memory is NULL.

  • SUNBRAID_ALLOCFAIL if allocating storage for the final solution fails.

  • SUNBRAID_BRAIDFAIL if an XBraid function fails. The return value can be retrieved with ARKBraid_GetLastBraidFlag().

3.4.2.5.3. A skeleton of the user’s main program with XBraid

In addition to the header files required for the integration of the ODE problem (see the section §3.4.1), to use the ARKBraid interace, the user’s program must include the header file arkode/arkode_xbraid.h which declares the needed function prototypes.

The following is a skeleton of the user’s main program (or calling program) for the integration of an ODE IVP using ARKStep with XBraid for parallel-in-time integration. Most steps are unchanged from the skeleton program presented in §3.4.2.1. New or updated steps are bold.

  1. Initialize MPI

    If parallelizing in space and time split the global communicator into communicators for space and time with braid_SplitCommworld().

  2. Set problem dimensions

  3. Set vector of initial values

  4. Create ARKStep object

  5. Specify integration tolerances

  6. Create matrix object

  7. Create linear solver object

  8. Set linear solver optional inputs

  9. Attach linear solver module

  10. Create nonlinear solver object

  11. Attach nonlinear solver module

  12. Set nonlinear solver optional inputs

  13. Set optional inputs

  14. Create ARKBraid interface

    Call the constructor ARKBraid_Create() to create the XBraid app structure.

  15. Set optional ARKBraid inputs

    See §3.4.2.5.2.2 for ARKBraid inputs.

  16. Initialize the ARKBraid interface

    Call the initialization function ARKBraid_BraidInit() to create the XBraid core memory structure and attach the ARKBraid interface app and functions.

  17. Set optional XBraid inputs

    See the XBraid documentation for available XBraid options.

  18. Evolve the problem

    Call braid_Drive() to evolve the problem with MGRIT.

  19. Get optional outputs

    See §3.4.2.5.2.3 for ARKBraid outputs.

  20. Deallocate memory for solution vector

  21. Free solver memory

  22. Free linear solver memory

  23. Free ARKBraid and XBraid memory

    Call ARKBraid_Free() and braid_Destroy to deallocate the ARKBraid interface and and XBraid core memory structures, respectively.

  24. Finalize MPI

3.4.2.5.4. Advanced ARKBraid Utility Functions

This section describes utility functions utilized in the ARKStep + XBraid interfacing. These functions are used internally by the above ARKBraid interface functions but are exposed to the user to assist in advanced usage of ARKODE and XBraid that requries defining a custom SUNBraidApp implementation.

int ARKBraid_TakeStep(void *arkode_mem, sunrealtype tstart, sunrealtype tstop, N_Vector y, int *ark_flag)

This function advances the vector y from tstart to tstop using a single ARKStep time step with step size h = tstop - start.

Arguments:
  • arkode_mem – input, the ARKStep memory structure pointer.

  • tstart – input, the step start time.

  • tstop – input, the step stop time.

  • y – input/output, on input the solution a tstop and on return, the solution at time tstop if the step was successful (ark_flag \(\geq 0\)) or the solution at time tstart if the step failed (ark_flag < 0).

  • ark_flag – output, the step status flag. If ark_flag is:

    \(= 0\) then the step succeeded and, if applicable, met the requested temporal accuracy.

    \(> 0\) then the step succeeded but failed to meet the requested temporal accuracy.

    \(< 0\) then the step failed e.g., a solver failure occurred.

Return value:

If all ARKStep function calls are successful the return value is ARK_SUCCESS, otherwise the return value is the error flag returned from the function that failed.