Reference

class leap.MethodBuilder[source]

An abstract base class for method implementations that generate code for dagrt.

generate(*solver_hooks) DAGCode[source]

Generate a method description.

Parameters:

solver_hooks – A list of callbacks that generate expressions for calling user-supplied implicit solvers

implicit_expression(expression_tag=None)[source]

Return a template that expressions in class:AssignImplicit instances will follow.

Parameters:

expression_tag – A name for the expression, if multiple expressions are present in the generated code.

Returns:

A tuple consisting of pymbolic expressions and the names of the free variables in the expressions.

Runge-Kutta Methods

leap.rk.ORDER_TO_RK_METHOD_BUILDER

A dictionary mapping desired order of accuracy to a corresponding RK method builder.

class leap.rk.ForwardEulerMethodBuilder(component_id, state_filter_name=None, rhs_func_name=None)[source]
__init__(component_id, state_filter_name=None, rhs_func_name=None)[source]
generate()[source]
Returns:

dagrt.language.DAGCode

class leap.rk.BackwardEulerMethodBuilder(component_id, state_filter_name=None, rhs_func_name=None)[source]
__init__(component_id, state_filter_name=None, rhs_func_name=None)[source]
generate()[source]
Returns:

dagrt.language.DAGCode

class leap.rk.MidpointMethodBuilder(component_id, state_filter_name=None, rhs_func_name=None)[source]
__init__(component_id, state_filter_name=None, rhs_func_name=None)[source]
generate()[source]
Returns:

dagrt.language.DAGCode

class leap.rk.HeunsMethodBuilder(component_id, state_filter_name=None, rhs_func_name=None)[source]
__init__(component_id, state_filter_name=None, rhs_func_name=None)[source]
generate()[source]
Returns:

dagrt.language.DAGCode

class leap.rk.RK3MethodBuilder(component_id, state_filter_name=None, rhs_func_name=None)[source]

Source: J. C. Butcher, Numerical MethodBuilders for Ordinary Differential Equations, 2nd ed., pages 94 - 99

__init__(component_id, state_filter_name=None, rhs_func_name=None)[source]
generate()[source]
Returns:

dagrt.language.DAGCode

class leap.rk.RK4MethodBuilder(component_id, state_filter_name=None, rhs_func_name=None)[source]
__init__(component_id, state_filter_name=None, rhs_func_name=None)[source]
generate()[source]
Returns:

dagrt.language.DAGCode

class leap.rk.RK5MethodBuilder(component_id, state_filter_name=None, rhs_func_name=None)[source]

Source: J. C. Butcher, Numerical MethodBuilders for Ordinary Differential Equations, 2nd ed., pages 94 - 99

Low-Storage Methods

class leap.rk.LSRK4MethodBuilder(component_id, state_filter_name=None, rhs_func_name=None)[source]

A low storage fourth-order Runge-Kutta method

See JSH, TW: Nodal Discontinuous Galerkin MethodBuilders p.64 or Carpenter, M.H., and Kennedy, C.A., Fourth-order-2N-storage Runge-Kutta schemes, NASA Langley Tech Report TM 109112, 1994

__init__(component_id, state_filter_name=None, rhs_func_name=None)[source]
Parameters:

component_id – an identifier to be used for the single state component supported.

generate()[source]
Returns:

dagrt.language.DAGCode

Adaptive/Embedded Methods

class leap.rk.ODE23MethodBuilder(component_id, use_high_order=True, state_filter_name=None, atol=0, rtol=0, max_dt_growth=None, min_dt_shrinkage=None)[source]

Bogacki-Shampine second/third-order Runge-Kutta.

(same as Matlab’s ode23)

Bogacki, Przemyslaw; Shampine, Lawrence F. (1989), “A 3(2) pair of Runge-Kutta formulas”, Applied Mathematics Letters 2 (4): 321-325, https://doi.org/10.1016/0893-9659(89)90079-7

__init__(component_id, use_high_order=True, state_filter_name=None, atol=0, rtol=0, max_dt_growth=None, min_dt_shrinkage=None)[source]
generate()[source]
Returns:

dagrt.language.DAGCode

class leap.rk.ODE45MethodBuilder(component_id, use_high_order=True, state_filter_name=None, atol=0, rtol=0, max_dt_growth=None, min_dt_shrinkage=None)[source]

Dormand-Prince fourth/fifth-order Runge-Kutta.

(same as Matlab’s ode45)

Dormand, J. R.; Prince, P. J. (1980), “A family of embedded Runge-Kutta formulae”, Journal of Computational and Applied Mathematics 6 (1): 19-26, http://dx.doi.org/10.1016/0771-050X(80)90013-3.

__init__(component_id, use_high_order=True, state_filter_name=None, atol=0, rtol=0, max_dt_growth=None, min_dt_shrinkage=None)[source]
generate()[source]
Returns:

dagrt.language.DAGCode

Strong Stability Preserving (SSP) Methods

class leap.rk.SSPRKMethodBuilder(component_id, state_filter_name=None, rhs_func_name=None)[source]

Explicit Strong Stability Preserving (SSP) Runge-Kutta Methods.

The methods are given in the now-standard Shu-Osher form

\[\begin{split}\begin{aligned} y^{(i)} =\,\, & \sum_{k = 0}^{i - 1} \alpha_{ik} y^{(i)} + \Delta t \beta_{ik} f(y^{(i)}), \\ y^{n + 1} =\,\, & y^{(n)}, \end{aligned}\end{split}\]

for \(i \in \{1, \dots, n\}\) and \(y^{(0)} = y^n\). For reference, see [gst-2001].

[gst-2001] (1,2,3)

S. Gottlieb, C.-W. Shu, E. Tadmor, Strong Stability Preserving High-Order Time Discretization Methods, SIAM, Vol. 43, pp. 89-112, 2001. https://doi.org/10.1137/S003614450036757X

class leap.rk.SSPRK22MethodBuilder(component_id, state_filter_name=None, rhs_func_name=None)[source]

Second-order SSP Runge-Kutta method from [gst-2001] Proposition 4.1.

generate()[source]
Returns:

a DAGCode that can be used to generate code for this method.

class leap.rk.SSPRK33MethodBuilder(component_id, state_filter_name=None, rhs_func_name=None)[source]

Third-order SSP Runge-Kutta method from [gst-2001] Proposition 4.1.

generate()[source]
Returns:

a DAGCode that can be used to generate code for this method.

Computing a stability function

leap.rk.stability_function(rk_a, rk_b)[source]

Given the matrix A and the ‘output coefficients’ b from a Butcher tableau, return the stability function of the method as a sympy.core.expr.Expr.

IMEX Method Builders

class leap.rk.imex.KennedyCarpenterIMEXRungeKuttaMethodBuilderBase(component_id, use_high_order=True, state_filter_name=None, use_explicit=True, use_implicit=True, atol=0, rtol=0, max_dt_growth=None, min_dt_shrinkage=None, implicit_rhs_name=None, explicit_rhs_name=None)[source]

Christopher A. Kennedy, Mark H. Carpenter. Additive Runge-Kutta schemes for convection-diffusion-reaction equations. Applied Numerical Mathematics Volume 44, Issues 1-2, January 2003, Pages 139-181 http://dx.doi.org/10.1016/S0168-9274(02)00138-1

User-supplied context:

<state> + component_id: The value that is integrated <func>rhs_expl_ + component_id: The explicit right hand side function <func>rhs_impl_ + component_id: The implicit right hand side function

__init__(component_id, use_high_order=True, state_filter_name=None, use_explicit=True, use_implicit=True, atol=0, rtol=0, max_dt_growth=None, min_dt_shrinkage=None, implicit_rhs_name=None, explicit_rhs_name=None)[source]
generate()[source]
Returns:

dagrt.language.DAGCode

Multi-Step Methods

class leap.multistep.AdamsIntegrationFunctionFamily[source]

An abstract interface for function families used for Adams-type time integration.

__len__()[source]
evaluate(func_idx, x)[source]
antiderivative(func_idx, x)[source]
class leap.multistep.AdamsMonomialIntegrationFunctionFamily(order)[source]

Implements AdamsMonomialIntegrationFunctionFamily.

class leap.multistep.AdamsBashforthMethodBuilder(component_id, function_family=None, state_filter_name=None, hist_length=None, static_dt=False, order=None)[source]
User-supplied context:

<state> + component_id: The value that is integrated <func> + component_id: The right hand side

__init__(component_id, function_family=None, state_filter_name=None, hist_length=None, static_dt=False, order=None)[source]
Parameters:
  • function_family – Accepts an instance of AdamsIntegrationFunctionFamily or an integer, in which case the classical monomial function family with the order given by the integer is used.

  • static_dt – If True, changing the timestep during time integration is not allowed.

generate()[source]
Returns:

dagrt.language.DAGCode

Multi-Rate Multi-Step Methods

class leap.multistep.multirate.rhs_policy[source]
late
early
early_and_late
class leap.multistep.multirate.MultiRateHistory(interval: int, func_name: str, arguments: Tuple[str, ...], order: int | None = None, rhs_policy: int = 0, invalidate_computed_state: bool = False, hist_length: int | None = None)[source]
arguments: Tuple[str, ...]

A tuple of component names which are passed to func_name.

func_name: str

The name of the function to call on history updates.

hist_length: int | None = None

The length of the history. If the length is greater than order, we use a least-squares solve rather than a linear sole to obtain the Adams coefficients for this history.

interval: int

An integer indicating the interval (relative to the smallest available time step) at which the history is to be update, where an update will typically involve a call to func_name.

invalidate_computed_state: bool = False

A flag that denotes whether evaluating the right-hand side func_name should force a recomputation of any state that depended upon now-superseded state.

order: int | None = None

The approximation order of the Adams method to be used for this history. If None, the default method is assumed.

rhs_policy: int = 0

One of the constants in rhs_policy.

class leap.multistep.multirate.MultiRateMultiStepMethodBuilder(default_order, system_description, state_filter_names=None, component_arg_names=None, static_dt=False, hist_consistency_threshold=None, early_hist_consistency_threshold=None)[source]

Simultaneously timesteps multiple parts of an ODE system, each with adjustable orders, rates, and dependencies.

Considerably generalizes [GearWells].

[GearWells]

C.W. Gear and D.R. Wells, “Multirate linear multistep methods,” BIT Numerical Mathematics, vol. 24, Dec. 1984,pg. 484-502.

__init__(default_order, system_description, state_filter_names=None, component_arg_names=None, static_dt=False, hist_consistency_threshold=None, early_hist_consistency_threshold=None)[source]
Parameters:
  • default_order – The order to be used for right-hand sides where no differing order is specified.

  • system_description

    A tuple of the form:

    (
        ("dt", "fast",
        "=", MultiRateHistory(1, "<func>f1", ("fast", "slow", "dep")),
        MultiRateHistory(3, "<func>f2", ("slow", "dep")),
        ),
        ("dt", "slow",
        "=", MultiRateHistory(3, "<func>f3", ("fast", "slow", "dep")),
        ),
        ("dep",
        "=", MultiRateHistory(3, "<func>f4", ("slow")),
        ),
    )
    

    i.e. the outermost tuple represents a list of state components. These can either be ODEs (in which case the first string in the tuple is "dt", followed by the component name, or computed/dependent state, in which case the first string in the tuple is just the name of the computed piece of state.

    The right-hand-side of each tuple (after the intervening "=" string) consists of one or more instances of MultiRateHistory describing the rate at which evaluations of the given functions should be stored (along with other parameters that can be configured in MultiRateHistory). The right-hand sides are combined additively.

    Computed state components may not (directly) or indirectly depend on themselves. The result will be undefined.

  • state_filter_names – a dictionary mapping (non-ODE) component names from the system description to the names of “state filter” functions (as in functions that receive this particular component state and return a potentially modified version thereof).

  • component_arg_names – A tuple of names of the components to be used as keywords for passing arguments to the right hand sides in rhss.

  • static_dt – If True, changing the timestep in between steps during time integration is not allowed.

generate(explainer=None)[source]
Parameters:

explainer – a subclass of SchemeExplainerBase, possibly TextualSchemeExplainer, or None.

Returns:

dagrt.language.DAGCode

class leap.multistep.multirate.SchemeExplainerBase[source]
class leap.multistep.multirate.TextualSchemeExplainer[source]
__init__()[source]
__str__()[source]

Return str(self).

Analysis

class leap.step_matrix.StepMatrixFinder(code, function_map, variables=None, exclude_variables=None)[source]

Constructs a step matrix on-the-fly while interpreting code.

Assumes that all function evaluations occur as the root node of a separate assignment statement.

leap.step_matrix.fast_evaluator(matrix, sparse=False)[source]

Generate a function to evaluate a step matrix quickly. The input comes from StepMatrixFinder.

leap.stability.find_stability_region(code, parallel=None, n_angles=100, prec=0.01, origin=-0.3)[source]

Transformation

leap.transform.strang_splitting(dag1, dag2, stepping_phase)[source]

Given two time advancement routines (in dag1 and dag2), returns a single second-order accurate time advancement routine representing the sum of both of those advancements.

Parameters:
Returns:

a dagrt.language.DAGCode