Struct ode_solvers::dopri5::Dopri5
source · pub struct Dopri5<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> Dopri5<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>,
impl<T, D: Dim, F> Dopri5<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>,
sourcepub fn new(
f: F,
x: T,
x_end: T,
dx: T,
y: OVector<T, D>,
rtol: T,
atol: T
) -> Self
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 Systemtrait x
- Initial value of the independent variable (usually time)x_end
- Final value of the independent variabledx
- Increment in the dense output. This argument has no effect if the output type is Sparsey
- Initial value of the dependent variable(s)rtol
- Relative tolerance used in the computation of the adaptive step sizeatol
- Absolute tolerance used in the computation of the adaptive step size
Examples found in repository?
examples/kepler_orbit.rs (line 30)
12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43
fn main() {
// Create the structure containing the ODEs.
let system = KeplerOrbit { mu: 398600.435436 };
let a: f64 = 20000.0;
let period = 2.0 * PI * (a.powi(3) / system.mu).sqrt();
// Orbit with: a = 20000km, e = 0.7, i = 35 deg, raan = 100 deg, arg_per = 65 deg, true_an = 30 deg.
let y0 = State::new(
-5007.248417988539,
-1444.918140151374,
3628.534606178356,
0.717716656891,
-10.224093784269,
0.748229399696,
);
// Create a stepper and run the integration.
let mut stepper = Dopri5::new(system, 0.0, 5.0 * period, 60.0, y0, 1.0e-10, 1.0e-10);
let res = stepper.integrate();
// Handle result.
match res {
Ok(stats) => {
println!("{}", stats);
let path = Path::new("./outputs/kepler_orbit_dopri5.dat");
save(stepper.x_out(), stepper.y_out(), path);
println!("Results saved in: {:?}", path);
}
Err(e) => println!("An error occured: {}", e),
}
}
More examples
examples/kepler_orbit_dvector.rs (line 30)
12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43
fn main() {
// Create the structure containing the ODEs.
let system = KeplerOrbit { mu: 398600.435436 };
let a: f64 = 20000.0;
let period = 2.0 * PI * (a.powi(3) / system.mu).sqrt();
// Orbit with: a = 20000km, e = 0.7, i = 35 deg, raan = 100 deg, arg_per = 65 deg, true_an = 30 deg.
let y0 = State::from(vec![
-5007.248417988539,
-1444.918140151374,
3628.534606178356,
0.717716656891,
-10.224093784269,
0.748229399696,
]);
// Create a stepper and run the integration.
let mut stepper = Dopri5::new(system, 0.0, 5.0 * period, 60.0, y0, 1.0e-10, 1.0e-10);
let res = stepper.integrate();
// Handle result.
match res {
Ok(stats) => {
println!("{}", stats);
let path = Path::new("./outputs/kepler_orbit_dopri5_dvector.dat");
save(stepper.x_out(), stepper.y_out(), path);
println!("Results saved in: {:?}", path);
}
Err(e) => println!("An error occured: {}", e),
}
}
sourcepub 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
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 Systemtrait x
- Initial value of the independent variable (usually time)x_end
- Final value of the independent variabledx
- Increment in the dense output. This argument has no effect if the output type is Sparsey
- Initial value of the dependent variable(s)rtol
- Relative tolerance used in the computation of the adaptive step sizeatol
- Absolute tolerance used in the computation of the adaptive step sizesafety_factor
- Safety factor used in the computation of the adaptive step sizebeta
- Value of the beta coefficient of the PI controller. Default is 0.04fac_min
- Minimum factor between two successive steps. Default is 0.2fac_max
- Maximum factor between two successive steps. Default is 10.0h_max
- Maximum step size. Default isx_end-x
h
- Initial value of the step size. If h = 0.0, the intial value of h is computed automaticallyn_max
- Maximum number of iterations. Default is 100000n_stiff
- Stifness is tested when the number of iterations is a multiple of n_stiff. Default is 1000out_type
- Type of the output. Must be a variant of the OutputType enum. Default is Dense
sourcepub fn integrate(&mut self) -> Result<Stats, IntegrationError>
pub fn integrate(&mut self) -> Result<Stats, IntegrationError>
Core integration method.
Examples found in repository?
examples/kepler_orbit.rs (line 31)
12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43
fn main() {
// Create the structure containing the ODEs.
let system = KeplerOrbit { mu: 398600.435436 };
let a: f64 = 20000.0;
let period = 2.0 * PI * (a.powi(3) / system.mu).sqrt();
// Orbit with: a = 20000km, e = 0.7, i = 35 deg, raan = 100 deg, arg_per = 65 deg, true_an = 30 deg.
let y0 = State::new(
-5007.248417988539,
-1444.918140151374,
3628.534606178356,
0.717716656891,
-10.224093784269,
0.748229399696,
);
// Create a stepper and run the integration.
let mut stepper = Dopri5::new(system, 0.0, 5.0 * period, 60.0, y0, 1.0e-10, 1.0e-10);
let res = stepper.integrate();
// Handle result.
match res {
Ok(stats) => {
println!("{}", stats);
let path = Path::new("./outputs/kepler_orbit_dopri5.dat");
save(stepper.x_out(), stepper.y_out(), path);
println!("Results saved in: {:?}", path);
}
Err(e) => println!("An error occured: {}", e),
}
}
More examples
examples/kepler_orbit_dvector.rs (line 31)
12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43
fn main() {
// Create the structure containing the ODEs.
let system = KeplerOrbit { mu: 398600.435436 };
let a: f64 = 20000.0;
let period = 2.0 * PI * (a.powi(3) / system.mu).sqrt();
// Orbit with: a = 20000km, e = 0.7, i = 35 deg, raan = 100 deg, arg_per = 65 deg, true_an = 30 deg.
let y0 = State::from(vec![
-5007.248417988539,
-1444.918140151374,
3628.534606178356,
0.717716656891,
-10.224093784269,
0.748229399696,
]);
// Create a stepper and run the integration.
let mut stepper = Dopri5::new(system, 0.0, 5.0 * period, 60.0, y0, 1.0e-10, 1.0e-10);
let res = stepper.integrate();
// Handle result.
match res {
Ok(stats) => {
println!("{}", stats);
let path = Path::new("./outputs/kepler_orbit_dopri5_dvector.dat");
save(stepper.x_out(), stepper.y_out(), path);
println!("Results saved in: {:?}", path);
}
Err(e) => println!("An error occured: {}", e),
}
}
sourcepub fn x_out(&self) -> &Vec<T>
pub fn x_out(&self) -> &Vec<T>
Getter for the independent variable’s output.
Examples found in repository?
examples/kepler_orbit.rs (line 38)
12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43
fn main() {
// Create the structure containing the ODEs.
let system = KeplerOrbit { mu: 398600.435436 };
let a: f64 = 20000.0;
let period = 2.0 * PI * (a.powi(3) / system.mu).sqrt();
// Orbit with: a = 20000km, e = 0.7, i = 35 deg, raan = 100 deg, arg_per = 65 deg, true_an = 30 deg.
let y0 = State::new(
-5007.248417988539,
-1444.918140151374,
3628.534606178356,
0.717716656891,
-10.224093784269,
0.748229399696,
);
// Create a stepper and run the integration.
let mut stepper = Dopri5::new(system, 0.0, 5.0 * period, 60.0, y0, 1.0e-10, 1.0e-10);
let res = stepper.integrate();
// Handle result.
match res {
Ok(stats) => {
println!("{}", stats);
let path = Path::new("./outputs/kepler_orbit_dopri5.dat");
save(stepper.x_out(), stepper.y_out(), path);
println!("Results saved in: {:?}", path);
}
Err(e) => println!("An error occured: {}", e),
}
}
More examples
examples/kepler_orbit_dvector.rs (line 38)
12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43
fn main() {
// Create the structure containing the ODEs.
let system = KeplerOrbit { mu: 398600.435436 };
let a: f64 = 20000.0;
let period = 2.0 * PI * (a.powi(3) / system.mu).sqrt();
// Orbit with: a = 20000km, e = 0.7, i = 35 deg, raan = 100 deg, arg_per = 65 deg, true_an = 30 deg.
let y0 = State::from(vec![
-5007.248417988539,
-1444.918140151374,
3628.534606178356,
0.717716656891,
-10.224093784269,
0.748229399696,
]);
// Create a stepper and run the integration.
let mut stepper = Dopri5::new(system, 0.0, 5.0 * period, 60.0, y0, 1.0e-10, 1.0e-10);
let res = stepper.integrate();
// Handle result.
match res {
Ok(stats) => {
println!("{}", stats);
let path = Path::new("./outputs/kepler_orbit_dopri5_dvector.dat");
save(stepper.x_out(), stepper.y_out(), path);
println!("Results saved in: {:?}", path);
}
Err(e) => println!("An error occured: {}", e),
}
}
sourcepub fn y_out(&self) -> &Vec<OVector<T, D>>
pub fn y_out(&self) -> &Vec<OVector<T, D>>
Getter for the dependent variables’ output.
Examples found in repository?
examples/kepler_orbit.rs (line 38)
12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43
fn main() {
// Create the structure containing the ODEs.
let system = KeplerOrbit { mu: 398600.435436 };
let a: f64 = 20000.0;
let period = 2.0 * PI * (a.powi(3) / system.mu).sqrt();
// Orbit with: a = 20000km, e = 0.7, i = 35 deg, raan = 100 deg, arg_per = 65 deg, true_an = 30 deg.
let y0 = State::new(
-5007.248417988539,
-1444.918140151374,
3628.534606178356,
0.717716656891,
-10.224093784269,
0.748229399696,
);
// Create a stepper and run the integration.
let mut stepper = Dopri5::new(system, 0.0, 5.0 * period, 60.0, y0, 1.0e-10, 1.0e-10);
let res = stepper.integrate();
// Handle result.
match res {
Ok(stats) => {
println!("{}", stats);
let path = Path::new("./outputs/kepler_orbit_dopri5.dat");
save(stepper.x_out(), stepper.y_out(), path);
println!("Results saved in: {:?}", path);
}
Err(e) => println!("An error occured: {}", e),
}
}
More examples
examples/kepler_orbit_dvector.rs (line 38)
12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43
fn main() {
// Create the structure containing the ODEs.
let system = KeplerOrbit { mu: 398600.435436 };
let a: f64 = 20000.0;
let period = 2.0 * PI * (a.powi(3) / system.mu).sqrt();
// Orbit with: a = 20000km, e = 0.7, i = 35 deg, raan = 100 deg, arg_per = 65 deg, true_an = 30 deg.
let y0 = State::from(vec![
-5007.248417988539,
-1444.918140151374,
3628.534606178356,
0.717716656891,
-10.224093784269,
0.748229399696,
]);
// Create a stepper and run the integration.
let mut stepper = Dopri5::new(system, 0.0, 5.0 * period, 60.0, y0, 1.0e-10, 1.0e-10);
let res = stepper.integrate();
// Handle result.
match res {
Ok(stats) => {
println!("{}", stats);
let path = Path::new("./outputs/kepler_orbit_dopri5_dvector.dat");
save(stepper.x_out(), stepper.y_out(), path);
println!("Results saved in: {:?}", path);
}
Err(e) => println!("An error occured: {}", e),
}
}
sourcepub fn results(&self) -> &SolverResult<T, OVector<T, D>>
pub fn results(&self) -> &SolverResult<T, OVector<T, D>>
Getter for the results type, a pair of independent and dependent variables
Trait Implementations§
Auto Trait Implementations§
impl<T, V, F> RefUnwindSafe for Dopri5<T, V, F>
impl<T, V, F> Send for Dopri5<T, V, F>
impl<T, V, F> Sync for Dopri5<T, V, F>
impl<T, V, F> Unpin for Dopri5<T, V, F>
impl<T, V, F> UnwindSafe for Dopri5<T, V, F>
Blanket Implementations§
source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Mutably borrows from an owned value. Read more
source§impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
source§fn to_subset(&self) -> Option<SS>
fn to_subset(&self) -> Option<SS>
The inverse inclusion map: attempts to construct
self
from the equivalent element of its
superset. Read moresource§fn is_in_subset(&self) -> bool
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
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
fn from_subset(element: &SS) -> SP
The inclusion map: converts
self
to the equivalent element of its superset.