pub trait ParametricOdeSystem<S: Scalar> {
Show 15 methods
// Required methods
fn n_states(&self) -> usize;
fn n_params(&self) -> usize;
fn params(&self) -> &[S];
fn rhs_with_params(&self, t: S, y: &[S], p: &[S], dydt: &mut [S]);
// Provided methods
fn rhs(&self, t: S, y: &[S], dydt: &mut [S]) { ... }
fn jacobian_y(&self, t: S, y: &[S], jac: &mut [S]) { ... }
fn jacobian_p(&self, t: S, y: &[S], jp: &mut [S]) { ... }
fn initial_sensitivity(&self, _y0: &[S], s0: &mut [S]) { ... }
fn has_analytical_jacobian_y(&self) -> bool { ... }
fn has_analytical_jacobian_p(&self) -> bool { ... }
fn is_autonomous(&self) -> bool { ... }
fn has_mass_matrix(&self) -> bool { ... }
fn mass_matrix(&self, mass: &mut [S]) { ... }
fn is_singular_mass(&self) -> bool { ... }
fn algebraic_indices(&self) -> Vec<usize> { ... }
}Expand description
An ODE system parameterised by a parameter vector p.
Implementors expose dimensions, the nominal parameter vector, the
parameterised right-hand side, and optionally analytical Jacobians.
Forward finite differences are provided as defaults for jacobian_y and
jacobian_p, so a minimal implementation supplies only n_states,
n_params, params, and rhs_with_params.
See the module-level documentation for the layout conventions
(J_y row-major, J_p and S column-major) and a worked 2-state,
2-parameter example.
Required Methods§
Sourcefn params(&self) -> &[S]
fn params(&self) -> &[S]
Nominal parameter vector. Length must equal Self::n_params.
Sourcefn rhs_with_params(&self, t: S, y: &[S], p: &[S], dydt: &mut [S])
fn rhs_with_params(&self, t: S, y: &[S], p: &[S], dydt: &mut [S])
Evaluate f(t, y, p) -> dy/dt for an arbitrary parameter vector p.
This is the primitive form: it lets the FD-default Jacobians perturb
p without mutating self.
Provided Methods§
Sourcefn rhs(&self, t: S, y: &[S], dydt: &mut [S])
fn rhs(&self, t: S, y: &[S], dydt: &mut [S])
Convenience: evaluate the RHS at the system’s nominal parameters.
Sourcefn jacobian_y(&self, t: S, y: &[S], jac: &mut [S])
fn jacobian_y(&self, t: S, y: &[S], jac: &mut [S])
Fill the state Jacobian J_y[i,j] = ∂f_i/∂y_j in row-major order
(jac[i*N + j], length N²).
Default: forward finite differences via Self::rhs_with_params.
Override to supply an analytical Jacobian — required for performance
on stiff problems with implicit solvers.
Sourcefn jacobian_p(&self, t: S, y: &[S], jp: &mut [S])
fn jacobian_p(&self, t: S, y: &[S], jp: &mut [S])
Fill the parameter Jacobian J_p[i,k] = ∂f_i/∂p_k in column-major
order (jp[k*N + i], length N · N_s).
Default: forward finite differences. Override for performance.
Sourcefn initial_sensitivity(&self, _y0: &[S], s0: &mut [S])
fn initial_sensitivity(&self, _y0: &[S], s0: &mut [S])
Initial sensitivity ∂y_0/∂p, column-major (s0[k*N + i], length
N · N_s).
Defaults to zero, which is correct when y_0 does not depend on p.
Override only if y_0 = y_0(p).
Sourcefn has_analytical_jacobian_y(&self) -> bool
fn has_analytical_jacobian_y(&self) -> bool
Returns true iff Self::jacobian_y has been overridden with an
analytical implementation. Default: false.
AugmentedSystem checks this flag inside its hot path. When false
(the default), it inlines its own forward-FD using reused scratch
buffers — bypassing both the trait default and any override.
When true, it calls system.jacobian_y(...) directly so the user’s
analytical override is honored.
This is the pattern used by CVODES (set_jacobian_user_supplied):
flagging analytical Jacobians lets the augmented-system path skip
FD scratch allocation in the common FD case while still respecting
analytical overrides in the stiff-problem case where they materially
matter. If you override jacobian_y, also override this method to
return true — otherwise your override is silently bypassed when
running through solve_forward_sensitivity.
§Debug-build safety net
To catch the silent-misconfiguration case (override the method,
forget the flag), AugmentedSystem runs a one-time consistency
check on the first RHS call: it evaluates Self::jacobian_y
(which dispatches to the user’s override or the FD default) and
compares against the inline FD result. If they disagree by more
than a Frobenius-norm relative threshold of 1e-3, the system
panics in debug builds with a message naming the missing flag
override. Release builds skip the check entirely (zero overhead).
This is a safety net, not a contract — fixing the panic by
returning true here is the canonical resolution.
Sourcefn has_analytical_jacobian_p(&self) -> bool
fn has_analytical_jacobian_p(&self) -> bool
Returns true iff Self::jacobian_p has been overridden with an
analytical implementation. Default: false. See
Self::has_analytical_jacobian_y for the rationale, contract,
and debug-build consistency check.
Sourcefn is_autonomous(&self) -> bool
fn is_autonomous(&self) -> bool
Is the system autonomous? (f does not depend on t explicitly)
Default: false. Override to true when f(t, y, p) = f(y, p)
independent of t. AugmentedSystem forwards this from the host
so the augmented system inherits the host’s time-dependence
structure — variational equations inherit autonomy from the
state equation.
Sourcefn has_mass_matrix(&self) -> bool
fn has_mass_matrix(&self) -> bool
Does this system have a mass matrix?
Default: false (the standard parameterised ODE dy/dt = f(t, y, p)).
Override to true for DAE-shaped parameterised systems
M · y' = f(t, y, p) where M is non-identity or singular.
See OdeSystem::has_mass_matrix for the equivalent surface on
the unparameterised side.
Sourcefn mass_matrix(&self, mass: &mut [S])
fn mass_matrix(&self, mass: &mut [S])
Get the mass matrix M for the DAE: M · y' = f(t, y, p).
Default returns identity (standard ODE), regardless of what
Self::has_mass_matrix reports — callers should consult
Self::has_mass_matrix before relying on this output.
Layout: row-major, mass[i * n + j] = M[i, j], length n² where
n = self.n_states(). Mirrors OdeSystem::mass_matrix.
Sourcefn is_singular_mass(&self) -> bool
fn is_singular_mass(&self) -> bool
Is the mass matrix singular? (i.e., is this a DAE?)
Default: false. Override to true for semi-explicit DAEs where
one or more rows of M are zero (algebraic constraints).
Sourcefn algebraic_indices(&self) -> Vec<usize>
fn algebraic_indices(&self) -> Vec<usize>
Indices of algebraic variables (rows i where M[i, i] = 0).
Default: empty (all rows differential). Override for DAE systems
to enumerate the algebraic rows so downstream solvers can apply
row-aware error scaling and the augmented-system lift can map
each host algebraic index i to its (N_s + 1) copies on the
augmented diagonal.
Implementations on Foreign Types§
Source§impl<S: Scalar, T: ParametricOdeSystem<S>> ParametricOdeSystem<S> for &T
Forwarding impl: a reference to a parametric system is itself a parametric
system. Lets solve_forward_sensitivity accept &Sys without taking
ownership of the user’s data.
impl<S: Scalar, T: ParametricOdeSystem<S>> ParametricOdeSystem<S> for &T
Forwarding impl: a reference to a parametric system is itself a parametric
system. Lets solve_forward_sensitivity accept &Sys without taking
ownership of the user’s data.