Module peroxide::numerical::ode

source ·
Expand description

§Ordinary Differential Equation (ODE) Solvers

This module provides traits and structs for solving ordinary differential equations (ODEs).

§Overview

  • ODEProblem: Trait for defining an ODE problem.
  • ODEIntegrator: Trait for ODE integrators.
  • ODESolver: Trait for ODE solvers.
  • ODEError: Enum for ODE errors.
    • ReachedMaxStepIter: Reached maximum number of steps per step. (internal error)
    • ConstraintViolation(f64, Vec<f64>, Vec<f64>): Constraint violation. (user-defined error)
    • ODE uses anyhow for error handling. So, you can customize your errors.

§Available integrators

  • Explicit
    • Ralston’s 3rd order (RALS3)
    • Runge-Kutta 4th order (RK4)
    • Ralston’s 4th order (RALS4)
    • Runge-Kutta 5th order (RK5)
  • Embedded
    • Bogacki-Shampine 2/3rd order (BS23)
    • Runge-Kutta-Fehlberg 4/5th order (RKF45)
    • Dormand-Prince 4/5th order (DP45)
    • Tsitouras 4/5th order (TSIT45)
  • Implicit
    • Gauss-Legendre 4th order (GL4)

§Available solvers

  • BasicODESolver: A basic ODE solver using a specified integrator.

You can implement your own ODE solver by implementing the ODESolver trait.

§Example

use peroxide::fuga::*;

fn main() -> Result<(), Box<dyn Error>> {
    // Same as : let rkf = RKF45::new(1e-4, 0.9, 1e-6, 1e-1, 100);
    let rkf = RKF45 {
        tol: 1e-6,
        safety_factor: 0.9,
        min_step_size: 1e-6,
        max_step_size: 1e-1,
        max_step_iter: 100,
    };
    let basic_ode_solver = BasicODESolver::new(rkf);
    let initial_conditions = vec![1f64];
    let (t_vec, y_vec) = basic_ode_solver.solve(
        &Test,
        (0f64, 10f64),
        0.01,
        &initial_conditions,
    )?;
    let y_vec: Vec<f64> = y_vec.into_iter().flatten().collect();
    println!("{}", y_vec.len());

    let mut plt = Plot2D::new();
    plt
        .set_domain(t_vec)
        .insert_image(y_vec)
        .set_xlabel(r"$t$")
        .set_ylabel(r"$y$")
        .set_style(PlotStyle::Nature)
        .tight_layout()
        .set_dpi(600)
        .set_path("example_data/rkf45_test.png")
        .savefig()?;
    Ok(())
}

// Extremely customizable struct
struct Test;

impl ODEProblem for Test {
    fn rhs(&self, t: f64, y: &[f64], dy: &mut [f64]) -> anyhow::Result<()> {
        Ok(dy[0] = (5f64 * t.powi(2) - y[0]) / (t + y[0]).exp())
    }
}

Structs§

  • Bogacki-Shampine 3(2) method
  • A basic ODE solver using a specified integrator.
  • Dormand-Prince 5(4) method
  • Gauss-Legendre 4th order integrator.
  • Ralston’s 3rd order integrator
  • Ralston’s 4th order integrator.
  • Runge-Kutta 4th order integrator.
  • Runge-Kutta 5th order integrator
  • Runge-Kutta-Fehlberg 4/5th order integrator.
  • Tsitouras 5(4) method

Enums§

Traits§