pub trait RungeKuttaSolver<N, const S: usize, const O: usize>: Sizedwhere
    N: ComplexField,{
    // Required methods
    fn t_coefficients() -> SVector<N::RealField, O>;
    fn k_coefficients() -> SMatrix<N::RealField, O, O>;
    fn avg_coefficients() -> SVector<N::RealField, O>;
    fn error_coefficients() -> SVector<N::RealField, O>;
    fn solve_ivp<T: Clone, F: FnMut(N::RealField, &[N], &mut T) -> Result<SVector<N, S>, String>>(
        self,
        f: F,
        params: &mut T
    ) -> Result<Vec<(N::RealField, SVector<N, S>)>, String>;
    fn with_tolerance(self, tol: N::RealField) -> Result<Self, String>;
    fn with_dt_max(self, max: N::RealField) -> Result<Self, String>;
    fn with_dt_min(self, min: N::RealField) -> Result<Self, String>;
    fn with_start(self, t_initial: N::RealField) -> Result<Self, String>;
    fn with_end(self, t_final: N::RealField) -> Result<Self, String>;
    fn with_initial_conditions(self, start: &[N]) -> Result<Self, String>;
    fn build(self) -> Self;
}
Expand description

This trait allows a struct to be used in the Runge-Kutta solver.

Things implementing RungeKuttaSolver should have an RKInfo to handle the actual IVP solving. It should also provide the with_* helper functions for convience.

Examples

See struct RK45 for an example of implementing this trait

Required Methods§

source

fn t_coefficients() -> SVector<N::RealField, O>

Returns a vec of coeffecients to multiply the time step by when getting intermediate results. Upper-left portion of Butch Tableaux

source

fn k_coefficients() -> SMatrix<N::RealField, O, O>

Returns the coefficients to use on the k_i’s when finding another k_i. Upper-right portion of the Butch Tableax. Should be an NxN-1 matrix, where N is the order of the Runge-Kutta Method (Or order+1 for adaptive methods)

source

fn avg_coefficients() -> SVector<N::RealField, O>

Coefficients to use when calculating the final step to take. These are the weights of the weighted average of k_i’s. Bottom portion of the Butch Tableaux. For adaptive methods, this is the first row of the bottom portion.

source

fn error_coefficients() -> SVector<N::RealField, O>

Coefficients to use on the k_i’s to find the error between the two orders of Runge-Kutta methods. In the Butch Tableaux, this is the first row of the bottom portion minus the second row.

source

fn solve_ivp<T: Clone, F: FnMut(N::RealField, &[N], &mut T) -> Result<SVector<N, S>, String>>( self, f: F, params: &mut T ) -> Result<Vec<(N::RealField, SVector<N, S>)>, String>

Ideally, call RKInfo.solve_ivp

source

fn with_tolerance(self, tol: N::RealField) -> Result<Self, String>

Set the error tolerance for this solver.

source

fn with_dt_max(self, max: N::RealField) -> Result<Self, String>

Set the maximum time step for this solver.

source

fn with_dt_min(self, min: N::RealField) -> Result<Self, String>

Set the minimum time step for this solver.

source

fn with_start(self, t_initial: N::RealField) -> Result<Self, String>

Set the initial time for this solver.

source

fn with_end(self, t_final: N::RealField) -> Result<Self, String>

Set the end time for this solver.

source

fn with_initial_conditions(self, start: &[N]) -> Result<Self, String>

Set the initial conditions for this solver.

source

fn build(self) -> Self

Build this solver.

Implementors§

source§

impl<N, const S: usize> RungeKuttaSolver<N, S, 4> for RK23<N, S>where N: ComplexField + FromPrimitive + Copy, <N as ComplexField>::RealField: FromPrimitive + Copy,

source§

impl<N, const S: usize> RungeKuttaSolver<N, S, 6> for RK45<N, S>where N: ComplexField + FromPrimitive + Copy, <N as ComplexField>::RealField: FromPrimitive + Copy,