Skip to content

Enforced constraint ansätze¤

These helpers construct ansätze that satisfy constraints by construction.

For composition/ordering (multi-field dependencies, applying several enforced constraints, etc.), see Enforced constraint pipelines.

Warning

Compatibility with coord-separable grids:

  • enforce_neumann, enforce_robin, enforce_traction, and enforce_sommerfeld rely on geometry boundary normals \(n(x)\) (via \(\partial/\partial n\)) and therefore do not support coord-separable (tuple-of-axes) evaluation. Phydrax raises a ValueError if you try to evaluate these ansätze on a CoordSeparableBatch.
  • enforce_dirichlet, enforce_initial, and enforce_blend do not require boundary normals and can be used in spectral/FNO-style interior grid evaluations.

phydrax.constraints.enforce_dirichlet(u: DomainFunction, component: DomainComponent, /, *, var: str = 'x', target: DomainFunction | ArrayLike | None = None) -> DomainFunction ¤

Enforced Dirichlet ansatz enforcing \(u=g\) on a component.

Constructs an ansatz \(u^*\) that satisfies the Dirichlet condition exactly on the selected component (boundary or fixed slice). For a geometry boundary with signed distance field \(\phi\) (so \(\phi=0\) on \(\partial\Omega\)), this uses

\[ u^*(x) = g(x) + \phi(x)\,(u(x) - g(x)), \]

which guarantees \(u^*(x)=g(x)\) on \(\partial\Omega\).

For scalar domains (e.g. time), an appropriate vanishing factor \(\phi(t)\) is constructed from \((t-t_0)\), \((t-t_1)\), etc.


phydrax.constraints.enforce_neumann(u: DomainFunction, component: DomainComponent, /, *, var: str = 'x', target: DomainFunction | ArrayLike | None = None, mode: Literal['reverse', 'forward'] = 'reverse') -> DomainFunction ¤

Enforced Neumann ansatz enforcing \(\partial u/\partial n = g\) on a boundary.

For a geometry boundary with signed distance field \(\phi\) and outward normal \(n\), this constructs

\[ u^* = u + \frac{\phi}{\partial\phi/\partial n}\,\bigl(g - \partial u/\partial n\bigr), \]

which yields \(\partial u^*/\partial n = g\) on \(\partial\Omega\) under mild regularity assumptions.


phydrax.constraints.enforce_robin(u: DomainFunction, component: DomainComponent, /, *, var: str = 'x', dirichlet_coeff: DomainFunction | ArrayLike | None = None, neumann_coeff: DomainFunction | ArrayLike | None = None, target: DomainFunction | ArrayLike | None = None, mode: Literal['reverse', 'forward'] = 'reverse') -> DomainFunction ¤

Enforced Robin ansatz enforcing \(a\,u + b\,\partial u/\partial n = g\).

On a geometry boundary, this constructs a corrected field that satisfies the Robin condition exactly by using a signed distance factor \(\phi\) and normal derivative of \(\phi\).

Let \(r = g - a\,u - b\,\partial u/\partial n\) be the residual. With signed distance field \(\phi\) and outward normal \(n\), define \(d\phi/dn = \partial\phi/\partial n\). The returned ansatz is

\[ u^* \;=\; u + \frac{\phi}{b\,(\partial\phi/\partial n)}\,r, \]

which yields \(a\,u^* + b\,\partial u^*/\partial n = g\) on \(\partial\Omega\) under mild regularity assumptions.


phydrax.constraints.enforce_sommerfeld(u: DomainFunction, component: DomainComponent, /, *, var: str = 'x', time_var: str = 't', wavespeed: DomainFunction | ArrayLike | None = None, target: DomainFunction | ArrayLike | None = None, mode: Literal['reverse', 'forward'] = 'reverse') -> DomainFunction ¤

Enforced Sommerfeld/absorbing boundary ansatz.

Enforces the first-order radiation condition

\[ \frac{\partial u}{\partial n} + \frac{1}{c}\frac{\partial u}{\partial t} = g \]

on the selected boundary component by correcting \(u\) using a signed-distance factor.

With signed distance field \(\phi\) and \(d\phi/dn = \partial\phi/\partial n\), let

\[ r = g - \frac{\partial u}{\partial n} - \frac{1}{c}\frac{\partial u}{\partial t}. \]

The returned ansatz is

\[ u^* \;=\; u + \frac{\phi}{\partial\phi/\partial n}\,r, \]

which yields the Sommerfeld condition on \(\partial\Omega\) under the same assumptions as enforce_neumann.


phydrax.constraints.enforce_traction(u: DomainFunction, component: DomainComponent, /, *, var: str = 'x', target: DomainFunction | ArrayLike | None = None, lambda_: DomainFunction | ArrayLike, mu: DomainFunction | ArrayLike, mode: Literal['reverse', 'forward'] = 'reverse') -> DomainFunction ¤

Enforced traction ansatz for linear elasticity on a boundary.

For isotropic linear elasticity, the traction is \(t=\sigma(u)\,n\). This constructs a corrected displacement \(u^*\) that aims to satisfy

\[ \sigma(u^*)\,n = t_{\text{target}} \]

on the boundary by adding a displacement correction proportional to the signed distance \(\phi\).

Let

\[ r \;=\; t_{\text{target}} - \sigma(u)\,n \]

be the traction residual, and decompose it into normal/tangential parts

\[ r_n = (r\cdot n)\,n,\qquad r_t = r - r_n. \]

Using Lamé parameters \(\lambda,\mu\), define the correction field

\[ v \;=\; \frac{r_t}{\mu} \;+\; \frac{r_n}{\lambda + 2\mu}, \]

and return the hard-constraint ansatz

\[ u^*(x) \;=\; u(x) + \phi(x)\,v(x). \]

In the idealized setting where \(\phi\) is a signed distance field near the boundary (so \(\partial\phi/\partial n \approx 1\)) and \(v\) varies slowly in the normal direction, the induced traction correction is approximately \(\sigma(\phi v)\,n \approx \mu\,v_t + (\lambda+2\mu)\,v_n\) (with \(v_n=(v\cdot n)\,n\) and \(v_t=v-v_n\)), making \(u^*\) enforce the target traction to first order.


phydrax.constraints.enforce_initial(u: DomainFunction, component: DomainComponent, /, *, var: str = 't', targets: Mapping[int, DomainFunction | ArrayLike], gate: Literal['rational', 'poly'] = 'rational', gate_eps: float = 0.01) -> DomainFunction ¤

Enforced initial-condition ansatz matching time derivatives at a fixed time.

Let \(t_0\) be the fixed time selected by component (e.g. FixedStart()). Given targets for derivatives \(u^{(k)}(t_0)\) for \(k=0,\dots,m\), this constructs a Taylor-like polynomial

\[ p(t) = \sum_{k=0}^{m} \frac{u^{(k)}(t_0)}{k!}\,(t-t_0)^k \]

and returns

\[ u^*(t) = p(t) + g(t)\,(u(t)-p(t)), \]

which matches all specified derivatives exactly at \(t=t_0\).

By default, g is a bounded rational gate based on normalized time-from-\(t_0\) that avoids unbounded polynomial growth away from the initial slice. Set gate="poly" to recover the original \((t-t_0)^{m+1}\) gate.


phydrax.constraints.enforce_blend(u: DomainFunction, pieces: Sequence[tuple[DomainComponent, DomainFunction]], /, *, var: str = 'x', include_identity_remainder: bool = True, num_reference: int = 3000000, sampler: str = 'latin_hypercube', key: Key[Array, ''] = jr.key(0)) -> DomainFunction ¤

Blend multiple enforced ansatz pieces via ported MLS/BVH weights.

Each piece contributes \(w_i\,u_i^*\) and the final output is:

\[ u_{\text{enforced}} = \frac{\sum_i w_i\,u_i^*}{\sum_i w_i}. \]

If include_identity_remainder=True, this also adds an identity term for the complement of the union of piece predicates (using component.where[var]), which prevents subset enforced constraints from leaking onto other boundary segments.