Struct ExplicitRungeKutta

Source
pub struct ExplicitRungeKutta<E, F, T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> {
    pub h0: T,
    pub rtol: T,
    pub atol: T,
    pub h_max: T,
    pub h_min: T,
    pub max_steps: usize,
    pub max_rejects: usize,
    pub safety_factor: T,
    pub min_scale: T,
    pub max_scale: T,
    /* private fields */
}
Expand description

Runge-Kutta solver that can handle:

  • Fixed-step methods with cubic Hermite interpolation
  • Adaptive step methods with embedded error estimation and cubic Hermite interpolation
  • Advanced methods with dense output interpolation using Butcher tableau coefficients

§Type Parameters

  • E: Equation type (e.g., Ordinary, Delay, Stochastic)
  • F: Family type (e.g., Adaptive, Fixed, DormandPrince, etc.)
  • T: Real number type (f32, f64)
  • Y: State vector type
  • D: Callback data type
  • const O: Order of the method
  • const S: Number of stages in the method
  • const I: Total number of stages including interpolation (equal to S for methods without dense output)

Fields§

§h0: T§rtol: T§atol: T§h_max: T§h_min: T§max_steps: usize§max_rejects: usize§safety_factor: T§min_scale: T§max_scale: T

Implementations§

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 5, 6, 6>

Source

pub fn rkf45() -> Self

Creates a Runge-Kutta-Fehlberg 4(5) method with error estimation.

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 5, 6, 6>

Source

pub fn cash_karp() -> Self

Creates a Cash-Karp 4(5) method with error estimation.

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 6, 9, 10>

Source

pub fn rkv655e() -> Self

Creates an Explicit Runge-Kutta 6(5) method with 9 stages and a 5th order interpolant.

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 6, 9, 12>

Source

pub fn rkv656e() -> Self

Creates an Explicit Runge-Kutta 6(5) method with 9 stages and a 6th order interpolant.

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 7, 10, 13>

Source

pub fn rkv766e() -> Self

Creates an Explicit Runge-Kutta 7(6) method with 10 stages and a 6th order interpolant.

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 7, 10, 16>

Source

pub fn rkv767e() -> Self

Creates an Explicit Runge-Kutta 7(6) method with 10 stages and a 7th order interpolant.

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 8, 13, 17>

Source

pub fn rkv877e() -> Self

Creates an Explicit Runge-Kutta 8(7) method with 13 stages and a 7th order interpolant.

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 8, 13, 21>

Source

pub fn rkv878e() -> Self

Creates an Explicit Runge-Kutta 8(7) method with 13 stages and an 8th order interpolant.

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 9, 16, 21>

Source

pub fn rkv988e() -> Self

Creates an Explicit Runge-Kutta 9(8) method with 16 stages and an 8th order interpolant.

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Adaptive, T, Y, D, 9, 16, 26>

Source

pub fn rkv989e() -> Self

Creates an Explicit Runge-Kutta 9(8) method with 16 stages and a 9th order interpolant.

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, DormandPrince, T, Y, D, 8, 12, 16>

Source

pub fn dop853() -> Self

Creates the DOP853 method (8th order, 12 stages, 4 dense output stages).

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, DormandPrince, T, Y, D, 5, 7, 7>

Source

pub fn dopri5() -> Self

Creates the DOPRI5 method (5th order, 7 stages).

Source§

impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> ExplicitRungeKutta<Delay, Fixed, T, Y, D, O, S, I>

Source

pub fn lagvals<const L: usize, H>( &mut self, t_stage: T, delays: &[T; L], y_delayed: &mut [Y; L], phi: &H, ) -> Result<(), Error<T, Y>>
where H: Fn(T) -> Y,

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 1, 1, 1>

Source

pub fn euler(h0: T) -> Self

Creates an Explicit Euler method (1st order, 1 stage).

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 2, 2, 2>

Source

pub fn midpoint(h0: T) -> Self

Creates an Explicit Midpoint method (2nd order, 2 stages).

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 2, 2, 2>

Source

pub fn heun(h0: T) -> Self

Creates an Explicit Heun method (2nd order, 2 stages).

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 2, 2, 2>

Source

pub fn ralston(h0: T) -> Self

Creates an Explicit Ralston method (2nd order, 2 stages).

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 4, 4, 4>

Source

pub fn rk4(h0: T) -> Self

Creates the classical 4th order Runge-Kutta method.

Source§

impl<E, T: Real, Y: State<T>, D: CallBackData> ExplicitRungeKutta<E, Fixed, T, Y, D, 4, 4, 4>

Source

pub fn three_eighths(h0: T) -> Self

Creates the three-eighths rule 4th order Runge-Kutta method.

Source§

impl<E, F, T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> ExplicitRungeKutta<E, F, T, Y, D, O, S, I>

Source

pub fn rtol(self, rtol: T) -> Self

Set the relative tolerance for error control

Source

pub fn atol(self, atol: T) -> Self

Set the absolute tolerance for error control

Source

pub fn h0(self, h0: T) -> Self

Set the initial step size

Source

pub fn h_min(self, h_min: T) -> Self

Set the minimum allowed step size

Source

pub fn h_max(self, h_max: T) -> Self

Set the maximum allowed step size

Source

pub fn max_steps(self, max_steps: usize) -> Self

Set the maximum number of steps allowed

Source

pub fn max_rejects(self, max_rejects: usize) -> Self

Set the maximum number of consecutive rejected steps before declaring stiffness

Source

pub fn safety_factor(self, safety_factor: T) -> Self

Set the safety factor for step size control (default: 0.9)

Source

pub fn min_scale(self, min_scale: T) -> Self

Set the minimum scale factor for step size changes (default: 0.2)

Source

pub fn max_scale(self, max_scale: T) -> Self

Set the maximum scale factor for step size changes (default: 10.0)

Source

pub fn order(&self) -> usize

Get the order of the method

Source

pub fn stages(&self) -> usize

Get the number of stages in the method

Source

pub fn dense_stages(&self) -> usize

Get the number of terms in the dense output interpolation polynomial

Source§

impl<F, T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> ExplicitRungeKutta<Delay, F, T, Y, D, O, S, I>

Source

pub fn max_delay(self, max_delay: T) -> Self

Set the maximum delay for DDEs

Trait Implementations§

Source§

impl<E, F, T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Default for ExplicitRungeKutta<E, F, T, Y, D, O, S, I>

Source§

fn default() -> Self

Returns the “default value” for a type. Read more
Source§

impl<const L: usize, T: Real, Y: State<T>, H: Fn(T) -> Y, D: CallBackData, const O: usize, const S: usize, const I: usize> DelayNumericalMethod<L, T, Y, H, D> for ExplicitRungeKutta<Delay, Adaptive, T, Y, D, O, S, I>

Source§

fn init<F>( &mut self, dde: &F, t0: T, tf: T, y0: &Y, phi: &H, ) -> Result<Evals, Error<T, Y>>
where F: DDE<L, T, Y, D>,

Initialize the solver before integration Read more
Source§

fn step<F>(&mut self, dde: &F, phi: &H) -> Result<Evals, Error<T, Y>>
where F: DDE<L, T, Y, D>,

Advance the solution by one step Read more
Source§

fn t(&self) -> T

Time of last accepted step
Source§

fn y(&self) -> &Y

State at last accepted step
Source§

fn t_prev(&self) -> T

Time of previous accepted step
Source§

fn y_prev(&self) -> &Y

State at previous accepted step
Source§

fn h(&self) -> T

Step size for next step
Source§

fn set_h(&mut self, h: T)

Set step size for next step
Source§

fn status(&self) -> &Status<T, Y, D>

Current solver status
Source§

fn set_status(&mut self, status: Status<T, Y, D>)

Set solver status
Source§

impl<const L: usize, T: Real, Y: State<T>, H: Fn(T) -> Y, D: CallBackData, const O: usize, const S: usize, const I: usize> DelayNumericalMethod<L, T, Y, H, D> for ExplicitRungeKutta<Delay, DormandPrince, T, Y, D, O, S, I>

Source§

fn init<F>( &mut self, dde: &F, t0: T, tf: T, y0: &Y, phi: &H, ) -> Result<Evals, Error<T, Y>>
where F: DDE<L, T, Y, D>,

Initialize the solver before integration Read more
Source§

fn step<F>(&mut self, dde: &F, phi: &H) -> Result<Evals, Error<T, Y>>
where F: DDE<L, T, Y, D>,

Advance the solution by one step Read more
Source§

fn t(&self) -> T

Time of last accepted step
Source§

fn y(&self) -> &Y

State at last accepted step
Source§

fn t_prev(&self) -> T

Time of previous accepted step
Source§

fn y_prev(&self) -> &Y

State at previous accepted step
Source§

fn h(&self) -> T

Step size for next step
Source§

fn set_h(&mut self, h: T)

Set step size for next step
Source§

fn status(&self) -> &Status<T, Y, D>

Current solver status
Source§

fn set_status(&mut self, status: Status<T, Y, D>)

Set solver status
Source§

impl<const L: usize, T: Real, Y: State<T>, H: Fn(T) -> Y, D: CallBackData, const O: usize, const S: usize, const I: usize> DelayNumericalMethod<L, T, Y, H, D> for ExplicitRungeKutta<Delay, Fixed, T, Y, D, O, S, I>

Source§

fn init<F>( &mut self, dde: &F, t0: T, tf: T, y0: &Y, phi: &H, ) -> Result<Evals, Error<T, Y>>
where F: DDE<L, T, Y, D>,

Initialize the solver before integration Read more
Source§

fn step<F>(&mut self, dde: &F, phi: &H) -> Result<Evals, Error<T, Y>>
where F: DDE<L, T, Y, D>,

Advance the solution by one step Read more
Source§

fn t(&self) -> T

Time of last accepted step
Source§

fn y(&self) -> &Y

State at last accepted step
Source§

fn t_prev(&self) -> T

Time of previous accepted step
Source§

fn y_prev(&self) -> &Y

State at previous accepted step
Source§

fn h(&self) -> T

Step size for next step
Source§

fn set_h(&mut self, h: T)

Set step size for next step
Source§

fn status(&self) -> &Status<T, Y, D>

Current solver status
Source§

fn set_status(&mut self, status: Status<T, Y, D>)

Set solver status
Source§

impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Delay, Adaptive, T, Y, D, O, S, I>

Source§

fn interpolate(&mut self, t_interp: T) -> Result<Y, Error<T, Y>>

Interpolates the solution at a given time t_interp.

Source§

impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Delay, DormandPrince, T, Y, D, O, S, I>

Source§

fn interpolate(&mut self, t_interp: T) -> Result<Y, Error<T, Y>>

Evaluate the step-local interpolant at the given time. Read more
Source§

impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Delay, Fixed, T, Y, D, O, S, I>

Source§

fn interpolate(&mut self, t_interp: T) -> Result<Y, Error<T, Y>>

Interpolates the solution at time t_interp within the last accepted step.

Source§

impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Ordinary, Adaptive, T, Y, D, O, S, I>

Source§

fn interpolate(&mut self, t_interp: T) -> Result<Y, Error<T, Y>>

Evaluate the step-local interpolant at the given time. Read more
Source§

impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Ordinary, DormandPrince, T, Y, D, O, S, I>

Source§

fn interpolate(&mut self, t_interp: T) -> Result<Y, Error<T, Y>>

Evaluate the step-local interpolant at the given time. Read more
Source§

impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Ordinary, Fixed, T, Y, D, O, S, I>

Source§

fn interpolate(&mut self, t_interp: T) -> Result<Y, Error<T, Y>>

Evaluate the step-local interpolant at the given time. Read more
Source§

impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ExplicitRungeKutta<Stochastic, Fixed, T, Y, D, O, S, I>

Source§

fn interpolate(&mut self, t_interp: T) -> Result<Y, Error<T, Y>>

Evaluate the step-local interpolant at the given time. Read more
Source§

impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> OrdinaryNumericalMethod<T, Y, D> for ExplicitRungeKutta<Ordinary, Adaptive, T, Y, D, O, S, I>

Source§

fn init<F>( &mut self, ode: &F, t0: T, tf: T, y0: &Y, ) -> Result<Evals, Error<T, Y>>
where F: ODE<T, Y, D>,

Initialize the solver before integration Read more
Source§

fn step<F>(&mut self, ode: &F) -> Result<Evals, Error<T, Y>>
where F: ODE<T, Y, D>,

Advance the solution by one step Read more
Source§

fn t(&self) -> T

Time of last accepted step
Source§

fn y(&self) -> &Y

State at last accepted step
Source§

fn t_prev(&self) -> T

Time of previous accepted step
Source§

fn y_prev(&self) -> &Y

State at previous accepted step
Source§

fn h(&self) -> T

Step size for next step
Source§

fn set_h(&mut self, h: T)

Set step size for next step
Source§

fn status(&self) -> &Status<T, Y, D>

Current solver status
Source§

fn set_status(&mut self, status: Status<T, Y, D>)

Set solver status
Source§

impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> OrdinaryNumericalMethod<T, Y, D> for ExplicitRungeKutta<Ordinary, DormandPrince, T, Y, D, O, S, I>

Source§

fn init<F>( &mut self, ode: &F, t0: T, tf: T, y0: &Y, ) -> Result<Evals, Error<T, Y>>
where F: ODE<T, Y, D>,

Initialize the solver before integration Read more
Source§

fn step<F>(&mut self, ode: &F) -> Result<Evals, Error<T, Y>>
where F: ODE<T, Y, D>,

Advance the solution by one step Read more
Source§

fn t(&self) -> T

Time of last accepted step
Source§

fn y(&self) -> &Y

State at last accepted step
Source§

fn t_prev(&self) -> T

Time of previous accepted step
Source§

fn y_prev(&self) -> &Y

State at previous accepted step
Source§

fn h(&self) -> T

Step size for next step
Source§

fn set_h(&mut self, h: T)

Set step size for next step
Source§

fn status(&self) -> &Status<T, Y, D>

Current solver status
Source§

fn set_status(&mut self, status: Status<T, Y, D>)

Set solver status
Source§

impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> OrdinaryNumericalMethod<T, Y, D> for ExplicitRungeKutta<Ordinary, Fixed, T, Y, D, O, S, I>

Source§

fn init<F>( &mut self, ode: &F, t0: T, tf: T, y0: &Y, ) -> Result<Evals, Error<T, Y>>
where F: ODE<T, Y, D>,

Initialize the solver before integration Read more
Source§

fn step<F>(&mut self, ode: &F) -> Result<Evals, Error<T, Y>>
where F: ODE<T, Y, D>,

Advance the solution by one step Read more
Source§

fn t(&self) -> T

Time of last accepted step
Source§

fn y(&self) -> &Y

State at last accepted step
Source§

fn t_prev(&self) -> T

Time of previous accepted step
Source§

fn y_prev(&self) -> &Y

State at previous accepted step
Source§

fn h(&self) -> T

Step size for next step
Source§

fn set_h(&mut self, h: T)

Set step size for next step
Source§

fn status(&self) -> &Status<T, Y, D>

Current solver status
Source§

fn set_status(&mut self, status: Status<T, Y, D>)

Set solver status
Source§

impl<T: Real, Y: State<T>, D: CallBackData, const O: usize, const S: usize, const I: usize> StochasticNumericalMethod<T, Y, D> for ExplicitRungeKutta<Stochastic, Fixed, T, Y, D, O, S, I>

Source§

fn init<F>( &mut self, sde: &mut F, t0: T, tf: T, y0: &Y, ) -> Result<Evals, Error<T, Y>>
where F: SDE<T, Y, D>,

Initialize the solver before integration Read more
Source§

fn step<F>(&mut self, sde: &mut F) -> Result<Evals, Error<T, Y>>
where F: SDE<T, Y, D>,

Advance the solution by one step Read more
Source§

fn t(&self) -> T

Time of last accepted step
Source§

fn y(&self) -> &Y

State at last accepted step
Source§

fn t_prev(&self) -> T

Time of previous accepted step
Source§

fn y_prev(&self) -> &Y

State at previous accepted step
Source§

fn h(&self) -> T

Step size for next step
Source§

fn set_h(&mut self, h: T)

Set step size for next step
Source§

fn status(&self) -> &Status<T, Y, D>

Current solver status
Source§

fn set_status(&mut self, status: Status<T, Y, D>)

Set solver status

Auto Trait Implementations§

§

impl<E, F, T, Y, D, const O: usize, const S: usize, const I: usize> Freeze for ExplicitRungeKutta<E, F, T, Y, D, O, S, I>
where T: Freeze, Y: Freeze, D: Freeze,

§

impl<E, F, T, Y, D, const O: usize, const S: usize, const I: usize> RefUnwindSafe for ExplicitRungeKutta<E, F, T, Y, D, O, S, I>

§

impl<E, F, T, Y, D, const O: usize, const S: usize, const I: usize> Send for ExplicitRungeKutta<E, F, T, Y, D, O, S, I>
where Y: Send, D: Send, F: Send, E: Send,

§

impl<E, F, T, Y, D, const O: usize, const S: usize, const I: usize> Sync for ExplicitRungeKutta<E, F, T, Y, D, O, S, I>
where Y: Sync, D: Sync, F: Sync, E: Sync,

§

impl<E, F, T, Y, D, const O: usize, const S: usize, const I: usize> Unpin for ExplicitRungeKutta<E, F, T, Y, D, O, S, I>
where T: Unpin, Y: Unpin, D: Unpin, F: Unpin, E: Unpin,

§

impl<E, F, T, Y, D, const O: usize, const S: usize, const I: usize> UnwindSafe for ExplicitRungeKutta<E, F, T, Y, D, O, S, I>

Blanket Implementations§

Source§

impl<T> Any for T
where T: 'static + ?Sized,

Source§

fn type_id(&self) -> TypeId

Gets the TypeId of self. Read more
Source§

impl<T> Borrow<T> for T
where T: ?Sized,

Source§

fn borrow(&self) -> &T

Immutably borrows from an owned value. Read more
Source§

impl<T> BorrowMut<T> for T
where T: ?Sized,

Source§

fn borrow_mut(&mut self) -> &mut T

Mutably borrows from an owned value. Read more
Source§

impl<T> From<T> for T

Source§

fn from(t: T) -> T

Returns the argument unchanged.

Source§

impl<T, U> Into<U> for T
where U: From<T>,

Source§

fn into(self) -> U

Calls U::from(self).

That is, this conversion is whatever the implementation of From<T> for U chooses to do.

Source§

impl<T> Same for T

Source§

type Output = T

Should always be Self
Source§

impl<SS, SP> SupersetOf<SS> for SP
where SS: SubsetOf<SP>,

Source§

fn to_subset(&self) -> Option<SS>

The inverse inclusion map: attempts to construct self from the equivalent element of its superset. Read more
Source§

fn is_in_subset(&self) -> bool

Checks if self is actually part of its subset T (and can be converted to it).
Source§

fn to_subset_unchecked(&self) -> SS

Use with care! Same as self.to_subset but without any property checks. Always succeeds.
Source§

fn from_subset(element: &SS) -> SP

The inclusion map: converts self to the equivalent element of its superset.
Source§

impl<T, U> TryFrom<U> for T
where U: Into<T>,

Source§

type Error = Infallible

The type returned in the event of a conversion error.
Source§

fn try_from(value: U) -> Result<T, <T as TryFrom<U>>::Error>

Performs the conversion.
Source§

impl<T, U> TryInto<U> for T
where U: TryFrom<T>,

Source§

type Error = <U as TryFrom<T>>::Error

The type returned in the event of a conversion error.
Source§

fn try_into(self) -> Result<U, <U as TryFrom<T>>::Error>

Performs the conversion.