ImplicitRungeKutta

Struct ImplicitRungeKutta 

Source
pub struct ImplicitRungeKutta<E, F, T: Real, Y: State<T>, const O: usize, const S: usize, const I: usize> {
    pub h0: T,
    pub newton_tol: T,
    pub max_newton_iter: usize,
    pub rtol: Tolerance<T>,
    pub atol: Tolerance<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

IRK solver core. Supports fixed/adaptive stepping and common IRK families (Gauss, Radau, Lobatto).

Type params: E (equation), F (family), T (scalar), Y (state), D (callback), O (order), S (stages), I (dense output terms).

Fields§

§h0: T§newton_tol: T§max_newton_iter: usize§rtol: Tolerance<T>§atol: Tolerance<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>> ImplicitRungeKutta<E, Adaptive, T, Y, 4, 2, 2>

Source

pub fn gauss_legendre_4() -> Self

Creates a new Gauss-Legendre 2-stage implicit Runge-Kutta method of order 4.

Source§

impl<E, T: Real, Y: State<T>> ImplicitRungeKutta<E, Adaptive, T, Y, 6, 3, 3>

Source

pub fn gauss_legendre_6() -> Self

Creates a new Gauss-Legendre 3-stage implicit Runge-Kutta method of order 6.

Source§

impl<E, T: Real, Y: State<T>> ImplicitRungeKutta<E, Adaptive, T, Y, 2, 2, 2>

Source

pub fn lobatto_iiic_2() -> Self

Creates a new Lobatto IIIC 2-stage implicit Runge-Kutta method of order 2.

Source§

impl<E, T: Real, Y: State<T>> ImplicitRungeKutta<E, Adaptive, T, Y, 4, 3, 3>

Source

pub fn lobatto_iiic_4() -> Self

Creates a new Lobatto IIIC 3-stage implicit Runge-Kutta method of order 4.

Source§

impl<E, T: Real, Y: State<T>> ImplicitRungeKutta<E, Fixed, T, Y, 1, 1, 1>

Source

pub fn backward_euler(h0: T) -> Self

Backward Euler, order 1, 1 stage. A-stable; stiff-suitable.

Source§

impl<E, T: Real, Y: State<T>> ImplicitRungeKutta<E, Fixed, T, Y, 2, 2, 2>

Source

pub fn crank_nicolson(h0: T) -> Self

Crank-Nicolson, order 2, 2 stages. A-stable; stiff-suitable.

Source§

impl<E, T: Real, Y: State<T>> ImplicitRungeKutta<E, Fixed, T, Y, 2, 2, 2>

Source

pub fn trapezoidal(h0: T) -> Self

Trapezoidal, order 2, 2 stages. A-stable; stiff-suitable.

Source§

impl<E, T: Real, Y: State<T>> ImplicitRungeKutta<E, Radau, T, Y, 5, 3, 3>

Constructor for Radau5

Source

pub fn radau5() -> Radau5<E, T, Y>

Creates a new Radau IIA 3-stage implicit Runge-Kutta method of order 5.

For full usage details, DAE index handling, tuning notes and examples, see the documentation on the [Radau5] type.

Source§

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

Source

pub fn rtol<V: Into<Tolerance<T>>>(self, rtol: V) -> Self

Set relative tolerance

Source

pub fn atol<V: Into<Tolerance<T>>>(self, atol: V) -> Self

Set absolute tolerance

Source

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

Set initial step size

Source

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

Set minimum step size

Source

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

Set maximum step size

Source

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

Set max steps

Source

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

Set max consecutive rejections

Source

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

Set step size safety factor (default: 0.9)

Source

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

Set minimum scale for step changes (default: 0.2)

Source

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

Set maximum scale for step changes (default: 10.0)

Source

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

Set Newton tolerance (default: 1e-10)

Source

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

Set max Newton iterations per stage (default: 50)

Source

pub fn order(&self) -> usize

Get method order

Source

pub fn stages(&self) -> usize

Get number of stages

Source

pub fn dense_stages(&self) -> usize

Get dense output terms

Trait Implementations§

Source§

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

Source§

fn default() -> Self

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

impl<T: Real, Y: State<T>, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ImplicitRungeKutta<Ordinary, Adaptive, T, Y, 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>, const O: usize, const S: usize, const I: usize> Interpolation<T, Y> for ImplicitRungeKutta<Ordinary, Fixed, T, Y, 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>, const O: usize, const S: usize, const I: usize> OrdinaryNumericalMethod<T, Y> for ImplicitRungeKutta<Ordinary, Adaptive, T, Y, 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>,

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>,

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>

Current solver status
Source§

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

Set solver status
Source§

impl<T: Real, Y: State<T>, const O: usize, const S: usize, const I: usize> OrdinaryNumericalMethod<T, Y> for ImplicitRungeKutta<Ordinary, Fixed, T, Y, 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>,

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>,

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>

Current solver status
Source§

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

Set solver status

Auto Trait Implementations§

§

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

§

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

§

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

§

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

§

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

§

impl<E, F, T, Y, const O: usize, const S: usize, const I: usize> UnwindSafe for ImplicitRungeKutta<E, F, T, Y, 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.