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>,
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>,
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/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
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),
}
}
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 size. Default is 0.9beta
- Value of the beta coefficient of the PI controller. Default is 0.0fac_min
- Minimum factor between two successive steps. Default is 0.333fac_max
- Maximum factor between two successive steps. Default is 6.0h_max
- Maximum step size. Default is `x_end-xh
- 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/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
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),
}
}
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/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
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),
}
}
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/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
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),
}
}
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 Dop853<T, V, F>
impl<T, V, F> Send for Dop853<T, V, F>
impl<T, V, F> Sync for Dop853<T, V, F>
impl<T, V, F> Unpin for Dop853<T, V, F>
impl<T, V, F> UnwindSafe for Dop853<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.