tidal.symbolic package#
Symbolic computation layer for Lagrangian-to-PDE pipeline.
This package provides the Python-side interface for loading and processing field equations derived symbolically from Lagrangians via Mathematica/xAct.
- class tidal.symbolic.BoundaryCondition(type, value=None, derivative=None, gamma=None)[source]#
Bases:
objectBoundary condition for one spatial axis.
- classmethod from_dict(data)[source]#
Create a BoundaryCondition from a dictionary.
- Raises:
ValueError – If the BC type is not recognized.
- Parameters:
data (Mapping[str, Any])
- Return type:
- to_side_bc()[source]#
Convert to a
SideBCSpecfor the operator layer.- Raises:
ValueError – If the BC type is “periodic” (not representable as a side BC).
- Return type:
- class tidal.symbolic.ComponentEquation(field_name, field_index, time_derivative_order, rhs_terms, constraint_solver=<factory>)[source]#
Bases:
objectEquation of motion for a single field component.
- For a wave-type equation:
d^2/dt^2 field = sum of OperatorTerms
- Parameters:
field_name (str)
field_index (int)
time_derivative_order (int)
rhs_terms (tuple[OperatorTerm, ...])
constraint_solver (ConstraintSolverConfig)
- rhs_terms#
Terms on the RHS of the equation.
- Type:
tuple[OperatorTerm, …]
- constraint_solver#
Configuration for elliptic constraint solving. Only meaningful when
time_derivative_order == 0.- Type:
- classmethod from_dict(data, fields_lookup)[source]#
Create a ComponentEquation from a dictionary.
- Raises:
ValueError – If the RHS type is not “linear_combination”, or if constraint_solver is enabled for a non-constraint equation.
- Parameters:
- Return type:
- rhs_terms: tuple[OperatorTerm, ...]#
- constraint_solver: ConstraintSolverConfig#
- class tidal.symbolic.ConstraintSolverConfig(enabled=False, method='auto', boundary_conditions=<factory>, max_iterations=20, tolerance=1e-08)[source]#
Bases:
objectConfiguration for elliptic constraint solving.
When
enabledis True, the constraint equation is solved at each timestep rather than remaining frozen at its initial value.Solver Methods#
“auto” (default): Automatically selects the best solver. Uses FFT for fully periodic grids (O(N log N)), sparse matrix for non-periodic grids (O(N) via LU). Handles all constraint types: Poisson, Helmholtz, algebraic, anisotropic, etc.
“fft”: Force FFT solver. Requires fully periodic grid.
“matrix”: Force sparse matrix solver. Works with any BCs.
“poisson”: Original py-pde Poisson solver. Requires exactly one
laplacian(self_field)term with no other self-referencing operators. Will warn if non-laplacian self-terms are present. Backward compatible.
Coupled Constraint Parameters#
When multiple constraints reference each other’s fields, the solver iterates using Gauss-Seidel until convergence or
max_iterations. For fully periodic grids, coupled constraints are solved exactly via Fourier-space block solve (no iteration needed).- enabled#
Whether to solve the constraint elliptically. Default False preserves existing frozen-constraint behavior.
- Type:
- boundary_conditions#
Per-axis boundary conditions (e.g.,
{"x": ..., "y": ...}).- Type:
- tolerance#
Convergence threshold for coupled constraint iteration. Iteration stops when max|field_new - field_old| < tolerance (scaled by field magnitude for robustness). Must be > 0.
- Type:
- classmethod from_dict(data)[source]#
Create from a dictionary or return default (disabled).
- Parameters:
data (Mapping[str, Any] | None) – Parsed
constraint_solverblock from JSON, or None.- Returns:
Configuration instance.
- Return type:
- Raises:
ValueError – If
methodis not one of the recognized solver methods.
- boundary_conditions: dict[str, BoundaryCondition]#
- class tidal.symbolic.EquationSystem(n_components, dimension, spatial_dimension, component_names, equations, mass_matrix, coupling_matrix, metadata, coordinates=(), mass_matrix_symbolic=(), coupling_matrix_symbolic=(), canonical=None)[source]#
Bases:
objectComplete system of field equations derived from a Lagrangian.
- Parameters:
n_components (int)
dimension (int)
spatial_dimension (int)
equations (tuple[ComponentEquation, ...])
coupling_matrix_symbolic (tuple[tuple[str | None, ...], ...])
canonical (CanonicalStructure | None)
- equations#
Equations for each component.
- Type:
- coordinates#
Coordinate names from JSON spacetime.coordinates (e.g., (“t”, “x”, “y”)). Defaults to empty tuple; use
effective_coordinatesfor a guaranteed non-empty result that infers names from dimension when not set.
- canonical#
Canonical momentum and Hamiltonian structure from Legendre transform. Present when the JSON spec includes a
"canonical"section (generated bytidal derivefor non-linearization theories). None for legacy specs or linearization theories.- Type:
CanonicalStructure | None
- canonical: CanonicalStructure | None = None#
- property effective_coordinates: tuple[str, ...]#
Coordinate names, inferred from dimension if not set explicitly.
- property equation_map: dict[str, int]#
Map from field name to equation index. Cached on frozen dataclass.
- classmethod from_dict(data)[source]#
Create an EquationSystem from a dictionary (parsed JSON).
- Raises:
ValueError – If the JSON data is invalid or component references are inconsistent.
TypeError – If a metadata parameter has an unsupported type.
- Parameters:
data (Mapping[str, Any])
- Return type:
- property has_constraint_velocity_terms: bool#
True if the Hamiltonian has kinetic coupling between constraint fields.
Detects
time_derivative(C_i) x time_derivative(C_j)terms where both C_i and C_j are constraint fields (time_derivative_order == 0). These indicate the naive Hamiltonian H = sum(pi * v) - L is not the correct conserved quantity for this theory (Dirac-Bergmann theory).See GitHub issue #178 for details.
- property spatial_coordinates: tuple[str, ...]#
Spatial coordinate names (all except first, which is time).
- property state_layout: tuple[tuple[str, str], ...]#
State vector layout as (field_name, slot_type) tuples.
slot_type is “field” or “momentum”. Second-order components produce two entries (field, momentum); first-order/constraint produce one (field).
- property state_size: int#
Total number of state fields.
Second-order components contribute 2 slots (field + momentum). First-order and constraint components contribute 1 slot (field only).
- equations: tuple[ComponentEquation, ...]#
- class tidal.symbolic.OperatorTerm(coefficient, operator, field, coefficient_symbolic=None, time_dependent=False, coordinate_dependent=())[source]#
Bases:
objectA single term in the RHS of a field equation.
Represents: coefficient * operator(field)
- Parameters:
- coefficient_symbolic#
Optional symbolic name for the coefficient (e.g., “m2”, “-kappa”). When present, the coefficient can be overridden at runtime by passing a parameters dict to the PDE constructor.
- Type:
str | None
- time_dependent#
Whether the coefficient depends on time. For curved spacetime, terms like -2H∂_t φ (Hubble friction) have time-dependent coefficients when the conformal factor Ω(t) varies with time. Default False for flat spacetime.
- Type:
- coordinate_dependent#
Coordinate names the coefficient depends on (e.g., (“x”, “y”) for position-dependent coefficients on curved spatial surfaces, or (“t”,) for time-dependent). Empty tuple for constant coefficients.
- classmethod from_dict(data)[source]#
Create an OperatorTerm from a dictionary.
- Raises:
ValueError – If required keys are missing or operator is unknown.
- Parameters:
data (Mapping[str, Any])
- Return type:
- property position_dependent: bool#
Whether the coefficient depends on spatial coordinates.
Returns
Truewhencoordinate_dependentis non-empty (explicit declaration), or whencoefficient_symboliccontains a spatial coordinate call pattern such asx[]ory[](auto-detection for JSON exports that predate thecoordinate_dependentfield). Time-only dependence (t[]) returnsFalse.
- tidal.symbolic.load_equation_system(json_path)[source]#
Load an equation system from a JSON file.
- Parameters:
json_path (Path | str) – Path to the JSON file exported from Mathematica.
- Returns:
The parsed equation system.
- Return type:
- Raises:
FileNotFoundError – If the JSON file does not exist.
Submodules#
tidal.symbolic.json_loader module#
Load and validate JSON equation specifications from Mathematica/xAct export.
This module provides the data structures and parsing logic for loading field equations that were derived symbolically from Lagrangians.
- tidal.symbolic.json_loader.AXIS_LETTERS: tuple[str, ...] = ('x', 'y', 'z', 'w', 'v', 'u')#
Canonical spatial axis letters, ordered by dimension index. Supports up to 6 spatial dimensions (sufficient for all foreseeable physics).
- tidal.symbolic.json_loader.is_known_operator(name)[source]#
Check whether an operator name is recognized.
Accepts static operators (identity, laplacian, gradient_x, …), user-registered custom operators (via
register_operator), dynamic patterns for generic Nth-order derivatives (derivative_3_x, derivative_5_y, derivative_2x_1y, …), mixed time-space derivative operators (mixed_T2_S2x, mixed_T_S1x, …), and pure higher-order time operators on RHS (d2_t, d3_t, d4_t).
- class tidal.symbolic.json_loader.LHSStructure(expression, time_order, space_order=0)[source]#
Bases:
objectStructure describing the left-hand side of a PDE.
Supports different PDE types: - Elliptic (time_order=0): ∇²φ = f (Poisson, Laplace) - Parabolic (time_order=1): ∂_t φ = … (heat, diffusion) - Hyperbolic (time_order=2): ∂²_t φ = … (wave) - Higher order: ∂^n_t φ = …
- classmethod from_dict(data)[source]#
Create LHSStructure from structured JSON data.
Expected format: {“expression”: “…”, “order”: {“time”: N, “space”: 0}}
- Parameters:
data (Mapping[str, Any]) – The structured LHS data from JSON.
- Returns:
Parsed LHS structure.
- Return type:
- Raises:
ValueError – If
orderdict does not contain a'time'key.
- class tidal.symbolic.json_loader.OperatorTerm(coefficient, operator, field, coefficient_symbolic=None, time_dependent=False, coordinate_dependent=())[source]#
Bases:
objectA single term in the RHS of a field equation.
Represents: coefficient * operator(field)
- Parameters:
- coefficient_symbolic#
Optional symbolic name for the coefficient (e.g., “m2”, “-kappa”). When present, the coefficient can be overridden at runtime by passing a parameters dict to the PDE constructor.
- Type:
str | None
- time_dependent#
Whether the coefficient depends on time. For curved spacetime, terms like -2H∂_t φ (Hubble friction) have time-dependent coefficients when the conformal factor Ω(t) varies with time. Default False for flat spacetime.
- Type:
- coordinate_dependent#
Coordinate names the coefficient depends on (e.g., (“x”, “y”) for position-dependent coefficients on curved spatial surfaces, or (“t”,) for time-dependent). Empty tuple for constant coefficients.
- property position_dependent: bool#
Whether the coefficient depends on spatial coordinates.
Returns
Truewhencoordinate_dependentis non-empty (explicit declaration), or whencoefficient_symboliccontains a spatial coordinate call pattern such asx[]ory[](auto-detection for JSON exports that predate thecoordinate_dependentfield). Time-only dependence (t[]) returnsFalse.
- classmethod from_dict(data)[source]#
Create an OperatorTerm from a dictionary.
- Raises:
ValueError – If required keys are missing or operator is unknown.
- Parameters:
data (Mapping[str, Any])
- Return type:
- class tidal.symbolic.json_loader.BoundaryCondition(type, value=None, derivative=None, gamma=None)[source]#
Bases:
objectBoundary condition for one spatial axis.
- classmethod from_dict(data)[source]#
Create a BoundaryCondition from a dictionary.
- Raises:
ValueError – If the BC type is not recognized.
- Parameters:
data (Mapping[str, Any])
- Return type:
- to_side_bc()[source]#
Convert to a
SideBCSpecfor the operator layer.- Raises:
ValueError – If the BC type is “periodic” (not representable as a side BC).
- Return type:
- class tidal.symbolic.json_loader.ConstraintSolverConfig(enabled=False, method='auto', boundary_conditions=<factory>, max_iterations=20, tolerance=1e-08)[source]#
Bases:
objectConfiguration for elliptic constraint solving.
When
enabledis True, the constraint equation is solved at each timestep rather than remaining frozen at its initial value.Solver Methods#
“auto” (default): Automatically selects the best solver. Uses FFT for fully periodic grids (O(N log N)), sparse matrix for non-periodic grids (O(N) via LU). Handles all constraint types: Poisson, Helmholtz, algebraic, anisotropic, etc.
“fft”: Force FFT solver. Requires fully periodic grid.
“matrix”: Force sparse matrix solver. Works with any BCs.
“poisson”: Original py-pde Poisson solver. Requires exactly one
laplacian(self_field)term with no other self-referencing operators. Will warn if non-laplacian self-terms are present. Backward compatible.
Coupled Constraint Parameters#
When multiple constraints reference each other’s fields, the solver iterates using Gauss-Seidel until convergence or
max_iterations. For fully periodic grids, coupled constraints are solved exactly via Fourier-space block solve (no iteration needed).- enabled#
Whether to solve the constraint elliptically. Default False preserves existing frozen-constraint behavior.
- Type:
- boundary_conditions#
Per-axis boundary conditions (e.g.,
{"x": ..., "y": ...}).- Type:
- tolerance#
Convergence threshold for coupled constraint iteration. Iteration stops when max|field_new - field_old| < tolerance (scaled by field magnitude for robustness). Must be > 0.
- Type:
- boundary_conditions: dict[str, BoundaryCondition]#
- classmethod from_dict(data)[source]#
Create from a dictionary or return default (disabled).
- Parameters:
data (Mapping[str, Any] | None) – Parsed
constraint_solverblock from JSON, or None.- Returns:
Configuration instance.
- Return type:
- Raises:
ValueError – If
methodis not one of the recognized solver methods.
- class tidal.symbolic.json_loader.HamiltonianFactor(field, operator)[source]#
Bases:
objectOne factor in a quadratic Hamiltonian term.
Represents a field (or its time/spatial derivative) that appears as a multiplicative factor in a term of the component-form Hamiltonian density.
- operator#
Differential operator: “identity” (bare field), “time_derivative” (∂_t field, maps to canonical momentum at evaluation), or a spatial operator (“gradient_x”, “laplacian”, etc.).
- Type:
- class tidal.symbolic.json_loader.HamiltonianTerm(coefficient, factor_a, factor_b, coefficient_symbolic=None, coordinate_dependent=(), term_class='unknown')[source]#
Bases:
objectA single quadratic term in the Hamiltonian density.
H = Σ coefficient * factor_a * factor_b
- Parameters:
coefficient (float)
factor_a (HamiltonianFactor)
factor_b (HamiltonianFactor)
coefficient_symbolic (str | None)
term_class (str)
- factor_a#
First field factor.
- Type:
- factor_b#
Second field factor (may equal factor_a for squared terms).
- Type:
- coordinate_dependent#
Spatial axes the coefficient depends on (e.g.
("x", "y")). When non-empty, the coefficient must be evaluated on the grid. Older JSON exports omit this field; auto-detection viaposition_dependentcovers those cases.
- term_class#
Classification:
"self"(both factors reference the same base field) or"interaction"(cross-field coupling). Defaults to"unknown"for older JSONs; useis_self_energyproperty which auto-classifies by comparing factor field names.- Type:
- factor_a: HamiltonianFactor#
- factor_b: HamiltonianFactor#
- property position_dependent: bool#
True if the coefficient is a function of spatial coordinates.
Returns
Truewhencoordinate_dependentis non-empty (explicit declaration), or whencoefficient_symboliccontains a coordinate call pattern such asx[]ory[](auto-detection for JSON exports that predate thecoordinate_dependentfield).
- property is_self_energy: bool#
True if both factors reference the same base field (self-energy).
Uses
term_classwhen available (from Wolfram export), otherwise classifies by comparing base field names (strippingv_prefix).
- class tidal.symbolic.json_loader.CanonicalStructure(hamiltonian_terms, volume_element=None)[source]#
Bases:
objectCanonical structure: Hamiltonian terms for energy measurement.
The Hamiltonian’s bilinear terms are used by energy measurement to compute H(q, v). The E-L equations are stored in the
equationsarray directly.- Parameters:
hamiltonian_terms (tuple[HamiltonianTerm, ...])
volume_element (str | None)
- hamiltonian_terms#
Quadratic terms in the component-form Hamiltonian density. Used by energy measurement to compute H(q, v).
- Type:
tuple[HamiltonianTerm, …]
- volume_element#
Symbolic expression for
sqrt|det(g_spatial)|, the spatial volume element.Nonefor flat (Minkowski) spacetimes where the volume element is 1. Used by energy measurement to weight the Hamiltonian density before spatial integration.- Type:
str or None
- hamiltonian_terms: tuple[HamiltonianTerm, ...]#
- class tidal.symbolic.json_loader.ComponentEquation(field_name, field_index, time_derivative_order, rhs_terms, constraint_solver=<factory>)[source]#
Bases:
objectEquation of motion for a single field component.
- For a wave-type equation:
d^2/dt^2 field = sum of OperatorTerms
- Parameters:
field_name (str)
field_index (int)
time_derivative_order (int)
rhs_terms (tuple[OperatorTerm, ...])
constraint_solver (ConstraintSolverConfig)
- rhs_terms#
Terms on the RHS of the equation.
- Type:
tuple[OperatorTerm, …]
- constraint_solver#
Configuration for elliptic constraint solving. Only meaningful when
time_derivative_order == 0.- Type:
- rhs_terms: tuple[OperatorTerm, ...]#
- constraint_solver: ConstraintSolverConfig#
- classmethod from_dict(data, fields_lookup)[source]#
Create a ComponentEquation from a dictionary.
- Raises:
ValueError – If the RHS type is not “linear_combination”, or if constraint_solver is enabled for a non-constraint equation.
- Parameters:
- Return type:
- class tidal.symbolic.json_loader.EquationSystem(n_components, dimension, spatial_dimension, component_names, equations, mass_matrix, coupling_matrix, metadata, coordinates=(), mass_matrix_symbolic=(), coupling_matrix_symbolic=(), canonical=None)[source]#
Bases:
objectComplete system of field equations derived from a Lagrangian.
- Parameters:
n_components (int)
dimension (int)
spatial_dimension (int)
equations (tuple[ComponentEquation, ...])
coupling_matrix_symbolic (tuple[tuple[str | None, ...], ...])
canonical (CanonicalStructure | None)
- equations#
Equations for each component.
- Type:
- coordinates#
Coordinate names from JSON spacetime.coordinates (e.g., (“t”, “x”, “y”)). Defaults to empty tuple; use
effective_coordinatesfor a guaranteed non-empty result that infers names from dimension when not set.
- canonical#
Canonical momentum and Hamiltonian structure from Legendre transform. Present when the JSON spec includes a
"canonical"section (generated bytidal derivefor non-linearization theories). None for legacy specs or linearization theories.- Type:
CanonicalStructure | None
- equations: tuple[ComponentEquation, ...]#
- canonical: CanonicalStructure | None = None#
- property state_size: int#
Total number of state fields.
Second-order components contribute 2 slots (field + momentum). First-order and constraint components contribute 1 slot (field only).
- property state_layout: tuple[tuple[str, str], ...]#
State vector layout as (field_name, slot_type) tuples.
slot_type is “field” or “momentum”. Second-order components produce two entries (field, momentum); first-order/constraint produce one (field).
- property effective_coordinates: tuple[str, ...]#
Coordinate names, inferred from dimension if not set explicitly.
- property spatial_coordinates: tuple[str, ...]#
Spatial coordinate names (all except first, which is time).
- property has_constraint_velocity_terms: bool#
True if the Hamiltonian has kinetic coupling between constraint fields.
Detects
time_derivative(C_i) x time_derivative(C_j)terms where both C_i and C_j are constraint fields (time_derivative_order == 0). These indicate the naive Hamiltonian H = sum(pi * v) - L is not the correct conserved quantity for this theory (Dirac-Bergmann theory).See GitHub issue #178 for details.
- property equation_map: dict[str, int]#
Map from field name to equation index. Cached on frozen dataclass.
- classmethod from_dict(data)[source]#
Create an EquationSystem from a dictionary (parsed JSON).
- Raises:
ValueError – If the JSON data is invalid or component references are inconsistent.
TypeError – If a metadata parameter has an unsupported type.
- Parameters:
data (Mapping[str, Any])
- Return type:
- tidal.symbolic.json_loader.validate_json_schema(data)[source]#
Validate that the JSON data matches the expected schema.
- Raises:
ValueError – If required fields are missing or have invalid types.
- Parameters:
data (Mapping[str, Any])
- Return type:
None
- tidal.symbolic.json_loader.load_equation_system(json_path)[source]#
Load an equation system from a JSON file.
- Parameters:
json_path (Path | str) – Path to the JSON file exported from Mathematica.
- Returns:
The parsed equation system.
- Return type:
- Raises:
FileNotFoundError – If the JSON file does not exist.
tidal.symbolic.latex module#
Convert JSON equation specifications to LaTeX math notation.
Provides functions to render TIDAL equation systems as publication-ready LaTeX, including:
Component PDEs with proper operator notation
Lagrangian expressions with tensor index notation (
\\tensor{}package)Hamiltonian density terms
Symbolic coefficients (Mathematica InputForm → LaTeX)
Primary public entry point:
system_to_latex(spec, ...)— full equation system
- tidal.symbolic.latex.coefficient_to_latex(expr)[source]#
Convert a Mathematica-style symbolic coefficient to LaTeX.
Examples
>>> coefficient_to_latex("-(B0^2*kappa^2)") '-B_0^{2} \\\\kappa^{2}' >>> coefficient_to_latex("1/2") '\\\\frac{1}{2}'
- tidal.symbolic.latex.field_to_latex(name, *, tensor_meta=None, coordinates=())[source]#
Convert a field component name to LaTeX.
- Parameters:
name (str) – Field name (e.g., “h_5”, “phi_0”, “v_phi_0”).
tensor_meta (dict, optional) – Tensor metadata from enriched JSON:
{"tensor_head": "h", "tensor_rank": 2, "tensor_indices": [2, 2]}.coordinates (tuple[str, ...], optional) – Coordinate names for resolving index labels (e.g., (“t”, “x”, “y”, “z”)).
- Returns:
LaTeX string for the field.
- Return type:
- tidal.symbolic.latex.operator_to_latex(operator, field_latex)[source]#
Render an operator applied to a field in LaTeX.
- tidal.symbolic.latex.equation_to_latex(eq, spec)[source]#
Convert a single component equation to LaTeX.
- Parameters:
eq (ComponentEquation) – The equation to render.
spec (EquationSystem) – The parent equation system (for coordinates and tensor metadata).
- Returns:
LaTeX string (without environment wrapping).
- Return type:
- tidal.symbolic.latex.hamiltonian_to_latex(terms, spec)[source]#
Render the Hamiltonian density as a LaTeX equation.
- Returns:
LaTeX for
\\mathcal{H} = ....- Return type:
- Parameters:
terms (list[HamiltonianTerm])
spec (EquationSystem)
- tidal.symbolic.latex.lagrangian_to_latex(expr)[source]#
Convert a Lagrangian expression from xAct notation to LaTeX.
This is a best-effort conversion. The xAct abstract index notation is rich and idiosyncratic; this handles the patterns found in the 33 example JSONs in this project.
- tidal.symbolic.latex.system_to_latex(spec, *, output_format='align', include_hamiltonian=True, include_lagrangian=True)[source]#
Convert a full equation system to LaTeX.
- Parameters:
spec (EquationSystem) – The equation system to render.
output_format ({"align", "document", "raw"}) – Output format.
include_hamiltonian (bool) – Whether to include the Hamiltonian density.
include_lagrangian (bool) – Whether to include the Lagrangian expression.
- Returns:
LaTeX output.
- Return type:
tidal.symbolic.ostrogradsky module#
Ostrogradsky reduction for higher-derivative field equations.
Converts 4th-order-in-time equations to 2nd-order by introducing
auxiliary fields. For d⁴h/dt⁴ = RHS, defines w = d²h/dt²
and rewrites as two coupled 2nd-order equations:
d²h/dt² = w (kinematic) d²w/dt² = RHS’ (dynamics, with substituted operators)
This is the standard Ostrogradsky reduction (Ostrogradsky 1850, “Mémoires sur les équations différentielles”). Modern review: Woodard (2015, arXiv:1506.02210).
Physics warning: Non-degenerate 4th-order Lagrangians generically have Ostrogradsky ghosts (unbounded Hamiltonian). Ghost-free parameter windows exist for PGT (Barker 2024, arXiv:2406.12826). Full ghost detection is tracked in issue #164.
The reduction is applied in-memory during JSON loading, transparent to all solvers. The JSON on disk is unchanged.
- tidal.symbolic.ostrogradsky.MAX_NATIVE_TIME_ORDER = 2#
Maximum time-derivative order that solvers support natively. Equations with higher order are reduced via Ostrogradsky.
- class tidal.symbolic.ostrogradsky.AuxiliaryField(original_field, auxiliary_name, original_time_order)[source]#
Bases:
objectMetadata for an Ostrogradsky auxiliary field.
- tidal.symbolic.ostrogradsky.apply_ostrogradsky_reduction(spec)[source]#
Apply Ostrogradsky reduction to higher-derivative equations.
For each equation with
time_derivative_order > 2, introduces auxiliary fields to reduce to 2nd-order. The reduction is applied in-memory — the JSON on disk is unchanged.- Parameters:
spec (EquationSystem) – Original equation system, possibly with higher-order equations.
- Returns:
Modified system where all equations have
time_order <= 2. Auxiliary fields are added to the component list.- Return type:
- Raises:
ValueError – If equations have odd time-derivative order (not physically meaningful for Euler-Lagrange equations from a Lagrangian).
References
Ostrogradsky M.V. (1850), Mem. Acad. St. Petersbourg VI 4, 385.
Woodard R.P. (2015), arXiv:1506.02210.
Barker W.E.V. et al. (2024), arXiv:2406.12826.
tidal.symbolic.reduction module#
Plane-wave dimensional reduction for JSON equation specs.
After the Wolfram pipeline exports a JSON spec with the original higher-dimensional metadata (dimension, coordinates, operator names), this module transforms it into a clean reduced-dimension spec suitable for efficient simulation.
For example, a 3+1D theory reduced along z produces a 1+1D spec with:
- dimension: 2, coordinates: ["t", "x"], signature: [-1, 1]
- Operators remapped: laplacian_z → laplacian_x, gradient_z → gradient_x
- Killed-axis operators (laplacian_x, gradient_y, etc.) removed
- coordinate_dependent arrays and coefficient_symbolic strings updated
- Provenance metadata recording the reduction
Curved-coordinate support: - Coefficients depending only on the surviving coordinate are preserved and remapped - Volume element is kept if it depends only on the surviving coordinate - If any surviving coefficient or volume element references a killed coordinate,
ValueErroris raised (incompatible reduction)
- tidal.symbolic.reduction.reduce_spec(spec_data, reduction_config)[source]#
Apply plane-wave dimensional reduction to a JSON equation spec.
Transforms a higher-dimensional spec into a clean 1+1D spec by:
Remapping operator names (
laplacian_z → laplacian_x)Removing terms with killed-axis operators
Updating spacetime metadata (dimension, signature, coordinates)
Remapping coordinate references in coefficient expressions
Handling volume element (keep if surviving-coord-only, error if not)
Adding provenance metadata