zoomy_core.model.models.ins_generator module

zoomy_core.model.models.ins_generator module#

General-purpose INS equation framework with composable projection operations.

Class hierarchy:

StateSpace — shared symbols (coordinates, fields, stress tensor, parameters) SymbolicBase — name + display ├── Expression — PDE terms with .project(), .ibp(), .apply(), .terms, [i] └── Relation — lhs = rhs substitution rules with .apply_to()

├── Assumption — physical conditions (kinematic BCs, hydrostatic) └── Material — constitutive models (Newtonian, inviscid)

FullINS(state) — builds INS equations from a StateSpace

Design: the user drives every step. No hidden assumptions.

class zoomy_core.model.models.ins_generator.StateSpace(dimension=2)#

Bases: object

Shared symbolic state: coordinates, velocity fields, pressure, stress tensor, bathymetry, and physical parameters.

dimension is the physical space dimension:

2 = xz plane (1D horizontal shallow water) 3 = xyz space (2D horizontal shallow water)

The vertical coordinate is always z. Horizontal coordinates are x (and y if dimension=3). The horizontal_dim property gives the number of horizontal directions (1 or 2).

property horizontal_dim#
property has_y#
class zoomy_core.model.models.ins_generator.SymbolicBase(name='')#

Bases: object

Base for all symbolic objects. Provides name and notebook display.

class zoomy_core.model.models.ins_generator.Expression(expr, name='', term_groups=None)#

Bases: SymbolicBase

Symbolic expression for PDE terms.

Supports: - Term access: .terms, [i], len(), iteration - Projection: .project(test, var, domain, …) - Integration by parts: .ibp(var, test_weight, domain) -> IBPResult - Apply conditions: .apply({old: new}, …) or .apply(relation) - SymPy: .subs(), .simplify(), .expand(), .doit() - Notebook: _repr_latex_

property terms#
apply_to_term(index, *operations)#

Apply operations to a specific term and return the full expression.

Usage:

xmom.apply_to_term(5, ProductRule())
xmom.apply_to_term(5, {old: new})
keep_groups(*group_names)#

Return a new Expression keeping only the specified term groups.

Usage:

zmom.keep_groups("pressure", "source")  # drops temporal, convection, stress
drop_groups(*group_names)#

Return a new Expression dropping the specified term groups.

Usage:

zmom.drop_groups("temporal", "convection", "stress")
project(test_func, var, domain=(0, 1), weight=1, scale=1, numerical=False, order=4)#

Galerkin projection: integral(expr * test * weight * scale, (var, a, b)).

ibp(var, test_weight, domain=(0, 1), scale=1)#

Integration by parts on the outermost Derivative w.r.t. var. Returns IBPResult(integrate, boundary_upper, boundary_lower).

apply(*conditions)#

Apply conditions: dicts, (old, new) tuples, or Relation objects. Relation objects use their .apply_to() method. Preserves term_groups if present.

depth_integrate(lower, upper, var, method='auto')#

Depth-integrate this expression over [lower, upper] w.r.t. var.

Parameters:
  • lower (sympy expressions for the integration bounds) – (typically b and b+H, both functions of t and x)

  • upper (sympy expressions for the integration bounds) – (typically b and b+H, both functions of t and x)

  • var (Symbol) – The vertical coordinate (z)

  • method (str) –

    ‘auto’ : detect derivative direction and choose method ‘leibniz’ : pull horizontal derivative outside:

    int df/dx dz = d/dx[int f dz] - f(upper)*d(upper)/dx + f(lower)*d(lower)/dx

    ’fundamental_theorem’for vertical derivatives:

    int df/dz dz = f(upper) - f(lower)

    ’direct’ : keep as Integral(expr, (var, lower, upper))

Returns:

  • DepthIntegralResult with .volume, .boundary_upper, .boundary_lower

  • or a plain Expression for ‘direct’ and ‘fundamental_theorem’.

subs(*args, **kwargs)#
simplify()#

Simplify: expand + cancel, preserving Integral and Derivative(Integral) terms.

expand()#
doit()#
has(*args)#
property free_symbols#
map(fn)#

Apply fn to each term, reassemble into a single Expression.

fn receives an Expression (single term) and must return either an Expression or a DepthIntegralResult. DepthIntegralResults are assembled (volume + boundaries) before summing.

Example

integrated = expr.map(lambda t: t.depth_integrate(b, eta, z))

map_with_bcs(fn, bcs)#

Like map(), but collects boundary terms and applies BCs globally.

This is the correct way to depth-integrate a full equation: boundary terms from ALL terms are combined first, then BCs are applied once (so cross-term cancellations happen properly).

Parameters:
  • fn (callable) – Applied to each term. Must return DepthIntegralResult or Expression.

  • bcs (list of Relation) – Kinematic BCs etc. applied to the combined boundary expression.

Returns:

The fully depth-integrated equation with BCs applied.

Return type:

Expression

classify(t=None, x=None, z=None)#

Classify each term by its role in the PDE.

Returns a dict: {role: Expression} where role is one of: ‘temporal’, ‘convective’, ‘diffusive’, ‘source’.

Detection rules: - Has d/dt → temporal - Has d/dx (first-order) of a product → convective flux - Has d²/dz² or d/dz of d/dz → diffusive - Otherwise → source (algebraic)

property temporal#

only temporal (d/dt) terms.

Type:

View

property convective#

only convective flux (d/dx) terms.

Type:

View

project_onto_basis(basis, level, field_map, z_var, lower=None, upper=None, test_mode=None)#

Project a depth-integrated equation onto a polynomial basis.

Replaces every Integral(f(u,...), (z, b, eta)) by substituting the basis expansion u(z) = sum alpha_k phi_k(zeta) and evaluating the resulting integrals using the SymbolicIntegrator.

If test_mode is an integer, multiplies each integral by the test function phi_{test_mode}(zeta) before evaluating (Galerkin projection for a specific mode). If test_mode=None, returns the scalar integral (e.g. for the mass equation where no test function is needed).

Parameters:
  • basis (Basisfunction class (e.g. Legendre_shifted)) –

  • level (int) –

  • field_map (dict) – Maps the original Function name to a list of SymPy Symbols for the basis coefficients. Example: {‘u’: [alpha_0, alpha_1, alpha_2]}

  • z_var (Symbol) – The vertical coordinate (z) that appears in the integrals.

  • lower (sympy expressions (optional)) – The integration bounds (b, eta). If None, detected from the first Integral found.

  • upper (sympy expressions (optional)) – The integration bounds (b, eta). If None, detected from the first Integral found.

  • test_mode (int or None) – If int, project onto test function phi_{test_mode}.

Returns:

With all depth integrals replaced by basis matrix products.

Return type:

Expression

latex(strip_args=False, multiline=False)#

LaTeX representation.

Parameters:
  • strip_args (bool) – u(t,x,z)u. Partial derivatives preserved.

  • multiline (bool) – Render as \begin{aligned} with one group per line (requires term_groups).

describe(header=True, final_equation=True, parameters=False, strip_args=False)#

Composable description of this expression.

Returns a Description that renders as markdown in Jupyter.

Parameters:
  • header (bool) – Show expression name + term count.

  • final_equation (bool) – Show the symbolic equation.

  • parameters (bool) – List free symbols.

  • strip_args (bool) – Display u instead of u(t, x, z).

class zoomy_core.model.models.ins_generator.DepthIntegralResult(volume, boundary_upper, boundary_lower)#

Bases: object

Result of depth-integrating a term over [lower, upper].

volume#

the volume integral (Expression)

boundary_upper#

boundary term at z=upper (Expression)

boundary_lower#

boundary term at z=lower (Expression)

The full integral = volume + boundary_upper - boundary_lower. Kinematic BCs can be applied to the boundary terms via .apply_bcs().

apply_bcs(bc_lower=None, bc_upper=None)#
assemble()#

Combine all terms: volume + boundary_upper + boundary_lower.

The sign convention is that each component already carries its correct sign. For Leibniz: upper = -f(eta)*d(eta)/dx, lower = +f(b)*db/dx. For fundamental theorem: upper = +f(eta), lower = -f(b).

class zoomy_core.model.models.ins_generator.IBPResult(integrate, boundary_upper, boundary_lower)#

Bases: object

Structured result from integration by parts.

integrate#

the volume integral (Expression)

boundary_upper#

boundary term at upper limit (Expression)

boundary_lower#

boundary term at lower limit (Expression)

apply_bcs(bc_lower=None, bc_upper=None)#
assemble()#
class zoomy_core.model.models.ins_generator.Relation(substitutions, name='')#

Bases: SymbolicBase

A symbolic relation: one or more substitution rules lhs_i = rhs_i.

Used as base for Assumption and Material. Can be applied to Expressions via expr.apply(relation), which calls relation.apply_to(expr).

Displays as a system of equations in notebooks.

apply_to(expr)#

Substitute all lhs -> rhs in the given SymPy expression.

class zoomy_core.model.models.ins_generator.Assumption(substitutions, name='')#

Bases: Relation

Physical assumption (kinematic BC, hydrostatic, etc.).

class zoomy_core.model.models.ins_generator.Material(substitutions, name='')#

Bases: Relation

Constitutive model (Newtonian, inviscid, etc.).

class zoomy_core.model.models.ins_generator.Operation(name='', description=None)#

Bases: SymbolicBase

An operation that transforms an Expression (e.g. depth integration).

Unlike a Relation (which substitutes symbols), an Operation applies a structural transformation. Works with both Expression.apply() and DerivedModel.apply().

Subclasses override apply_to(expr) which receives and returns a sympy expression, or apply_to_expression(expression) which receives and returns an Expression object.

apply_to(expr)#

Transform a sympy expression. Override for simple operations.

apply_to_expression(expression)#

Transform an Expression object. Override for operations that need access to .terms, .map(), etc.

class zoomy_core.model.models.ins_generator.DepthIntegrate(state)#

Bases: Operation

Depth-integrate all equations over [b, b+h] w.r.t. z.

Applies Leibniz rule and fundamental theorem term-by-term. Boundary values (w at z=b, u at z=b+h, etc.) remain as Subs objects. Apply ApplyKinematicBCs to see the cancellations.

apply_to_expression(eq)#

Transform an Expression object. Override for operations that need access to .terms, .map(), etc.

class zoomy_core.model.models.ins_generator.ApplyKinematicBCs(state)#

Bases: Operation

Apply kinematic BCs globally to combined boundary terms.

Evaluates Subs boundary values, applies kinematic BCs at surface and bottom, and simplifies. The Leibniz boundary u-terms cancel with the fundamental theorem w-terms.

Must be applied immediately after DepthIntegrate.

apply_to(expr)#

Transform a sympy expression. Override for simple operations.

class zoomy_core.model.models.ins_generator.SimplifyIntegrals(state)#

Bases: Operation

Evaluate integrals with constant integrand and remove zero integrals.

  • 0 dz 0

  • c dz c·h (if integrand has no z-dependence)

  • ∂/∂x 0 dz 0

apply_to(expr)#

Transform a sympy expression. Override for simple operations.

class zoomy_core.model.models.ins_generator.Integrate(var, lower, upper)#

Bases: Operation

Integrate the expression w.r.t. a variable over given bounds.

Sympy performs the integration analytically. For an equation like ∂p/∂z/ρ + g = 0, integrating from η to z gives p(z)/ρ - p(η)/ρ + g(z - η) = 0, i.e. the hydrostatic pressure.

Usage:

zmom.apply(Integrate(z, lower=eta, upper=z))
apply_to(expr)#

Transform a sympy expression. Override for simple operations.

class zoomy_core.model.models.ins_generator.ZetaTransform(state)#

Bases: Operation

Transform vertical coordinate: z = ζ·h + b, dz = h·dζ.

Transforms integrals: ∫_b^{b+h} f(z) dz → h·∫_0^1 f(ζ·h+b) dζ.

apply_to(expr)#

Transform a sympy expression. Override for simple operations.

class zoomy_core.model.models.ins_generator.ProductRule#

Bases: Operation

Inverse product rule: combine coeff · d(f)/dx into d(F)/dx.

Applied to a single term (use expr.terms[i].apply(ProductRule())). Detects the derivative in the term, checks if the remaining coefficient can be integrated w.r.t. the differentiated variable, and combines.

Example: g·h·∂h/∂x ∂/∂x(g·h²/2)

apply_to(expr)#

Transform a sympy expression. Override for simple operations.

zoomy_core.model.models.ins_generator.FullINS(state, equations=None)#

Build the incompressible Navier-Stokes as a System.

Returns a mutable system with .apply() for in-place transformations.

Parameters:
  • state (StateSpace) –

  • equations (list of str, optional) – Which equations to include. Default: all. Options: “continuity”, “x_momentum”, “y_momentum”, “z_momentum”.

  • Usage:: – ins = FullINS(state) ins.z_momentum.apply({state.w: 0}) ins.x_momentum.apply(HydrostaticPressure(state)) ins.describe()

class zoomy_core.model.models.ins_generator.Newtonian(state, nu=None)#

Bases: Material

Newtonian fluid: tau_ij = mu * (du_i/dx_j + du_j/dx_i), mu = rho * nu.

Parameters:

state (StateSpace) –

class zoomy_core.model.models.ins_generator.Inviscid(state)#

Bases: Material

Inviscid fluid: all tau_ij = 0.

Parameters:

state (StateSpace) –

class zoomy_core.model.models.ins_generator.materials#

Bases: object

Material model library. Usage: materials.newtonian(state)

newtonian#

alias of Newtonian

inviscid#

alias of Inviscid

class zoomy_core.model.models.ins_generator.KinematicBCBottom(state)#

Bases: Assumption

w|_{z=b} = db/dt + u_b * db/dx [+ v_b * db/dy]

Parameters:

state (StateSpace) –

class zoomy_core.model.models.ins_generator.KinematicBCSurface(state)#

Bases: Assumption

w|_{z=eta} = d(eta)/dt + u_s * d(eta)/dx [+ v_s * d(eta)/dy]

Parameters:

state (StateSpace) –

class zoomy_core.model.models.ins_generator.HydrostaticPressure(state)#

Bases: Assumption

p = p_atm + rho * g * (eta - z)

Parameters:

state (StateSpace) –

class zoomy_core.model.models.ins_generator.StressFreeSurface(state)#

Bases: Assumption

Stress-free surface: τ·n|_{z=η} = 0.

For a free surface with normal n ≈ (0,0,1), this gives:

τ_xz|_{z=η} = 0, τ_zz|_{z=η} = 0 (and τ_yz|_{z=η} = 0 in 3D)

Parameters:

state (StateSpace) –

class zoomy_core.model.models.ins_generator.ZeroAtmosphericPressure(state)#

Bases: Assumption

p_atm = 0 (no atmospheric pressure).

Parameters:

state (StateSpace) –

class zoomy_core.model.models.ins_generator.assumptions#

Bases: object

Assumptions library.

kinematic_bc_bottom#

alias of KinematicBCBottom

kinematic_bc_surface#

alias of KinematicBCSurface

hydrostatic_pressure#

alias of HydrostaticPressure

stress_free_surface#

alias of StressFreeSurface

zero_atmospheric_pressure#

alias of ZeroAtmosphericPressure

zoomy_core.model.models.ins_generator.integrate_by_parts(f, g, var, domain=(0, 1))#

Standalone IBP: integral d(f)/dvar * g dvar = [f*g]_a^b - integral f * dg/dvar dvar Returns IBPResult(integrate, boundary_upper, boundary_lower).

zoomy_core.model.models.ins_generator.gauss_legendre_integrate(expr, var, a, b, order=4)#