Time Derivatives

Time derivative operators approximate \(D\phi/Dt\) or \(DF/Dt\) for transient solvers. All operators share a common interface: update_pre_solve(dt) before each timestep, bdf() for the BDF approximation in the weak form, and update_post_solve(dt) after the solve completes.

History is initialised automatically on the first solve call, and BDF order ramps from 1 up to the requested order over the first few timesteps.

Base Class

Symbolic

class underworld3.systems.ddt.Symbolic[source]

Bases: uw_object

Symbolic history manager for time derivative approximations.

Manages the update of a variable \(\psi\) across timesteps. The history operator stores \(\psi\) over several timesteps (given by order) so that it can compute backward differentiation (BDF) or Adams-Moulton expressions.

The history operator is defined as:

\[\begin{split}\psi_p^{t-n\Delta t} &\leftarrow \psi_p^{t-(n-1)\Delta t} \\ \psi_p^{t-(n-1)\Delta t} &\leftarrow \psi_p^{t-(n-2)\Delta t} \cdots \\ \psi_p^{t-\Delta t} &\leftarrow \psi_p^{t}\end{split}\]

This is a purely symbolic history manager that operates on sympy expressions without mesh or swarm storage. It is useful for tracking symbolic expressions through time-stepping algorithms.

Parameters:
  • psi_fn (sympy.Basic) – The sympy expression to track. Can be scalar or matrix form.

  • theta (float, optional) – Implicitness parameter for Adams-Moulton order 1 (default 0.5). Values: 0 = explicit, 1 = implicit, 0.5 = Crank-Nicolson.

  • varsymbol (str, optional) – LaTeX symbol for display (default r"\\psi").

  • verbose (bool, optional) – Enable verbose output (default False).

  • bcs (list, optional) – Boundary conditions (default []).

  • order (int, optional) – Order of time integration (1-3) (default 1).

  • smoothing (float, optional) – Smoothing parameter (default 0.0).

Notes

The Symbolic class is the base for understanding BDF and Adams-Moulton formulas without the complexity of mesh or swarm storage. It is primarily useful for:

  • Understanding time-stepping algorithm behavior

  • Debugging symbolic expressions in time-dependent problems

  • Prototyping before implementing with mesh/swarm storage

For actual simulations, use Eulerian, SemiLagrangian, or Lagrangian which store history on computational meshes or swarms.

See also

Eulerian

Mesh-based history with BDF time-stepping.

SemiLagrangian

Nodal-swarm approach for advection-dominated problems.

Lagrangian

Swarm-based material tracking.

__init__(psi_fn, theta=0.5, varsymbol='\\psi', verbose=False, bcs=[], order=1, smoothing=0.0)[source]
Parameters:
property psi_fn

Current symbolic expression \(\psi\) being tracked.

property effective_order

Current effective BDF order, accounting for history startup.

For BDF order k, k distinct history values are needed. During startup, effective_order ramps from 1 to self.order as successive solves populate the history slots with distinct values.

update_history_fn()[source]

Copy current \(\psi\) to the first history slot psi_star[0].

initialise_history()[source]

Initialize all history slots to the current value of \(\psi\).

Called automatically on the first update_pre_solve. Can also be called manually after setting initial conditions.

initiate_history_fn()[source]

Deprecated: use initialise_history instead.

update(dt, evalf=False, verbose=False)[source]

Update history (alias for update_pre_solve).

Parameters:
  • evalf (bool | None)

  • verbose (bool | None)

update_pre_solve(dt, evalf=False, verbose=False)[source]

Pre-solve update hook. Auto-initialises history on first call.

Parameters:
  • evalf (bool | None)

  • verbose (bool | None)

update_post_solve(dt, evalf=False, verbose=False)[source]

Shift history chain after solve: \(\psi^{*n} \leftarrow \psi^{*(n-1)}\).

Parameters:
  • evalf (bool | None)

  • verbose (bool | None)

bdf(order=None)[source]

Compute the backward differentiation approximation of the time-derivative of ψ. For order 1: bdf ≡ ψ - psi_star[0]

When order is not specified, uses effective_order which ramps up from 1 during startup to avoid using higher-order BDF coefficients before distinct history values are available.

For order > 1 with variable timesteps, the BDF coefficients are adjusted using the ratio of consecutive timesteps (see _bdf_coefficients).

Parameters:

order (int | None)

adams_moulton_flux(order=None)[source]

Adams-Moulton flux approximation for implicit time integration.

When order is not specified, uses effective_order which ramps up from 1 during startup.

Parameters:

order (int, optional) – Order of the approximation (1-3). Defaults to effective_order.

Returns:

Weighted average of \(\psi\) and history terms.

Return type:

sympy.Matrix

Notes

The Adams-Moulton formulas for order 1-3 are:

  • Order 1: \(\theta \psi + (1-\theta) \psi^*\)

  • Order 2: \(\frac{5\psi + 8\psi^* - \psi^{**}}{12}\)

  • Order 3: \(\frac{9\psi + 19\psi^* - 5\psi^{**} + \psi^{***}}{24}\)

Eulerian Derivatives

Eulerian

class underworld3.systems.ddt.Eulerian[source]

Bases: uw_object

Eulerian (mesh-based) history manager for time derivatives.

Manages the update of a variable \(\psi\) on the mesh across timesteps, storing history values on mesh variables for backward differentiation.

\[\begin{split}\psi_p^{t-n\Delta t} &\leftarrow \psi_p^{t-(n-1)\Delta t} \\ \psi_p^{t-(n-1)\Delta t} &\leftarrow \psi_p^{t-(n-2)\Delta t} \cdots \\ \psi_p^{t-\Delta t} &\leftarrow \psi_p^{t}\end{split}\]

When V_fn is provided, the update_pre_solve method applies an explicit advection correction so that bdf() approximates the full material derivative \(D\psi/Dt = \partial\psi/\partial t + \mathbf{u}\cdot\nabla\psi\) rather than the partial time derivative alone.

Note

The optional advection capability (V_fn parameter) is suitable for problems where the advection is weak or where a purely grid-based approach is desired (e.g., Richards equation with no transport). For advection-dominated problems, SemiLagrangian is more mature and generally preferred.

Parameters:
  • mesh (underworld3.discretisation.Mesh) – The computational mesh.

  • psi_fn (MeshVariable or sympy.Basic) – The quantity to track. Can be a mesh variable or symbolic expression.

  • vtype (VarType) – Variable type (SCALAR, VECTOR, etc.) for history storage.

  • degree (int) – Polynomial degree for history mesh variables.

  • continuous (bool) – Whether history variables are continuous across element boundaries.

  • V_fn (sympy.Basic, optional) – Velocity field for grid-based advection correction. If None (default), computes pure ∂ψ/∂t. If set, computes D/Dt = ∂/∂t + u·∇ via operator splitting.

  • evalf (bool, default=False) – If True, evaluate expressions numerically during updates.

  • theta (float, default=0.5) – Time-stepping parameter for implicit/explicit blending. theta=0 is fully explicit, theta=1 is fully implicit.

  • varsymbol (str, default=r"u") – LaTeX symbol for display.

  • verbose (bool, default=False) – Enable verbose output during updates.

  • bcs (list, default=[]) – Boundary conditions to apply to projections.

  • order (int, default=1) – Number of history timesteps to store (for multi-step methods).

  • smoothing (float, default=0.0) – Smoothing parameter for projections.

See also

SemiLagrangian

For advection-dominated problems with nodal swarm.

Lagrangian

For full Lagrangian tracking on swarms.

Symbolic

For purely symbolic history (no mesh storage).

__init__(mesh, psi_fn, vtype, degree, continuous, V_fn=None, evalf=False, theta=0.5, varsymbol='u', verbose=False, bcs=[], order=1, smoothing=0.0)[source]
Parameters:
  • mesh (Mesh)

  • psi_fn (EnhancedMeshVariable | Basic)

  • vtype (VarType)

  • degree (int)

  • continuous (bool)

  • evalf (bool | None)

  • theta (float | None)

  • varsymbol (str | None)

  • verbose (bool | None)

property psi_fn

Current symbolic expression \(\psi\) being tracked.

property effective_order

Current effective BDF order, accounting for history startup.

For BDF order k, k distinct history values are needed. During startup, effective_order ramps from 1 to self.order as successive solves populate the history slots with distinct values.

update_history_fn()[source]

Copy current \(\psi\) to psi_star[0] via evaluation or projection.

initialise_history()[source]

Initialize all history slots to the current value of \(\psi\).

Called automatically on the first update_pre_solve. Can also be called manually after setting initial conditions to ensure the history chain starts from the correct state.

initiate_history_fn()[source]

Deprecated: use initialise_history instead.

update(dt, evalf=False, verbose=False)[source]

Update history (alias for update_pre_solve).

Parameters:
  • evalf (bool | None)

  • verbose (bool | None)

update_pre_solve(dt, evalf=False, verbose=False)[source]

Pre-solve: auto-initialise history and apply advection correction.

On the first call, automatically initialises history from the current field values. If V_fn is set, also applies an explicit grid-based advection correction so that bdf() approximates the material derivative Dφ/Dt rather than ∂φ/∂t.

Parameters:
  • evalf (bool | None)

  • verbose (bool | None)

update_post_solve(dt, evalf=False, verbose=False)[source]

Shift history chain after solve: \(\psi^{*n} \leftarrow \psi^{*(n-1)}\).

Parameters:
  • evalf (bool | None)

  • verbose (bool | None)

bdf(order=None)[source]

Backwards differentiation form for calculating DuDt.

Note that you will need bdf / \(\delta t\) in computing derivatives.

When order is not specified, uses effective_order which ramps up from 1 during startup to avoid using higher-order BDF coefficients before distinct history values are available.

For order > 1 with variable timesteps, the BDF coefficients are adjusted using the ratio of consecutive timesteps (see _bdf_coefficients).

adams_moulton_flux(order=None)[source]

Adams-Moulton flux approximation for implicit time integration.

When order is not specified, uses effective_order which ramps up from 1 during startup.

Parameters:

order (int, optional) – Order of the approximation (0-3). Defaults to effective_order.

Returns:

  • sympy.Basic – Weighted average of \(\psi\) and history terms.

  • .. warning:: – For order > 1, the coefficients assume uniform timestep spacing. Variable Δt with order > 1 requires modified coefficients. Use order=1 when timestep varies significantly.

Lagrangian Derivatives

SemiLagrangian

class underworld3.systems.ddt.SemiLagrangian[source]

Bases: uw_object

Semi-Lagrangian history manager using nodal swarm.

Manages the semi-Lagrangian update of a mesh variable \(\psi\) across timesteps. Uses a nodal swarm to track departure points and interpolate values back to the mesh.

\[\begin{split}\psi_p^{t-n\Delta t} &\leftarrow \psi_p^{t-(n-1)\Delta t} \\ \psi_p^{t-(n-1)\Delta t} &\leftarrow \psi_p^{t-(n-2)\Delta t} \cdots \\ \psi_p^{t-\Delta t} &\leftarrow \psi_p^{t}\end{split}\]

The semi-Lagrangian method traces characteristics backward in time to find departure points, providing stable advection without CFL restrictions while maintaining accuracy.

Parameters:
  • mesh (underworld3.discretisation.Mesh) – The computational mesh.

  • psi_fn (sympy.Function) – The quantity to advect (typically a mesh variable’s symbolic form).

  • V_fn (sympy.Function) – Velocity field for advection (e.g., stokes.u.sym).

  • vtype (VarType) – Variable type (SCALAR, VECTOR, SYM_TENSOR, etc.).

  • degree (int) – Polynomial degree for mesh variable storage.

  • continuous (bool) – Whether variables are continuous across element boundaries.

  • swarm_degree (int, optional) – Polynomial degree for swarm interpolation. Defaults to degree.

  • swarm_continuous (bool, optional) – Continuity for swarm variables. Defaults to continuous.

  • varsymbol (str, optional) – LaTeX symbol for display.

  • verbose (bool, default=False) – Enable verbose output during updates.

  • bcs (list, default=[]) – Boundary conditions for projections.

  • order (int, default=1) – Number of history timesteps (1 for first-order, 2 for second-order).

  • smoothing (float, default=0.0) – Smoothing parameter for projections.

  • preserve_moments (bool, default=False) – Use moment-preserving projection (experimental).

Notes

The semi-Lagrangian method is particularly useful for:

  • Advection-dominated problems (high Péclet number)

  • Problems where CFL stability is restrictive

  • Viscoelastic stress advection

See also

Eulerian

For fixed-mesh time derivatives without advection.

Lagrangian

For full particle-following Lagrangian tracking.

__init__(mesh, psi_fn, V_fn, vtype, degree, continuous, swarm_degree=None, swarm_continuous=None, varsymbol=None, verbose=False, bcs=[], order=1, smoothing=0.0, preserve_moments=False)[source]
Parameters:
property psi_fn

Current symbolic expression \(\psi\) being tracked.

property effective_order

Current effective BDF order, accounting for history startup.

For BDF order k, k distinct history values are needed. During startup, effective_order ramps from 1 to self.order as successive solves populate the history slots with distinct values.

initialise_history()[source]

Initialize all history slots to the current value of \(\psi\).

Called automatically on the first update_pre_solve. Can also be called manually after setting initial conditions.

initiate_history_fn()[source]

Deprecated: use initialise_history instead.

update(dt, evalf=False, verbose=False, dt_physical=None)[source]

Update history (alias for update_pre_solve).

Parameters:
update_post_solve(dt, evalf=False, verbose=False, dt_physical=None)[source]

Post-solve: record timestep and increment solve counter.

Parameters:
update_pre_solve(dt, evalf=False, verbose=False, dt_physical=None)[source]

Sample upstream values along characteristics before solve.

On the first call, automatically initialises history from the current field values so that bdf() returns zero on the first step.

Parameters:
bdf(order=None)[source]

Backwards differentiation form for calculating DuDt.

Note that you will need bdf / \(\delta t\) in computing derivatives.

When order is not specified, uses effective_order which ramps up from 1 during startup. For order > 1 with variable timesteps, coefficients are adjusted automatically.

adams_moulton_flux(order=None)[source]

Adams-Moulton flux approximation for implicit time integration.

When order is not specified, uses effective_order which ramps up from 1 during startup.

Parameters:

order (int, optional) – Order of the approximation (0-3). Defaults to effective_order.

Returns:

Weighted average of \(\psi\) and history terms.

Return type:

sympy.Basic

Lagrangian

class underworld3.systems.ddt.Lagrangian[source]

Bases: uw_object

Swarm-based Lagrangian history manager for material tracking.

Manages the update of a Lagrangian variable \(\psi\) on a swarm across timesteps. Creates and manages its own internal swarm for tracking material properties through the flow.

\[ \begin{align}\begin{aligned}\psi_p^{t-n\Delta t} \leftarrow \psi_p^{t-(n-1)\Delta t}\\\psi_p^{t-(n-1)\Delta t} \leftarrow \psi_p^{t-(n-2)\Delta t} \cdots\\\psi_p^{t-\Delta t} \leftarrow \psi_p^{t}\end{aligned}\end{align} \]

The Lagrangian approach follows material points through the flow, avoiding numerical diffusion in advection. History values are stored on swarm variables with proxy mesh variables for solver integration.

Parameters:
  • mesh (underworld3.discretisation.Mesh) – The computational mesh.

  • psi_fn (sympy.Function) – The quantity to track (e.g., stress tensor stokes.stress).

  • V_fn (sympy.Function) – Velocity field for particle advection.

  • vtype (VarType) – Variable type (SCALAR, VECTOR, SYM_TENSOR, etc.).

  • degree (int) – Polynomial degree for proxy mesh variables.

  • continuous (bool) – Whether proxy variables are continuous across elements.

  • varsymbol (str, default=r"u") – LaTeX symbol for display.

  • verbose (bool, default=False) – Enable verbose output during updates.

  • bcs (list, default=[]) – Boundary conditions (currently unused for swarm variables).

  • order (int, default=1) – Number of history timesteps to store.

  • smoothing (float, default=0.0) – Smoothing parameter for projections.

  • fill_param (int, default=3) – Fill parameter for swarm population density.

Notes

The Lagrangian method is ideal for:

  • Viscoelastic stress history tracking

  • Material property advection without diffusion

  • Tracking compositional fields

The internal swarm is automatically advected during updates.

See also

SemiLagrangian

Mesh-based semi-Lagrangian with departure points.

Lagrangian_Swarm

For user-provided swarms.

instances = 0
__init__(mesh, psi_fn, V_fn, vtype, degree, continuous, varsymbol='u', verbose=False, bcs=[], order=1, smoothing=0.0, fill_param=3)[source]
Parameters:
property effective_order

Current effective BDF order, accounting for history startup.

initialise_history()[source]

Initialize all history slots to the current value of \(\psi\).

Called automatically on the first update_pre_solve. Can also be called manually after setting initial conditions.

initiate_history_fn()[source]

Deprecated: use initialise_history instead.

update(dt, evalf=False, verbose=False)[source]

Update history (alias for update_post_solve).

Parameters:
update_pre_solve(dt, evalf=False, verbose=False)[source]

Pre-solve: auto-initialise history on first call.

Parameters:
update_post_solve(dt, evalf=False, verbose=False)[source]

Shift history chain and advect swarm after solve.

Parameters:
bdf(order=None)[source]

Backwards differentiation form for calculating DuDt.

Note that you will need bdf / \(\delta t\) in computing derivatives.

When order is not specified, uses effective_order which ramps up from 1 during startup. For order > 1 with variable timesteps, coefficients are adjusted automatically.

adams_moulton_flux(order=None)[source]

Adams-Moulton flux approximation for implicit time integration.

When order is not specified, uses effective_order which ramps up from 1 during startup.

Parameters:

order (int, optional) – Order of the approximation (0-3). Defaults to effective_order.

Returns:

Weighted average of \(\psi\) and history terms.

Return type:

sympy.Basic

Lagrangian_Swarm

class underworld3.systems.ddt.Lagrangian_Swarm[source]

Bases: uw_object

Swarm-based Lagrangian history manager (user-provided swarm).

Manages the update of a Lagrangian variable \(\psi\) on a user-supplied swarm across timesteps:

\[ \begin{align}\begin{aligned}\psi_p^{t-n\Delta t} \leftarrow \psi_p^{t-(n-1)\Delta t}\\\psi_p^{t-(n-1)\Delta t} \leftarrow \psi_p^{t-(n-2)\Delta t} \cdots\\\psi_p^{t-\Delta t} \leftarrow \psi_p^{t}\end{aligned}\end{align} \]

Unlike Lagrangian, this class uses a user-provided swarm rather than creating its own. The swarm should already be populated and configured for tracking material points.

Parameters:
  • swarm (underworld3.swarm.Swarm) – User-provided swarm for material point tracking.

  • psi_fn (sympy.Function) – The quantity to track (e.g., stress tensor from a solver).

  • vtype (underworld3.VarType) – Variable type (SCALAR, VECTOR, SYM_TENSOR, etc.).

  • degree (int) – Interpolation degree for proxy mesh variables.

  • continuous (bool) – Whether proxy mesh variables should be continuous.

  • varsymbol (str, optional) – LaTeX symbol for display (default "u").

  • verbose (bool, optional) – Enable verbose output (default False).

  • bcs (list, optional) – Boundary conditions (currently unused, default []).

  • order (int, optional) – Order of time integration (1-3) (default 1).

  • smoothing (float, optional) – Smoothing parameter (default 0.0).

  • step_averaging (int, optional) – Number of steps for history averaging (default 2).

mesh

Reference to the computational mesh (from swarm).

Type:

underworld3.discretisation.Mesh

swarm

The user-provided swarm.

Type:

underworld3.swarm.Swarm

psi_fn

Symbolic expression for the tracked quantity.

Type:

sympy.Function

order

Order of BDF integration.

Type:

int

step_averaging

Number of steps for averaging (affects BDF scaling).

Type:

int

psi_star

History values \(\psi^*, \psi^{**}, \ldots\) as swarm variables.

Type:

list

Notes

Key differences from Lagrangian class:

  • Uses user-provided swarm (not internally created)

  • Swarm advection is NOT performed (user’s responsibility)

  • Step averaging for smoothing history updates

  • Suitable when swarm is shared between multiple history managers

The step_averaging parameter scales the BDF formula to account for updates that occur over multiple sub-steps within a timestep.

See also

Lagrangian

Creates and manages its own swarm with advection.

SemiLagrangian

Nodal-swarm approach for advection-dominated problems.

Eulerian

Pure mesh-based history (no particle tracking).

instances = 0
__init__(swarm, psi_fn, vtype, degree, continuous, varsymbol='u', verbose=False, bcs=[], order=1, smoothing=0.0, step_averaging=2)[source]
Parameters:
property effective_order

Current effective BDF order, accounting for history startup.

initialise_history()[source]

Initialize all history slots to the current value of \(\psi\).

Called automatically on the first update_pre_solve. Can also be called manually after setting initial conditions.

initiate_history_fn()[source]

Deprecated: use initialise_history instead.

update(dt, evalf=False, verbose=False)[source]

Update history (alias for update_post_solve).

Parameters:
update_pre_solve(dt, evalf=False, verbose=False)[source]

Pre-solve: auto-initialise history on first call.

Parameters:
update_post_solve(dt, evalf=False, verbose=False)[source]

Shift history chain and evaluate current \(\psi\) on swarm.

Parameters:
bdf(order=None)[source]

Backwards differentiation form for calculating DuDt.

Note that you will need bdf / \(\delta t\) in computing derivatives.

When order is not specified, uses effective_order which ramps up from 1 during startup. For order > 1 with variable timesteps, coefficients are adjusted automatically.

adams_moulton_flux(order=None)[source]

Adams-Moulton flux approximation for implicit time integration.

When order is not specified, uses effective_order which ramps up from 1 during startup.

Parameters:

order (int, optional) – Order of the approximation (1-3). Defaults to effective_order.

Returns:

Weighted average of \(\psi\) and history terms.

Return type:

sympy.Basic

Convenience Aliases

The following aliases are available via underworld3.systems: