Struct AdamsPredictorCorrector

Source
pub struct AdamsPredictorCorrector<E, F, T: Real, V: State<T>, D: CallBackData, const S: usize> {
    pub h0: T,
    pub tol: T,
    pub h_max: T,
    pub h_min: T,
    pub max_steps: usize,
    pub stages: usize,
    /* private fields */
}

Fields§

§h0: T§tol: T§h_max: T§h_min: T§max_steps: usize§stages: usize

Implementations§

Source§

impl<T: Real, V: State<T>, D: CallBackData> AdamsPredictorCorrector<Ordinary, Adaptive, T, V, D, 4>

Source

pub fn v4() -> Self

The Adams-Predictor-Corrector method is an explicit method that uses the previous states to predict the next state. This implementation uses a variable step size to maintain a desired accuracy. It is recommended to start with a small step size so that tolerance can be quickly met and the algorithm can adjust the step size accordingly.

The First 3 steps are calculated using the Runge-Kutta method of order 4(5) and then the Adams-Predictor-Corrector method is used to calculate the remaining steps until the final time./ Create a Adams-Predictor-Corrector 4th Order Variable Step Size Method instance.

§Example
use differential_equations::prelude::*;
use nalgebra::{SVector, vector};

struct HarmonicOscillator {
    k: f64,
}

impl ODE<f64, SVector<f64, 2>> for HarmonicOscillator {
    fn diff(&self, _t: f64, y: &SVector<f64, 2>, dydt: &mut SVector<f64, 2>) {
        dydt[0] = y[1];
        dydt[1] = -self.k * y[0];
    }
}
let mut apcv4 = AdamsPredictorCorrector::v4();
let t0 = 0.0;
let tf = 10.0;
let y0 = vector![1.0, 0.0];
let system = HarmonicOscillator { k: 1.0 };
let results = ODEProblem::new(system, t0, tf, y0).solve(&mut apcv4).unwrap();
let expected = vector![-0.83907153, 0.54402111];
assert!((results.y.last().unwrap()[0] - expected[0]).abs() < 1e-6);
assert!((results.y.last().unwrap()[1] - expected[1]).abs() < 1e-6);
§Warning

This method is not suitable for stiff problems and can results in extremely small step sizes and long computation times.```

Source§

impl<T: Real, V: State<T>, D: CallBackData> AdamsPredictorCorrector<Ordinary, Fixed, T, V, D, 4>

Source

pub fn f4(h: T) -> Self

Adams-Predictor-Corrector 4th Order Fixed Step Size Method.

The Adams-Predictor-Corrector method is an explicit method that uses the previous states to predict the next state.

The First 3 steps, of fixed step size h, are calculated using the Runge-Kutta method of order 4(5) and then the Adams-Predictor-Corrector method is used to calculate the remaining steps until the final time.

§Example
use differential_equations::prelude::*;
use nalgebra::{SVector, vector};

struct HarmonicOscillator {
    k: f64,
}

impl ODE<f64, SVector<f64, 2>> for HarmonicOscillator {
    fn diff(&self, _t: f64, y: &SVector<f64, 2>, dydt: &mut SVector<f64, 2>) {
        dydt[0] = y[1];
        dydt[1] = -self.k * y[0];
    }
}
let mut apcf4 = AdamsPredictorCorrector::f4(0.01);
let t0 = 0.0;
let tf = 10.0;
let y0 = vector![1.0, 0.0];
let system = HarmonicOscillator { k: 1.0 };
let results = ODEProblem::new(system, t0, tf, y0).solve(&mut apcf4).unwrap();
let expected = vector![-0.83907153, 0.54402111];
assert!((results.y.last().unwrap()[0] - expected[0]).abs() < 1e-2);
assert!((results.y.last().unwrap()[1] - expected[1]).abs() < 1e-2);
§Settings
  • h - Step Size
Source§

impl<E, F, T: Real, V: State<T>, D: CallBackData, const S: usize> AdamsPredictorCorrector<E, F, T, V, D, S>

Source

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

Set the 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 stages(&self) -> usize

Get the number of stages in the method

Trait Implementations§

Source§

impl<E, F, T: Real, V: State<T>, D: CallBackData, const S: usize> Default for AdamsPredictorCorrector<E, F, T, V, D, S>

Source§

fn default() -> Self

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

impl<T: Real, V: State<T>, D: CallBackData> Interpolation<T, V> for AdamsPredictorCorrector<Ordinary, Adaptive, T, V, D, 4>

Source§

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

Interpolate between previous and current step Read more
Source§

impl<T: Real, V: State<T>, D: CallBackData> Interpolation<T, V> for AdamsPredictorCorrector<Ordinary, Fixed, T, V, D, 4>

Source§

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

Interpolate between previous and current step Read more
Source§

impl<T: Real, V: State<T>, D: CallBackData> OrdinaryNumericalMethod<T, V, D> for AdamsPredictorCorrector<Ordinary, Adaptive, T, V, D, 4>

Source§

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

Initialize OrdinaryNumericalMethod before solving ODE Read more
Source§

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

Step through solving the ODE by one step Read more
Source§

fn t(&self) -> T

Access time of last accepted step
Source§

fn y(&self) -> &V

Access solution of last accepted step
Source§

fn t_prev(&self) -> T

Access time of previous accepted step
Source§

fn y_prev(&self) -> &V

Access solution of previous accepted step
Source§

fn h(&self) -> T

Access step size of next step
Source§

fn set_h(&mut self, h: T)

Set step size of next step
Source§

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

Status of solver
Source§

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

Set status of solver
Source§

impl<T: Real, V: State<T>, D: CallBackData> OrdinaryNumericalMethod<T, V, D> for AdamsPredictorCorrector<Ordinary, Fixed, T, V, D, 4>

Source§

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

Initialize OrdinaryNumericalMethod before solving ODE Read more
Source§

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

Step through solving the ODE by one step Read more
Source§

fn t(&self) -> T

Access time of last accepted step
Source§

fn y(&self) -> &V

Access solution of last accepted step
Source§

fn t_prev(&self) -> T

Access time of previous accepted step
Source§

fn y_prev(&self) -> &V

Access solution of previous accepted step
Source§

fn h(&self) -> T

Access step size of next step
Source§

fn set_h(&mut self, h: T)

Set step size of next step
Source§

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

Status of solver
Source§

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

Set status of solver

Auto Trait Implementations§

§

impl<E, F, T, V, D, const S: usize> Freeze for AdamsPredictorCorrector<E, F, T, V, D, S>
where T: Freeze, V: Freeze, D: Freeze,

§

impl<E, F, T, V, D, const S: usize> RefUnwindSafe for AdamsPredictorCorrector<E, F, T, V, D, S>

§

impl<E, F, T, V, D, const S: usize> Send for AdamsPredictorCorrector<E, F, T, V, D, S>
where V: Send, D: Send, F: Send, E: Send,

§

impl<E, F, T, V, D, const S: usize> Sync for AdamsPredictorCorrector<E, F, T, V, D, S>
where V: Sync, D: Sync, F: Sync, E: Sync,

§

impl<E, F, T, V, D, const S: usize> Unpin for AdamsPredictorCorrector<E, F, T, V, D, S>
where T: Unpin, V: Unpin, D: Unpin, F: Unpin, E: Unpin,

§

impl<E, F, T, V, D, const S: usize> UnwindSafe for AdamsPredictorCorrector<E, F, T, V, D, S>

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.