pub struct EM<T: Real, V: State<T>, D: CallBackData> {
pub h: T,
/* private fields */
}Expand description
Euler-Maruyama Method for solving stochastic differential equations.
The Euler-Maruyama method is the simplest numerical method for SDEs, essentially extending Euler’s method to stochastic equations.
For an SDE of the form: dY = a(t,Y)dt + b(t,Y)dW
The Euler-Maruyama update is: Y_{n+1} = Y_n + a(t_n, Y_n)Δt + b(t_n, Y_n)ΔW_n
where ΔW_n is a Wiener process increment, produced by the SDE’s noise method.
§Example
use differential_equations::prelude::*;
use nalgebra::SVector;
use rand::SeedableRng;
use rand_distr::{Distribution, Normal};
struct GBM {
rng: rand::rngs::StdRng,
}
impl GBM {
fn new(seed: u64) -> Self {
Self {
rng: rand::rngs::StdRng::seed_from_u64(seed),
}
}
}
impl SDE<f64, SVector<f64, 1>> for GBM {
fn drift(&self, _t: f64, y: &SVector<f64, 1>, dydt: &mut SVector<f64, 1>) {
dydt[0] = 0.1 * y[0]; // μS
}
fn diffusion(&self, _t: f64, y: &SVector<f64, 1>, dydw: &mut SVector<f64, 1>) {
dydw[0] = 0.2 * y[0]; // σS
}
fn noise(&self, dt: f64, dw: &mut SVector<f64, 1>) {
let normal = Normal::new(0.0, dt.sqrt()).unwrap();
dw[0] = normal.sample(&mut self.rng.clone());
}
}
let t0 = 0.0;
let tf = 1.0;
let y0 = SVector::<f64, 1>::new(100.0);
let mut solver = EM::new(0.01);
let gbm = GBM::new(42);
let gbm_problem = SDEProblem::new(gbm, t0, tf, y0);
// Solve the SDE
let result = gbm_problem.solve(&mut solver);Fields§
§h: TImplementations§
Trait Implementations§
Source§impl<T: Real, V: State<T>, D: CallBackData> Interpolation<T, V> for EM<T, V, D>
impl<T: Real, V: State<T>, D: CallBackData> Interpolation<T, V> for EM<T, V, D>
Source§impl<T: Real, V: State<T>, D: CallBackData> SDENumericalMethod<T, V, D> for EM<T, V, D>
impl<T: Real, V: State<T>, D: CallBackData> SDENumericalMethod<T, V, D> for EM<T, V, D>
Source§fn init<F>(
&mut self,
sde: &F,
t0: T,
tf: T,
y0: &V,
) -> Result<Evals, Error<T, V>>where
F: SDE<T, V, D>,
fn init<F>(
&mut self,
sde: &F,
t0: T,
tf: T,
y0: &V,
) -> Result<Evals, Error<T, V>>where
F: SDE<T, V, D>,
Initialize SDENumericalMethod before solving SDE Read more
Source§fn step<F>(&mut self, sde: &F) -> Result<Evals, Error<T, V>>where
F: SDE<T, V, D>,
fn step<F>(&mut self, sde: &F) -> Result<Evals, Error<T, V>>where
F: SDE<T, V, D>,
Step through solving the SDE by one step Read more
Source§fn set_status(&mut self, status: Status<T, V, D>)
fn set_status(&mut self, status: Status<T, V, D>)
Set status of solver
Auto Trait Implementations§
impl<T, V, D> Freeze for EM<T, V, D>
impl<T, V, D> RefUnwindSafe for EM<T, V, D>
impl<T, V, D> Send for EM<T, V, D>
impl<T, V, D> Sync for EM<T, V, D>
impl<T, V, D> Unpin for EM<T, V, D>
impl<T, V, D> UnwindSafe for EM<T, V, D>
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.