Integral operators¤
Note
For a more detailed mathematical guide (measures, sampling, quadrature), see Guides → Integrals and measures.
Note
When integrating over a CoordSeparableBatch, Phydrax uses per-axis quadrature weights
from batch.axis_discretization_by_axis when available (e.g. Gauss–Legendre weights
from LegendreAxisSpec). Otherwise it falls back to uniform per-axis weights based on
the factor's axis-aligned bounding box.
phydrax.operators.build_quadrature(component: DomainComponent, batch: PointsBatch) -> QuadratureBatch
¤
Build a uniform tensor-product quadrature for a PointsBatch.
For each sampling axis \(a\) with \(n_a\) points, this constructs weights
where \(\mu_a\) is the measure associated with the labels in the corresponding
sampling block (as determined by component). The full weight is the product
\(w=\prod_a w_a\).
Arguments:
component: TheDomainComponentwhose measure is being integrated.batch: APointsBatchsampled with a canonicalizedProductStructure.
Returns:
- A
QuadratureBatchwith per-axis weights compatible withbatch.
phydrax.operators.build_ball_quadrature(*, radius: float, dim: int, num_points: int = 2048, seed: int = 0, method: typing.Literal['fibonacci', 'grid'] = 'fibonacci', angular_design: str | None = None, radial_strata: int | list[float] | None = None) -> dict
¤
Precompute isotropic quadrature offsets for a ball.
Constructs a rule on the Euclidean ball
returning a dict with:
- "offsets": array of shape (N, dim) containing offsets \(\xi_i\);
- "weights": array of shape (N,) with equal weights
\(w_i = |B_r(0)|/N\).
This is intended for local/nonlocal operators of the form:
phydrax.operators.integral._batch_ops.integral(f: DomainFunction | ArrayLike | None, batch: PointsBatch | CoordSeparableBatch | tuple[PointsBatch, ...], /, *, component: DomainComponent | DomainComponentUnion, quadrature: QuadratureBatch | tuple[QuadratureBatch | None, ...] | None = None, over: str | tuple[str, ...] | None = None, key: Key[Array, ''] = jr.key(0), **kwargs: Any) -> cx.Field
¤
Estimate an integral over a DomainComponent.
Given an integrand \(f\) and a component \(\Omega_{\text{comp}}\) with measure \(\mu\), this computes a Monte Carlo / quadrature estimate of
Sampling structure is provided by batch:
- PointsBatch: paired sampling according to a ProductStructure;
- CoordSeparableBatch: coordinate-separable sampling for some geometry labels.
Filtering and weighting:
- component.where / component.where_all act as indicator functions;
- component.weight_all multiplies the integrand.
Arguments:
f: Integrand as aDomainFunctionor array-like constant.batch: Sampled points (or a tuple of batches forDomainComponentUnion).component: Component (or union of components) to integrate over.quadrature: Optional explicit per-axis weights (QuadratureBatchforPointsBatchonly).over: Which label/block to integrate over;Noneintegrates over all axes implied bybatch.key: PRNG key forwarded towhere/weightcallables andf(when needed).kwargs: Extra keyword arguments forwarded tofand component callables.
Returns:
- A
coordax.Fieldcontaining the reduced integral value, with remaining named axes corresponding to any non-integrated sampling axes.
phydrax.operators.mean(f: DomainFunction | ArrayLike | None, batch: PointsBatch | CoordSeparableBatch | tuple[PointsBatch, ...], /, *, component: DomainComponent | DomainComponentUnion, quadrature: QuadratureBatch | tuple[QuadratureBatch | None, ...] | None = None, over: str | tuple[str, ...] | None = None, key: Key[Array, ''] = jr.key(0), **kwargs: Any) -> cx.Field
¤
Estimate the mean value of an integrand over a component.
Computes
Arguments:
f,batch,component,quadrature,over,key,kwargs: As inintegral.
Returns:
- A
coordax.Fieldcontaining the mean value.
phydrax.operators.integrate_interior(f: DomainFunction | ArrayLike | None, batch: PointsBatch | CoordSeparableBatch | tuple[PointsBatch, ...], /, *, component: DomainComponent | DomainComponentUnion, quadrature: QuadratureBatch | tuple[QuadratureBatch | None, ...] | None = None, over: str | tuple[str, ...] | None = None, key: Key[Array, ''] = jr.key(0), **kwargs: Any) -> cx.Field
¤
Alias for integral.
Interior/boundary semantics are encoded by component (e.g. its ComponentSpec).
phydrax.operators.integrate_boundary(f: DomainFunction | ArrayLike | None, batch: PointsBatch | CoordSeparableBatch | tuple[PointsBatch, ...], /, *, component: DomainComponent | DomainComponentUnion, quadrature: QuadratureBatch | tuple[QuadratureBatch | None, ...] | None = None, over: str | tuple[str, ...] | None = None, key: Key[Array, ''] = jr.key(0), **kwargs: Any) -> cx.Field
¤
Alias for integral.
Interior/boundary semantics are encoded by component (e.g. its ComponentSpec).
phydrax.operators.spatial_integral(u: DomainFunction, /, *, quad: dict, kernel: Callable | None = None, nonlinearity: Callable | None = None, importance_weight: Callable[[jax.Array, jax.Array], jax.Array] | None = None, var: str = 'x', time_var: str | None = None) -> DomainFunction
¤
Spatial integral operator with a fixed quadrature rule.
Constructs a new function \(v\) defined by
approximated using the provided quadrature nodes/weights quad = {"points": y_j, "weights": w_j}:
Arguments:
u: Input field \(u\).quad: A dict with keys"points"(shape(N_y, d)) and"weights"(shape(N_y,)).kernel: Optional kernel \(K\). If provided, it is evaluated askernel(concat([x, y]))with input in \(\mathbb{R}^{2d}\).nonlinearity: Optional nonlinearity \(g\) applied to \(u(y)\) before integration.importance_weight: Optional extra factor \(M(x,y)\) multiplied into the quadrature weights.var: Label of the spatial variable (default"x").time_var: Optional time label to include as an additional dependency.
Notes:
- This operator does not support coord-separable inputs for
var.
phydrax.operators.local_integral(u: DomainFunction, /, *, integrand: Callable, ball_quad: dict, var: str = 'x', time_var: str | None = None) -> DomainFunction
¤
Local (ball) integral operator around each query point.
Constructs a function of the form
approximated using a fixed quadrature rule on offsets \(\xi_j\):
The integrand is provided as integrand(ctx) and receives a context dictionary
similar to nonlocal_integral, including keys "x", "y", "ux", "uy", "du",
and "xi".
Notes:
- This operator does not support coord-separable inputs for
var.
phydrax.operators.local_integral_ball(u: DomainFunction, f_bond: Callable, *, ball_quad: dict, var: str = 'x', time_var: str | None = None) -> DomainFunction
¤
Convenience wrapper for a bond-based local integral.
Uses an integrand of the form \(\mathcal{I} = f(\Delta u, \xi)\) with \(\Delta u = u(y)-u(x)\) and \(\xi=y-x\).
phydrax.operators.nonlocal_integral(u: DomainFunction, /, *, integrand: Callable, quad: dict, importance_weight: Callable[[jax.Array, jax.Array], jax.Array] | None = None, var: str = 'x', time_var: str | None = None) -> DomainFunction
¤
Nonlocal integral operator with a context-aware integrand.
Constructs a function
where the integrand is provided as integrand(ctx) and is evaluated on a context
dictionary containing (at least) the keys:
"x": Full coordinate (including time iftime_varis provided)."y": Full coordinate (including time iftime_varis provided)."x_space": Spatial part ofx."y_space": Spatial part ofy."t": Scalar time orNone."ux": Value \(u(x)\)."uy": Value \(u(y)\)."du": \(u(y)-u(x)\)."xi": displacement \(y-x\) in space.
Arguments:
u: Input function \(u\).integrand: A callable returning the integrand value given actxdict.quad: A dict with"points"and"weights"as inspatial_integral.importance_weight: Optional extra factor \(M(x,y)\) multiplied into weights.var: Spatial domain label used for the integral variable.time_var: Optional time label included in"x"/"y"when present.
phydrax.operators.time_convolution(k: Callable[[Array], ArrayLike], u: DomainFunction, /, *, time_var: str = 't', order: int = 48, cluster_exponent: float = 1.0, mode: Literal['gl', 'qmc', 'gmc_1d'] = 'gl', sampler: str = 'sobol_scrambled', kernel_exponent: float | None = None) -> DomainFunction
¤
Time convolution operator on a labeled time coordinate.
Constructs the causal convolution
where \(t_0\) is the start of the time interval when time_var is a TimeInterval
factor (otherwise \(t_0=0\)).
The integral is approximated using one of:
- mode="gl": Gauss–Legendre quadrature on \([t_0,t]\) (optionally clustered near \(t_0\));
- mode="qmc": quasi Monte Carlo on \([t_0,t]\);
- mode="gmc_1d": importance sampling for weakly singular kernels
(requires kernel_exponent \(\gamma\)).
Arguments:
k: Kernel \(k(\tau)\) evaluated at \(\tau=t-s\).u: Input function \(u(s)\).time_var: Label for the time coordinate (default"t").order: Number of quadrature/MC samples.cluster_exponent: When using Gauss–Legendre, applies a power-law transform to cluster nodes toward the start of the interval.sampler: QMC sampler name used formode="qmc"/mode="gmc_1d".kernel_exponent: Exponent \(\gamma\) used inmode="gmc_1d"for sampling \(\tau^{-\gamma}\)-type singularities.