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.

For analytical-IC benchmarks (no startup transient) or checkpoint restarts, set_initial_history(values, dt=...) plants the BDF history directly and bypasses the order ramp, so the first solve runs at full BDF order.

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 state: DDtSymbolicState

Return a snapshot-of-state dataclass for this DDt instance.

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)

property bdf_coefficients

Current BDF coefficients [c0, c1, …] accounting for variable timesteps.

bdf(order=None)[source]

Backward differentiation approximation of the time-derivative of ψ.

Returns a fixed-structure symbolic expression using UWexpression coefficients. The coefficient values are updated each step in update_pre_solve — no JIT recompilation needed when the order ramps up or the timestep changes.

Parameters:

order (int, optional) – Ignored (kept for API compatibility). The effective order is controlled by the coefficient values.

adams_moulton_flux(order=None)[source]

Adams-Moulton flux approximation for implicit time integration.

Returns a fixed-structure symbolic expression using UWexpression coefficients. Values are updated each step in update_pre_solve.

Parameters:

order (int, optional) – Ignored (kept for API compatibility).

update_exp_coefficients(dt, tau_eff)[source]

Update the ETD-2 (exponential) coefficient values for this step.

Sets self._exp_coeffs[0].sym = α and self._exp_coeffs[1].sym = φ from current dt and tau_eff (Maxwell relaxation time \(\tau = \eta_\mathrm{eff}/\mu\)). Called by the constitutive model (which owns τ_eff) before each solve, peer to the BDF/AM coefficient updates that happen automatically in update_pre_solve.

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 state: DDtEulerianState

Return a snapshot-of-state dataclass for this Eulerian DDt.

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.

set_initial_history(values, dt=None)[source]

Plant history values for BDF restart or analytical IC.

Bypasses the automatic effective_order ramp so the very first solve runs at the full BDF order rather than starting at BDF-1. Use this when you have known values at \(t\) and past times — e.g. an analytical periodic solution, or a checkpointed history loaded from disk.

Parameters:
  • values (sequence of length self.order) – values[k] is \(\psi\) at \(t - k\,\Delta t\), i.e. values[0] is the current state. Each entry must be assignable to psi_star[k].array — either an array of matching shape, or a scalar that broadcasts.

  • dt (float, optional) – Uniform timestep assumed between history slots. Required for order >= 2 to seed correct multistep coefficients on the first solve. Ignored for order = 1.

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)

property bdf_coefficients

Current BDF coefficients [c0, c1, …] accounting for variable timesteps.

bdf(order=None)[source]

Backward differentiation approximation of the time-derivative of \(\psi\).

Returns a fixed-structure symbolic expression using UWexpression coefficients. Values are updated each step in update_pre_solve.

Parameters:

order (int, optional) – Ignored (kept for API compatibility).

adams_moulton_flux(order=None)[source]

Adams-Moulton flux approximation for implicit time integration.

Returns a fixed-structure symbolic expression using UWexpression coefficients. Values are updated each step in update_pre_solve.

Parameters:

order (int, optional) – Ignored (kept for API compatibility).

update_exp_coefficients(dt, tau_eff)[source]

Update the ETD-2 (exponential) coefficient values for this step.

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

  • with_forcing_history (bool, default=False) – When True, allocate an additional forcing_star MeshVariable (matching psi_star[0]’s shape, vtype, degree, continuity) to store one history slot for the strain-rate forcing. Required by ETD-2 exponential integration of the Maxwell relaxation operator; ignored for BDF/AM. Populated each step via update_forcing_history() (direct nodal evaluation of forcing_fn — typically the constitutive model’s strain-rate symbol).

  • theta (float, default=0.5) –

    Adams-Moulton θ for the implicit flux integrator at order 1. The order-1 AM coefficients are [θ, 1-θ]:

    • θ = 0.5 → Crank-Nicolson (trapezoidal, second-order accurate, A-stable). Default, matches legacy SLCN behaviour.

    • θ = 1.0 → Backward Euler (L-stable, monotone for diffusion, first-order accurate). Use for stiff parabolic terms (under-resolved sharp gradients on deformed cells) where CN’s lack of L-stability causes sign-flip ringing on stiff modes.

    Settable after construction as a property: adv_diff.DuDt.theta = 1.0.

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, with_forcing_history=False, monotone_mode=None, theta=0.5)[source]
Parameters:
  • mesh (Mesh)

  • psi_fn (Function)

  • V_fn (Function)

  • vtype (VarType)

  • degree (int)

  • continuous (bool)

  • swarm_degree (int | None)

  • swarm_continuous (bool | None)

  • varsymbol (str | None)

  • verbose (bool | None)

  • with_forcing_history (bool)

  • monotone_mode (str | None)

  • theta (float)

on_remesh(ctx)[source]

Adapt-time hook: ALE for the SL history stack.

Two branches:

  • Standard ALE (smooth adapt). The SL-owned vars (psi_star[i], forcing_star, _workVar, the flattening view _psi_star_flat_var) are CARRY + operator-managed — the generic per-variable pass already skipped them, and we leave their .data untouched here. Accumulate ctx.total_disp onto self._pending_v_mesh_disp so the next update_pre_solve() runs the SL trace-back along (V_fn v_mesh) with v_mesh = Δx / dt — that subtraction is exactly what compensates for the arbitrary mesh motion when reading the CARRY’d history. One-step pulse: the next solve consumes Δx and clears it.

  • Opt-out (e.g. OT_adapt reset). When the adapt is a discrete jump rather than a smooth displacement (ctx.scratch.get("ale_opt_out")), the linear Δx/dt v_mesh interpretation breaks down. Fall back to Phase-1 REMAP for this DDt’s managed vars: call remap_var_set() with the pre-move snapshot in ctx.managed_snapshot. The pending v_mesh is cleared because REMAP already brought the history onto the new positions.

Accumulation across multiple adapts before one solve is linear: v_mesh_disp += ctx.total_disp. The trace-back uses the SUM divided by the next dt, which is the correct node-frame velocity for that step.

property state: DDtSemiLagrangianState
property psi_fn

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

enable_source_snapshot()[source]

Enable snapshot substitution in the projection’s source field.

Call this once when the source expression (psi_fn) references psi_star[0] itself — without it the projection’s residual (target flux(psi_star[0]))·weight is implicit in the target because target and source share the same data field. With Min-mode plasticity at the yield kink, the implicit projection admits two fixed points (elastic and yield branches); under timestep change the iteration drifts to the elastic-branch fixed point and σ violates the yield surface.

The snapshot is a separate mesh variable matching psi_star[0]’s shape/vtype/degree. Each call to update_pre_solve copies psi_star[0].array psi_snapshot.array, freezing the source’s input for the upcoming projection. Substitution makes the projection’s compiled C code read from psi_snapshot.array instead of psi_star[0].array — there’s no recompile per step, just a memcpy.

Idempotent: safe to call more than once.

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.

set_initial_history(values, dt=None)[source]

Plant history values for BDF restart or analytical IC.

Bypasses the automatic effective_order ramp so the very first solve runs at the full BDF order rather than starting at BDF-1. Use this when you have known values at \(t\) and past times — e.g. an analytical periodic solution, or a checkpointed history loaded from disk.

Parameters:
  • values (sequence of length self.order) – values[k] is \(\psi\) at \(t - k\,\Delta t\), i.e. values[0] is the current state. Each entry must be assignable to psi_star[k].array — either an array of matching shape, or a scalar that broadcasts.

  • dt (float, optional) – Uniform timestep assumed between history slots. Required for order >= 2 to seed correct multistep coefficients on the first solve. Ignored for order = 1.

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, store_result=True, monotone_mode='__instance__')[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:
  • store_result (bool, optional) – If True (default), evaluate psi_fn at current positions and store in psi_star[0] before advection — the standard DDt behaviour. If False, skip this step and the history shift: only advect the existing psi_star levels upstream. Used by VE_Stokes where psi_star[0] already contains the projected actual stress from the previous solve.

  • monotone_mode (str or None, optional) – Override the instance self.monotone_mode for this call. Default "__instance__" (sentinel) means “use whatever is on the instance”. Pass None, "clamp", or "pick" to force a particular mode for one call.

  • dt (float)

  • evalf (bool | None)

  • verbose (bool | None)

  • dt_physical (float | None)

property bdf_coefficients

Current BDF coefficients [c0, c1, …] accounting for variable timesteps.

bdf(order=None)[source]

Backward differentiation approximation of the time-derivative of \(\psi\).

Returns a fixed-structure symbolic expression using UWexpression coefficients. Values are updated each step in update_pre_solve.

Parameters:

order (int, optional) – Ignored (kept for API compatibility).

adams_moulton_flux(order=None)[source]

Adams-Moulton flux approximation for implicit time integration.

Returns a fixed-structure symbolic expression using UWexpression coefficients. Values are updated each step in update_pre_solve.

Parameters:

order (int, optional) – Ignored (kept for API compatibility).

update_exp_coefficients(dt, tau_eff)[source]

Update the scalar ETD-2 (exponential) coefficient UWexpressions.

Sets self._exp_coeffs[0].sym = α = exp(-Δt/τ_eff) and self._exp_coeffs[1].sym = φ = (1-α)/(Δt/τ_eff) so the next solve uses the correct exponential coefficients via PetscDSSetConstants on the next _update_constants call.

update_forcing_history(forcing_fn=None, evalf=False, verbose=False)[source]

Refresh forcing_star from forcing_fn via direct nodal evaluation.

Used by ETD-2 exponential integration to store the current strain-rate field as \(\dot\varepsilon^{n}\) for the next step’s history term. Called by the constitutive model from the solver’s post-solve hook. Direct nodal evaluation (rather than an L2 projection through SNES) is sufficient because the strain rate is \(\nabla\mathbf{u}\), well-defined at nodes, with no history-coupled term that would make a projection implicit.

Unit handling

forcing_star is allocated with units=None (see __init__). When the model is unit-aware, forcing_fn is a symbolic expression of the velocity field whose evaluation returns a UnitAwareArray carrying strain-rate units (1/time). We non-dimensionalise that result via the active scaling system before assigning to forcing_star.array, which keeps the stored values consistent with the rest of the variable storage (codebase convention: variable storage is non-dimensional; units are re-attached at the .data interface). When the model is not unit-aware, the evaluation returns a plain ndarray and assignment is a straight numpy copy.

param forcing_fn:

Symbolic strain-rate field to evaluate at each node. If None, falls back to self._forcing_fn (set by the constitutive model at solve-attach time). No-op if neither is set or with_forcing_history=False.

type forcing_fn:

sympy expression, optional

param evalf:

Forwarded to uw.function.evaluate (forces numerical evaluation when True).

type evalf:

bool, optional

param verbose:

Enable verbose output.

type verbose:

bool, optional

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:
on_remesh(ctx)[source]

Adapt-time hook (Phase 1 no-op).

Lagrangian history is carried by the underlying swarm (each particle holds its own psi_star values), so there is no mesh-side state to transfer on an adapt — the swarm’s particle positions stay put and the new cells re-claim them. Method is defined so the registration shim has a target.

property state: DDtLagrangianState
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:
property bdf_coefficients

Current BDF coefficients [c0, c1, …] accounting for variable timesteps.

bdf(order=None)[source]

Backward differentiation approximation of the time-derivative of \(\psi\).

Returns a fixed-structure symbolic expression using UWexpression coefficients. Values are updated each step in update_pre_solve.

Parameters:

order (int, optional) – Ignored (kept for API compatibility).

adams_moulton_flux(order=None)[source]

Adams-Moulton flux approximation for implicit time integration.

Returns a fixed-structure symbolic expression using UWexpression coefficients. Values are updated each step in update_pre_solve.

Parameters:

order (int, optional) – Ignored (kept for API compatibility).

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 state: DDtLagrangianSwarmState
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:
property bdf_coefficients

Current BDF coefficients [c0, c1, …] accounting for variable timesteps.

bdf(order=None)[source]

Backward differentiation approximation of the time-derivative of \(\psi\).

Returns a fixed-structure symbolic expression using UWexpression coefficients. Values are updated each step in update_pre_solve.

Parameters:

order (int, optional) – Ignored (kept for API compatibility).

adams_moulton_flux(order=None)[source]

Adams-Moulton flux approximation for implicit time integration.

Returns a fixed-structure symbolic expression using UWexpression coefficients. Values are updated each step in update_pre_solve.

Parameters:

order (int, optional) – Ignored (kept for API compatibility).

Convenience Aliases

The following aliases are available via underworld3.systems: