pub trait ODE<T = f64, V = f64, D = String>{
// Required method
fn diff(&self, t: T, y: &V, dydt: &mut V);
// Provided methods
fn event(&self, t: T, y: &V) -> ControlFlag<T, V, D> { ... }
fn jacobian(&self, t: T, y: &V, j: &mut DMatrix<T>) { ... }
}Expand description
ODE Trait for Differential Equations
ODE trait defines the differential equation dydt = f(t, y) for the solver. The differential equation is used to solve the ordinary differential equation. The trait also includes a solout function to interupt the solver when a condition is met or event occurs.
§Impl
diff- Differential Equation dydt = f(t, y) in form f(t, &y, &mut dydt).event- Solout function to interupt solver when condition is met or event occurs.
Note that the solout function is optional and can be left out when implementing in which can it will be set to return false by default.
Required Methods§
Sourcefn diff(&self, t: T, y: &V, dydt: &mut V)
fn diff(&self, t: T, y: &V, dydt: &mut V)
Differential Equation dydt = f(t, y)
An ordinary differential equation (ODE) takes a independent variable which in this case is ‘t’ as it is typically time and a dependent variable which is a vector of values ‘y’. The ODE returns the derivative of the dependent variable ‘y’ with respect to the independent variable ‘t’ as dydt = f(t, y).
For efficiency and ergonomics the derivative is calculated from an argument of a mutable reference to the derivative vector dydt. This allows for a derivatives to be calculated in place which is more efficient as iterative ODE solvers require the derivative to be calculated at each step without regard to the previous value.
§Arguments
t- Independent variable point.y- Dependent variable point.dydt- Derivative point.
Provided Methods§
Sourcefn event(&self, t: T, y: &V) -> ControlFlag<T, V, D>
fn event(&self, t: T, y: &V) -> ControlFlag<T, V, D>
Event function to detect significant conditions during integration.
Called after each step to detect events like threshold crossings, singularities, or other mathematically/physically significant conditions. Can be used to terminate integration when conditions occur.
§Arguments
t- Current independent variable point.y- Current dependent variable point.
§Returns
ControlFlag- Command to continue or stop solver.
Sourcefn jacobian(&self, t: T, y: &V, j: &mut DMatrix<T>)
fn jacobian(&self, t: T, y: &V, j: &mut DMatrix<T>)
Jacobian matrix J = df/dy
The Jacobian matrix is a matrix of partial derivatives of a vector-valued function. It describes the local behavior of the system of equations and can be used to improve the efficiency of certain solvers by providing information about the local behavior of the system of equations.
By default, this method uses a finite difference approximation. Users can override this with an analytical implementation for better efficiency.
§Arguments
t- Independent variable grid point.y- Dependent variable vector.j- Jacobian matrix. This matrix should be pre-sized by the caller todim x dimwheredim = y.len().