use crate::StrError;
use crate::{OdeSolverTrait, Params, System, Workspace};
use russell_lab::{vec_add, vec_copy, Vector};
use russell_sparse::CooMatrix;
pub(crate) struct EulerForward<'a, F, J, A>
where
F: Fn(&mut Vector, f64, &Vector, &mut A) -> Result<(), StrError>,
J: Fn(&mut CooMatrix, f64, f64, &Vector, &mut A) -> Result<(), StrError>,
{
system: &'a System<F, J, A>,
k: Vector,
w: Vector,
}
impl<'a, F, J, A> EulerForward<'a, F, J, A>
where
F: Fn(&mut Vector, f64, &Vector, &mut A) -> Result<(), StrError>,
J: Fn(&mut CooMatrix, f64, f64, &Vector, &mut A) -> Result<(), StrError>,
{
pub fn new(system: &'a System<F, J, A>) -> Self {
let ndim = system.ndim;
EulerForward {
system,
k: Vector::new(ndim),
w: Vector::new(ndim),
}
}
}
impl<'a, F, J, A> OdeSolverTrait<A> for EulerForward<'a, F, J, A>
where
F: Fn(&mut Vector, f64, &Vector, &mut A) -> Result<(), StrError>,
J: Fn(&mut CooMatrix, f64, f64, &Vector, &mut A) -> Result<(), StrError>,
{
fn enable_dense_output(&mut self) -> Result<(), StrError> {
Err("dense output is not available for the FwEuler method")
}
fn step(&mut self, work: &mut Workspace, x: f64, y: &Vector, h: f64, args: &mut A) -> Result<(), StrError> {
work.stats.n_function += 1;
(self.system.function)(&mut self.k, x, y, args)?; vec_add(&mut self.w, 1.0, &y, h, &self.k).unwrap(); Ok(())
}
fn accept(
&mut self,
_work: &mut Workspace,
x: &mut f64,
y: &mut Vector,
h: f64,
_args: &mut A,
) -> Result<(), StrError> {
*x += h;
vec_copy(y, &self.w).unwrap();
Ok(())
}
fn reject(&mut self, _work: &mut Workspace, _h: f64) {}
fn dense_output(&self, _y_out: &mut Vector, _x_out: f64, _x: f64, _y: &Vector, _h: f64) {}
fn update_params(&mut self, _params: Params) {}
}
#[cfg(test)]
mod tests {
use super::EulerForward;
use crate::{no_jacobian, HasJacobian, Method, NoArgs, OdeSolverTrait, Params, Samples, System, Workspace};
use russell_lab::{array_approx_eq, Vector};
#[test]
fn euler_forward_works() {
let (system, x0, y0, mut args, y_fn_x) = Samples::kreyszig_eq6_page902();
let ndim = system.ndim;
let mut solver = EulerForward::new(&system);
let mut work = Workspace::new(Method::FwEuler);
assert_eq!(
solver.enable_dense_output().err(),
Some("dense output is not available for the FwEuler method")
);
let h = 0.2;
let mut x = x0;
let mut y = y0.clone();
let mut y_ana = Vector::new(ndim);
y_fn_x(&mut y_ana, x, &mut args);
let mut xx = vec![x];
let mut yy_num = vec![y[0]];
let mut yy_ana = vec![y_ana[0]];
let mut errors = vec![f64::abs(yy_num[0] - yy_ana[0])];
for n in 0..5 {
solver.step(&mut work, x, &y, h, &mut args).unwrap();
assert_eq!(work.stats.n_function, n + 1);
work.stats.n_accepted += 1; solver.accept(&mut work, &mut x, &mut y, h, &mut args).unwrap();
xx.push(x);
yy_num.push(y[0]);
y_fn_x(&mut y_ana, x, &mut args);
yy_ana.push(y_ana[0]);
errors.push(f64::abs(yy_num.last().unwrap() - yy_ana.last().unwrap()));
}
let xx_correct = &[0.0, 0.2, 0.4, 0.6, 0.8, 1.0];
let yy_correct = &[0.0, 0.0, 0.04000000000000001, 0.128, 0.2736000000000001, 0.48832];
let errors_correct = &[
0.0,
0.0214027581601699,
0.05182469764127042,
0.094118800390509,
0.1519409284924678,
0.229961828459045,
];
array_approx_eq(&xx, xx_correct, 1e-15);
array_approx_eq(&yy_num, yy_correct, 1e-15);
array_approx_eq(&errors, errors_correct, 1e-15);
}
#[test]
fn euler_forward_handles_errors() {
let system = System::new(
1,
|f, _, _, _: &mut NoArgs| {
f[0] = 1.0;
Err("stop")
},
no_jacobian,
HasJacobian::No,
None,
None,
);
let mut solver = EulerForward::new(&system);
let mut work = Workspace::new(Method::FwEuler);
let x = 0.0;
let y = Vector::from(&[0.0]);
let h = 0.1;
let mut args = 0;
assert_eq!(solver.step(&mut work, x, &y, h, &mut args).err(), Some("stop"));
let mut y_out = Vector::new(1);
let x_out = 0.1;
solver.reject(&mut work, h);
solver.dense_output(&mut y_out, x_out, x, &y, h);
let params = Params::new(Method::FwEuler);
solver.update_params(params);
}
}