Struct ode_solvers::dop853::Dop853

source ·
pub struct Dop853<T, V, F>
where T: FloatNumber, F: System<T, V>,
{ /* private fields */ }
Expand description

Structure containing the parameters for the numerical integration.

Implementations§

source§

impl<T, D: Dim, F> Dop853<T, OVector<T, D>, F>
where f64: From<T>, T: FloatNumber, F: System<T, OVector<T, D>>, OVector<T, D>: Mul<T, Output = OVector<T, D>>, DefaultAllocator: Allocator<T, D>,

source

pub fn new( f: F, x: T, x_end: T, dx: T, y: OVector<T, D>, rtol: T, atol: T ) -> Self

Default initializer for the structure.

Arguments
  • f - Structure implementing the System trait
  • x - Initial value of the independent variable (usually time)
  • x_end - Final value of the independent variable
  • dx - Increment in the dense output. This argument has no effect if the output type is Sparse
  • y - Initial value of the dependent variable(s)
  • rtol - Relative tolerance used in the computation of the adaptive step size
  • atol - Absolute tolerance used in the computation of the adaptive step size
Examples found in repository?
examples/chemical_reaction.rs (line 18)
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
fn main() {
    // Initial state.
    let y0 = State::new(1.0, 0.0, 0.0);

    // Create the structure containing the ODEs.
    let system = ChemicalReaction;

    // Create a stepper and run the integration.
    let mut stepper = Dop853::new(system, 0., 0.3, 0.3, y0, 1.0e-2, 1.0e-6);
    let res = stepper.integrate();

    // Handle result.
    match res {
        Ok(stats) => println!("{}", stats),
        Err(e) => println!("An error occured: {}", e),
    }
}
More examples
Hide additional examples
examples/lorenz.rs (line 22)
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
fn main() {
    // Initial state
    let y0 = State::new(1.0, 1.0, 1.0);

    // Define problem specific constants
    let system = LorenzAttractor {
        sigma: 10.,
        beta: 8. / 3.,
        rho: 28.,
    };
    // Create stepper and integrate
    let mut stepper = Dop853::new(system, 0.0, 100.0, 1e-3, y0, 1e-4, 1e-4);
    let res = stepper.integrate();

    // Handle result
    match res {
        Ok(stats) => {
            println!("{}", stats);
            let path = Path::new("./outputs/lorenz_dop853.dat");
            save(stepper.x_out(), stepper.y_out(), path);
            println!("Results saved in: {:?}", path);
        }
        Err(e) => println!("An error occured: {}", e),
    }
}
examples/three_body_system.rs (line 20)
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
fn main() {
    // Create the structure containing the problem specific constant and equations.
    let system = ThreeBodyProblem {
        mu: 0.012300118882173,
    };

    // Initial state.
    let y0 = State::new(-0.271, -0.42, 0.0, 0.3, -1.0, 0.0);

    // Create a stepper and run the integration.
    let mut stepper = Dop853::new(system, 0.0, 150.0, 0.002, y0, 1.0e-14, 1.0e-14);
    let res = stepper.integrate();

    // Handle result
    match res {
        Ok(stats) => {
            println!("{}", stats);
            let path = Path::new("./outputs/three_body_dop853.dat");
            save(stepper.x_out(), stepper.y_out(), path);
            println!("Results saved in: {:?}", path);
        }
        Err(e) => println!("An error occured: {}", e),
    }
}
examples/three_body_system_dvector.rs (line 20)
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
fn main() {
    // Create the structure containing the problem specific constant and equations.
    let system = ThreeBodyProblem {
        mu: 0.012300118882173,
    };

    // Initial state.
    let y0 = State::from(vec![-0.271, -0.42, 0.0, 0.3, -1.0, 0.0]);

    // Create a stepper and run the integration.
    let mut stepper = Dop853::new(system, 0.0, 150.0, 0.002, y0, 1.0e-14, 1.0e-14);
    let res = stepper.integrate();

    // Handle result
    match res {
        Ok(stats) => {
            println!("{}", stats);
            let path = Path::new("./outputs/three_body_dop853_dvector.dat");
            save(stepper.x_out(), stepper.y_out(), path);
            println!("Results saved in: {:?}", path);
        }
        Err(e) => println!("An error occured: {}", e),
    }
}
source

pub fn from_param( f: F, x: T, x_end: T, dx: T, y: OVector<T, D>, rtol: T, atol: T, safety_factor: T, beta: T, fac_min: T, fac_max: T, h_max: T, h: T, n_max: u32, n_stiff: u32, out_type: OutputType ) -> Self

Advanced initializer for the structure.

Arguments
  • f - Structure implementing the System trait
  • x - Initial value of the independent variable (usually time)
  • x_end - Final value of the independent variable
  • dx - Increment in the dense output. This argument has no effect if the output type is Sparse
  • y - Initial value of the dependent variable(s)
  • rtol - Relative tolerance used in the computation of the adaptive step size
  • atol - Absolute tolerance used in the computation of the adaptive step size
  • safety_factor - Safety factor used in the computation of the adaptive step size. Default is 0.9
  • beta - Value of the beta coefficient of the PI controller. Default is 0.0
  • fac_min - Minimum factor between two successive steps. Default is 0.333
  • fac_max - Maximum factor between two successive steps. Default is 6.0
  • h_max - Maximum step size. Default is `x_end-x
  • h - Initial value of the step size. If h = 0.0, the intial value of h is computed automatically
  • n_max - Maximum number of iterations. Default is 100000
  • n_stiff - Stifness is tested when the number of iterations is a multiple of n_stiff. Default is 1000
  • out_type - Type of the output. Must be a variant of the OutputType enum. Default is Dense
source

pub fn integrate(&mut self) -> Result<Stats, IntegrationError>

Core integration method.

Examples found in repository?
examples/chemical_reaction.rs (line 19)
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
fn main() {
    // Initial state.
    let y0 = State::new(1.0, 0.0, 0.0);

    // Create the structure containing the ODEs.
    let system = ChemicalReaction;

    // Create a stepper and run the integration.
    let mut stepper = Dop853::new(system, 0., 0.3, 0.3, y0, 1.0e-2, 1.0e-6);
    let res = stepper.integrate();

    // Handle result.
    match res {
        Ok(stats) => println!("{}", stats),
        Err(e) => println!("An error occured: {}", e),
    }
}
More examples
Hide additional examples
examples/lorenz.rs (line 23)
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
fn main() {
    // Initial state
    let y0 = State::new(1.0, 1.0, 1.0);

    // Define problem specific constants
    let system = LorenzAttractor {
        sigma: 10.,
        beta: 8. / 3.,
        rho: 28.,
    };
    // Create stepper and integrate
    let mut stepper = Dop853::new(system, 0.0, 100.0, 1e-3, y0, 1e-4, 1e-4);
    let res = stepper.integrate();

    // Handle result
    match res {
        Ok(stats) => {
            println!("{}", stats);
            let path = Path::new("./outputs/lorenz_dop853.dat");
            save(stepper.x_out(), stepper.y_out(), path);
            println!("Results saved in: {:?}", path);
        }
        Err(e) => println!("An error occured: {}", e),
    }
}
examples/three_body_system.rs (line 21)
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
fn main() {
    // Create the structure containing the problem specific constant and equations.
    let system = ThreeBodyProblem {
        mu: 0.012300118882173,
    };

    // Initial state.
    let y0 = State::new(-0.271, -0.42, 0.0, 0.3, -1.0, 0.0);

    // Create a stepper and run the integration.
    let mut stepper = Dop853::new(system, 0.0, 150.0, 0.002, y0, 1.0e-14, 1.0e-14);
    let res = stepper.integrate();

    // Handle result
    match res {
        Ok(stats) => {
            println!("{}", stats);
            let path = Path::new("./outputs/three_body_dop853.dat");
            save(stepper.x_out(), stepper.y_out(), path);
            println!("Results saved in: {:?}", path);
        }
        Err(e) => println!("An error occured: {}", e),
    }
}
examples/three_body_system_dvector.rs (line 21)
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
fn main() {
    // Create the structure containing the problem specific constant and equations.
    let system = ThreeBodyProblem {
        mu: 0.012300118882173,
    };

    // Initial state.
    let y0 = State::from(vec![-0.271, -0.42, 0.0, 0.3, -1.0, 0.0]);

    // Create a stepper and run the integration.
    let mut stepper = Dop853::new(system, 0.0, 150.0, 0.002, y0, 1.0e-14, 1.0e-14);
    let res = stepper.integrate();

    // Handle result
    match res {
        Ok(stats) => {
            println!("{}", stats);
            let path = Path::new("./outputs/three_body_dop853_dvector.dat");
            save(stepper.x_out(), stepper.y_out(), path);
            println!("Results saved in: {:?}", path);
        }
        Err(e) => println!("An error occured: {}", e),
    }
}
source

pub fn x_out(&self) -> &Vec<T>

Getter for the independent variable’s output.

Examples found in repository?
examples/lorenz.rs (line 30)
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
fn main() {
    // Initial state
    let y0 = State::new(1.0, 1.0, 1.0);

    // Define problem specific constants
    let system = LorenzAttractor {
        sigma: 10.,
        beta: 8. / 3.,
        rho: 28.,
    };
    // Create stepper and integrate
    let mut stepper = Dop853::new(system, 0.0, 100.0, 1e-3, y0, 1e-4, 1e-4);
    let res = stepper.integrate();

    // Handle result
    match res {
        Ok(stats) => {
            println!("{}", stats);
            let path = Path::new("./outputs/lorenz_dop853.dat");
            save(stepper.x_out(), stepper.y_out(), path);
            println!("Results saved in: {:?}", path);
        }
        Err(e) => println!("An error occured: {}", e),
    }
}
More examples
Hide additional examples
examples/three_body_system.rs (line 28)
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
fn main() {
    // Create the structure containing the problem specific constant and equations.
    let system = ThreeBodyProblem {
        mu: 0.012300118882173,
    };

    // Initial state.
    let y0 = State::new(-0.271, -0.42, 0.0, 0.3, -1.0, 0.0);

    // Create a stepper and run the integration.
    let mut stepper = Dop853::new(system, 0.0, 150.0, 0.002, y0, 1.0e-14, 1.0e-14);
    let res = stepper.integrate();

    // Handle result
    match res {
        Ok(stats) => {
            println!("{}", stats);
            let path = Path::new("./outputs/three_body_dop853.dat");
            save(stepper.x_out(), stepper.y_out(), path);
            println!("Results saved in: {:?}", path);
        }
        Err(e) => println!("An error occured: {}", e),
    }
}
examples/three_body_system_dvector.rs (line 28)
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
fn main() {
    // Create the structure containing the problem specific constant and equations.
    let system = ThreeBodyProblem {
        mu: 0.012300118882173,
    };

    // Initial state.
    let y0 = State::from(vec![-0.271, -0.42, 0.0, 0.3, -1.0, 0.0]);

    // Create a stepper and run the integration.
    let mut stepper = Dop853::new(system, 0.0, 150.0, 0.002, y0, 1.0e-14, 1.0e-14);
    let res = stepper.integrate();

    // Handle result
    match res {
        Ok(stats) => {
            println!("{}", stats);
            let path = Path::new("./outputs/three_body_dop853_dvector.dat");
            save(stepper.x_out(), stepper.y_out(), path);
            println!("Results saved in: {:?}", path);
        }
        Err(e) => println!("An error occured: {}", e),
    }
}
source

pub fn y_out(&self) -> &Vec<OVector<T, D>>

Getter for the dependent variables’ output.

Examples found in repository?
examples/lorenz.rs (line 30)
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
fn main() {
    // Initial state
    let y0 = State::new(1.0, 1.0, 1.0);

    // Define problem specific constants
    let system = LorenzAttractor {
        sigma: 10.,
        beta: 8. / 3.,
        rho: 28.,
    };
    // Create stepper and integrate
    let mut stepper = Dop853::new(system, 0.0, 100.0, 1e-3, y0, 1e-4, 1e-4);
    let res = stepper.integrate();

    // Handle result
    match res {
        Ok(stats) => {
            println!("{}", stats);
            let path = Path::new("./outputs/lorenz_dop853.dat");
            save(stepper.x_out(), stepper.y_out(), path);
            println!("Results saved in: {:?}", path);
        }
        Err(e) => println!("An error occured: {}", e),
    }
}
More examples
Hide additional examples
examples/three_body_system.rs (line 28)
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
fn main() {
    // Create the structure containing the problem specific constant and equations.
    let system = ThreeBodyProblem {
        mu: 0.012300118882173,
    };

    // Initial state.
    let y0 = State::new(-0.271, -0.42, 0.0, 0.3, -1.0, 0.0);

    // Create a stepper and run the integration.
    let mut stepper = Dop853::new(system, 0.0, 150.0, 0.002, y0, 1.0e-14, 1.0e-14);
    let res = stepper.integrate();

    // Handle result
    match res {
        Ok(stats) => {
            println!("{}", stats);
            let path = Path::new("./outputs/three_body_dop853.dat");
            save(stepper.x_out(), stepper.y_out(), path);
            println!("Results saved in: {:?}", path);
        }
        Err(e) => println!("An error occured: {}", e),
    }
}
examples/three_body_system_dvector.rs (line 28)
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
fn main() {
    // Create the structure containing the problem specific constant and equations.
    let system = ThreeBodyProblem {
        mu: 0.012300118882173,
    };

    // Initial state.
    let y0 = State::from(vec![-0.271, -0.42, 0.0, 0.3, -1.0, 0.0]);

    // Create a stepper and run the integration.
    let mut stepper = Dop853::new(system, 0.0, 150.0, 0.002, y0, 1.0e-14, 1.0e-14);
    let res = stepper.integrate();

    // Handle result
    match res {
        Ok(stats) => {
            println!("{}", stats);
            let path = Path::new("./outputs/three_body_dop853_dvector.dat");
            save(stepper.x_out(), stepper.y_out(), path);
            println!("Results saved in: {:?}", path);
        }
        Err(e) => println!("An error occured: {}", e),
    }
}
source

pub fn results(&self) -> &SolverResult<T, OVector<T, D>>

Getter for the results type, a pair of independent and dependent variables

Trait Implementations§

source§

impl<T, D: Dim, F> Into<SolverResult<T, Matrix<T, D, Const<1>, <DefaultAllocator as Allocator<T, D>>::Buffer>>> for Dop853<T, OVector<T, D>, F>
where T: FloatNumber, F: System<T, OVector<T, D>>, DefaultAllocator: Allocator<T, D>,

source§

fn into(self) -> SolverResult<T, OVector<T, D>>

Converts this type into the (usually inferred) input type.

Auto Trait Implementations§

§

impl<T, V, F> RefUnwindSafe for Dop853<T, V, F>

§

impl<T, V, F> Send for Dop853<T, V, F>
where F: Send, T: Send, V: Send,

§

impl<T, V, F> Sync for Dop853<T, V, F>
where F: Sync, T: Sync, V: Sync,

§

impl<T, V, F> Unpin for Dop853<T, V, F>
where F: Unpin, T: Unpin, V: Unpin,

§

impl<T, V, F> UnwindSafe for Dop853<T, V, F>
where F: UnwindSafe, T: UnwindSafe, V: UnwindSafe,

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

§

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

§

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

§

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.