Continuous constraints¤
These helpers construct sampled constraints over domain components. For details on reduction, measures, and filtering, see Guides → Constraints and objectives.
For residual-style continuous constraints, weight can be either a scalar global
multiplier or a pointwise DomainFunction.
Interior / initial sampling constraints¤
phydrax.constraints.ContinuousPointwiseInteriorConstraint(constraint_vars: str | Sequence[str], domain, operator, *, num_points: SamplingNumPoints, structure: ProductStructure, dense_structure: ProductStructure | None = None, sampler: str = 'latin_hypercube', weight: DomainFunction | ArrayLike = 1.0, label: str | None = None, over: str | tuple[str, ...] | None = None, reduction: Literal[mean, integral] = 'mean', sampling_mode: Literal[resample, fixed] = 'resample', fixed_batch: PointsBatch | CoordSeparableBatch | tuple[PointsBatch, ...] | None = None, fixed_batch_key: Key[Array, ''] = jr.key(0), where: Mapping[str, Any] | None = None, where_all: DomainFunction | None = None) -> FunctionalConstraint
¤
Pointwise residual constraint over an interior domain component.
Let operator define a residual DomainFunction \(r\) from one or more fields:
where \(r(z)\) is evaluated at sampled interior points \(z\in\Omega\).
The resulting objective is the mean or integral of the squared residual.
For reduction="mean":
For reduction="integral":
Arguments:
constraint_vars: Name (or names) of the field functions passed intooperator.domain: Aphydrax.domainobject to sample from.operator: Callable mapping one or moreDomainFunctionobjects to a residualDomainFunction.num_points: Sampling spec. Accepts dense counts (int/tuple[int, ...]), coord-separable mappings (e.g.{"x": 64}), or mixed forms(dense_num_points, coord_map).structure: AProductStructuredescribing how variables are sampled/blocked.dense_structure: Optional dense structure used when sampling produces dense batches.sampler: Sampling scheme (e.g."latin_hypercube").weight: Scalar multiplier \(w\) applied to this term.over: Optional subset of labels to reduce/integrate over.reduction: Reduction mode:"mean"(measure-normalized) or"integral"(unnormalized).sampling_mode:"resample"(default; sample every call) or"fixed"(reuse one batch).fixed_batch: Optional prebuilt batch used whensampling_mode="fixed".fixed_batch_key: Key used to sample the fixed batch whensampling_mode="fixed"andfixed_batchis not provided.where: Optional per-label filters, treated as indicator functions.where_all: Optional global filter, evaluated on the full point tuple.
phydrax.constraints.ContinuousInitialConstraint(constraint_var: str, component: DomainComponent | DomainComponentUnion, /, *, evolution_var: str = 't', func: DomainFunction | ArrayLike | None = None, time_derivative_order: int = 0, mode: Literal[reverse, forward] = 'reverse', time_derivative_backend: Literal[ad, jet] = 'ad', num_points: SamplingNumPoints, structure: ProductStructure, dense_structure: ProductStructure | None = None, sampler: str = 'latin_hypercube', weight: DomainFunction | ArrayLike = 1.0, label: str | None = None, over: str | tuple[str, ...] | None = None, reduction: Literal[mean, integral] = 'mean') -> FunctionalConstraint
¤
Continuous initial-condition constraint on a fixed-start time slice.
Enforces
by sampling the initial slice, where n = time_derivative_order.
Arguments:
constraint_var: Name of the constrained field.component: Domain component on the initial surface (or a union of such components).evolution_var: Time label used for the initial slice (default:"t").func: Target value/function \(g\) (defaults to 0).time_derivative_order: Derivative order \(n\) for \(\partial^n/\partial t^n\).mode: Differentiation mode ("reverse"or"forward").time_derivative_backend: Backend for time derivatives ("ad"or"jet").num_points: Sampling spec. Accepts dense counts (int/tuple[int, ...]), coord-separable mappings (e.g.{"x": 64}), or mixed forms(dense_num_points, coord_map).structure: AProductStructuredescribing how variables are sampled/blocked.dense_structure: Optional dense structure used when sampling produces dense batches.sampler: Sampling scheme (e.g."latin_hypercube").weight: Scalar multiplier applied to this term.label: Optional label for logging.over: Optional subset of labels to reduce/integrate over.reduction:"mean"or"integral".
phydrax.constraints.ContinuousInitialFunctionConstraint(constraint_vars: str | Sequence[str], domain, /, *, func: DomainFunction | ArrayLike | None = None, evolution_var: str = 't', time_derivative_order: int = 0, mode: Literal[reverse, forward] = 'reverse', time_derivative_backend: Literal[ad, jet] = 'ad', num_points: SamplingNumPoints, structure: ProductStructure, dense_structure: ProductStructure | None = None, sampler: str = 'latin_hypercube', weight: DomainFunction | ArrayLike = 1.0, label: str | None = None, over: str | tuple[str, ...] | None = None, reduction: Literal[mean, integral] = 'mean', where: Mapping[str, Any] | None = None, where_all: DomainFunction | None = None) -> FunctionalConstraint
¤
Initial-condition constraint matching time derivatives on the initial slice.
Builds a constraint on the component slice evolution_var = FixedStart() and
enforces
where n = time_derivative_order, \(t_0\) is the start of the time interval, and
g is provided by func (defaults to \(0\)).
The loss is formed by sampling points on the initial slice and reducing the
squared residual as in FunctionalConstraint.
Arguments:
constraint_vars: Name (or names) of the field functions passed into the residual.domain: Aphydrax.domainobject to sample from.func: Target value \(g\) as aDomainFunction, callable, or array-like. Defaults to0.0.evolution_var: Name of the time-like label used for the initial slice (default"t").time_derivative_order: Derivative order \(n\) for \(\partial^n/\partial t^n\).mode: Differentiation mode ("reverse"or"forward").time_derivative_backend: Backend for time derivatives ("ad"or"jet").num_points: Sampling spec. Accepts dense counts (int/tuple[int, ...]), coord-separable mappings (e.g.{"x": 64}), or mixed forms(dense_num_points, coord_map).structure: AProductStructuredescribing how variables are sampled/blocked.dense_structure: Optional dense structure used when sampling produces dense batches.sampler: Sampling scheme (e.g."latin_hypercube").weight: Scalar multiplier applied to this term.over: Optional subset of labels to reduce/integrate over.reduction:"mean"or"integral".where: Optional per-label filters, treated as indicator functions.where_all: Optional global filter, evaluated on the full point tuple.
Integral constraints¤
phydrax.constraints.ContinuousIntegralInteriorConstraint(constraint_vars: str | Sequence[str], domain, operator: Callable[..., DomainFunction] | DomainFunction, /, *, num_points: NumPoints | tuple[Any, ...], structure: ProductStructure | None = None, sampler: str = 'latin_hypercube', weight: ArrayLike = 1.0, label: str | None = None, over: str | tuple[str, ...] | None = None, equal_to: ArrayLike | None = None, where: Mapping[str, Callable] | None = None, where_all: DomainFunction | None = None) -> IntegralEqualityConstraint
¤
Integral equality constraint over an interior component.
Enforces an interior integral condition of the form
where the integrand \(f\) is produced by operator applied to the named
constraint_vars.
phydrax.constraints.ContinuousIntegralBoundaryConstraint(constraint_vars: str | Sequence[str], domain, operator: Callable[..., DomainFunction] | DomainFunction, /, *, num_points: NumPoints | tuple[Any, ...], structure: ProductStructure | None = None, sampler: str = 'latin_hypercube', weight: ArrayLike = 1.0, label: str | None = None, over: str | tuple[str, ...] | None = None, equal_to: ArrayLike | None = None, where: Mapping[str, Callable] | None = None, where_all: DomainFunction | None = None) -> IntegralEqualityConstraint
¤
Integral equality constraint over a boundary component.
Enforces a boundary integral condition
where operator receives the boundary normal field \(n(x)\) as an additional
argument.
phydrax.constraints.ContinuousIntegralInitialConstraint(constraint_vars: str | Sequence[str], domain, operator: Callable[..., DomainFunction] | DomainFunction, /, *, evolution_var: str = 't', num_points: NumPoints | tuple[Any, ...], structure: ProductStructure | None = None, sampler: str = 'latin_hypercube', weight: ArrayLike = 1.0, label: str | None = None, over: str | tuple[str, ...] | None = None, equal_to: ArrayLike | None = None, where: Mapping[str, Callable] | None = None, where_all: DomainFunction | None = None) -> IntegralEqualityConstraint
¤
Integral equality constraint over the initial surface.
Constructs a component slice evolution_var = FixedStart() (e.g. \(t=t_0\)) and
enforces an integral condition on that slice.
ODE constraints¤
phydrax.constraints.ContinuousODEConstraint(constraint_var: str, domain: TimeInterval, operator: Callable[[DomainFunction], DomainFunction], /, *, num_points: NumPoints | tuple[Any, ...], structure: ProductStructure | None = None, sampler: str = 'latin_hypercube', weight: DomainFunction | ArrayLike = 1.0, label: str | None = None, over: str | tuple[str, ...] | None = None, reduction: Literal[mean, integral] = 'mean') -> FunctionalConstraint
¤
Continuous ODE residual constraint on a TimeInterval.
Given an ODE residual operator \(\mathcal{R}\) (provided by operator), this
enforces \(\mathcal{R}(u)(t)=0\) in a sampled/weak sense by minimizing
(or the unnormalized integral when reduction="integral").
phydrax.constraints.InitialODEConstraint(constraint_var: str, domain: TimeInterval, /, *, func: DomainFunction | ArrayLike | Callable | None = None, time_derivative_order: int = 0, time_derivative_backend: Literal[ad, jet] = 'ad', weight: DomainFunction | ArrayLike = 1.0, label: str | None = None, structure: ProductStructure | None = None, sampler: str = 'latin_hypercube', reduction: Literal[mean, integral] = 'mean') -> FunctionalConstraint
¤
Initial-condition constraint on the slice \(t=t_0\).
Enforces
where n = time_derivative_order and g is given by func (defaults to \(0\)).