Struct ode_solvers::rk4::Rk4
source · pub struct Rk4<T, V, F>where
F: System<T, V>,
T: FloatNumber,{ /* private fields */ }
Expand description
Structure containing the parameters for the numerical integration.
Implementations§
source§impl<T, D: Dim, F> Rk4<T, OVector<T, D>, F>where
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> Rk4<T, OVector<T, D>, F>where
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, y: OVector<T, D>, x_end: T, step_size: T) -> Self
pub fn new(f: F, x: T, y: OVector<T, D>, x_end: T, step_size: T) -> Self
Default initializer for the structure
Arguments
f
- Structure implementing the Systemtrait x
- Initial value of the independent variable (usually time)y
- Initial value of the dependent variable(s)x_end
- Final value of the independent variablestep_size
- Step size used in the method
Examples found in repository?
examples/bouncing_ball.rs (line 31)
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 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68
fn main() {
// Initial state: At 10m with zero velocity
let mut y0 = State::new(10.0, 0., 0.);
let mut num_bounces = 0;
let mut combined_solver_results = Result::default();
while num_bounces < MAX_BOUNCES && y0[0] >= 0.001 {
// Create the structure containing the ODEs.
let system = BouncingBall;
// Create a stepper and run the integration.
// Use comments to see differences with Dopri
//let mut stepper = Dopri5::new(system, 0., 10.0, 0.01, y0, 1.0e-2, 1.0e-6);
let mut stepper = Rk4::new(system, 0f32, y0, 10f32, 0.01f32);
let res = stepper.integrate();
// Handle result.
match res {
Ok(stats) => println!("{}", stats),
Err(e) => println!("An error occured: {}", e),
}
num_bounces = num_bounces + 1;
// solout may not be called and therefore end not "smooth" when observing dense values with dopri5 or dop853
// Therefore we seach for the point where the results turn zero
let (_, y_out) = stepper.results().get();
let f = y_out.iter().find(|y| y[0] <= 0.);
if f.is_none() {
// that should not happen...
break;
}
let last_state = f.unwrap();
println!("Last state: {:?}", last_state);
y0[0] = last_state[0].abs();
y0[1] = -1. * last_state[1] * BOUNCE;
// beware in the case of dopri5 or dop853 the results contain a lot of invalid data with y[0] < 0
combined_solver_results.append(stepper.into());
}
let path = Path::new("./outputs/bouncing_ball.dat");
save(
combined_solver_results.get().0,
combined_solver_results.get().1,
path,
);
}
sourcepub fn integrate(&mut self) -> Result<Stats, IntegrationError>
pub fn integrate(&mut self) -> Result<Stats, IntegrationError>
Core integration method.
Examples found in repository?
examples/bouncing_ball.rs (line 32)
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 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68
fn main() {
// Initial state: At 10m with zero velocity
let mut y0 = State::new(10.0, 0., 0.);
let mut num_bounces = 0;
let mut combined_solver_results = Result::default();
while num_bounces < MAX_BOUNCES && y0[0] >= 0.001 {
// Create the structure containing the ODEs.
let system = BouncingBall;
// Create a stepper and run the integration.
// Use comments to see differences with Dopri
//let mut stepper = Dopri5::new(system, 0., 10.0, 0.01, y0, 1.0e-2, 1.0e-6);
let mut stepper = Rk4::new(system, 0f32, y0, 10f32, 0.01f32);
let res = stepper.integrate();
// Handle result.
match res {
Ok(stats) => println!("{}", stats),
Err(e) => println!("An error occured: {}", e),
}
num_bounces = num_bounces + 1;
// solout may not be called and therefore end not "smooth" when observing dense values with dopri5 or dop853
// Therefore we seach for the point where the results turn zero
let (_, y_out) = stepper.results().get();
let f = y_out.iter().find(|y| y[0] <= 0.);
if f.is_none() {
// that should not happen...
break;
}
let last_state = f.unwrap();
println!("Last state: {:?}", last_state);
y0[0] = last_state[0].abs();
y0[1] = -1. * last_state[1] * BOUNCE;
// beware in the case of dopri5 or dop853 the results contain a lot of invalid data with y[0] < 0
combined_solver_results.append(stepper.into());
}
let path = Path::new("./outputs/bouncing_ball.dat");
save(
combined_solver_results.get().0,
combined_solver_results.get().1,
path,
);
}
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
Examples found in repository?
examples/bouncing_ball.rs (line 44)
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 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68
fn main() {
// Initial state: At 10m with zero velocity
let mut y0 = State::new(10.0, 0., 0.);
let mut num_bounces = 0;
let mut combined_solver_results = Result::default();
while num_bounces < MAX_BOUNCES && y0[0] >= 0.001 {
// Create the structure containing the ODEs.
let system = BouncingBall;
// Create a stepper and run the integration.
// Use comments to see differences with Dopri
//let mut stepper = Dopri5::new(system, 0., 10.0, 0.01, y0, 1.0e-2, 1.0e-6);
let mut stepper = Rk4::new(system, 0f32, y0, 10f32, 0.01f32);
let res = stepper.integrate();
// Handle result.
match res {
Ok(stats) => println!("{}", stats),
Err(e) => println!("An error occured: {}", e),
}
num_bounces = num_bounces + 1;
// solout may not be called and therefore end not "smooth" when observing dense values with dopri5 or dop853
// Therefore we seach for the point where the results turn zero
let (_, y_out) = stepper.results().get();
let f = y_out.iter().find(|y| y[0] <= 0.);
if f.is_none() {
// that should not happen...
break;
}
let last_state = f.unwrap();
println!("Last state: {:?}", last_state);
y0[0] = last_state[0].abs();
y0[1] = -1. * last_state[1] * BOUNCE;
// beware in the case of dopri5 or dop853 the results contain a lot of invalid data with y[0] < 0
combined_solver_results.append(stepper.into());
}
let path = Path::new("./outputs/bouncing_ball.dat");
save(
combined_solver_results.get().0,
combined_solver_results.get().1,
path,
);
}
Trait Implementations§
Auto Trait Implementations§
impl<T, V, F> RefUnwindSafe for Rk4<T, V, F>
impl<T, V, F> Send for Rk4<T, V, F>
impl<T, V, F> Sync for Rk4<T, V, F>
impl<T, V, F> Unpin for Rk4<T, V, F>
impl<T, V, F> UnwindSafe for Rk4<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.