Skip to main content

Crate ivp

Crate ivp 

Source
Expand description

ivp: Initial value problem solvers for ODEs.

This crate provides explicit Runge–Kutta methods (RK4, RK23, DOPRI5, DOP853), LSODA automatic-switching multistep methods, implicit methods (Radau, BDF), and structured fixed-step symplectic solvers for separable Hamiltonian and second-order systems.

Highlights

  • Methods: RK4, RK23, DOPRI5, DOP853, LSODA, Radau, BDF
  • Controls: rtol, atol, first_step, min_step, max_step, nmax
  • Sampling: internal accepted steps by default, or exact t_eval times
  • Dense output: sol(t), sol_many(&ts), sol_span() on the returned Solution
  • Iteration: iterate stored samples via solution.iter()
  • Typed Rust API: Ivp::first_order, Ivp::second_order, and Ivp::hamiltonian

Quick start

use ivp::prelude::*;
use std::f64::consts::PI;

struct SHO;
impl FirstOrderSystem for SHO {
    fn derivative(&self, _x: f64, y: &[f64], dydx: &mut [f64]) {
        dydx[0] = y[1];
        dydx[1] = -y[0];
    }
}

let f = SHO;
let x0 = 0.0;
let xend = 2.0 * PI; // one period
let y0 = [1.0, 0.0];

let sol = Ivp::first_order(&f, x0, xend, &y0)
    .method(Method::DOP853)
    .rtol(1e-9).atol(1e-9)
    .dense_output(true)
    .solve()
    .unwrap();

// Discrete samples
println!("Discrete output at accepted steps:");
for (t, y) in sol.iter() {
    println!("x = {:>8.5}, y = {:?}", t, y);
}

// Continuous evaluation within the solution span
if let Some((t0, t1)) = sol.sol_span() {
    let ts = [t0, 0.5*(t0+t1), t1];
    let ys = sol.sol_many(&ts).unwrap();
    println!("\nDense output at t0, (t0+t1)/2, t1:");
    for (t, y) in ts.iter().zip(ys.iter()) {
        println!("x = {:>8.5}, y = {:?}", t, y);
    }
}

Structured symplectic example

use ivp::prelude::*;

struct HarmonicOscillator;
impl SecondOrderSystem for HarmonicOscillator {
    fn acceleration(&self, _t: f64, q: &[f64], a: &mut [f64]) {
        a[0] = -q[0];
    }
}

let sol = Ivp::second_order(&HarmonicOscillator, 0.0, 1.0, &[1.0], &[0.0])
    .method(SymplecticMethod::VelocityVerlet)
    .step_size(0.05)
    .solve()
    .unwrap();
assert_eq!(sol.y.last().unwrap().len(), 2);

§License

Copyright 2025 Ryan D. Gast

Licensed under the Apache License, Version 2.0 (the "License");
you may not use this file except in compliance with the License.
You may obtain a copy of the License at

    http://www.apache.org/licenses/LICENSE-2.0

Unless required by applicable law or agreed to in writing, software
distributed under the License is distributed on an "AS IS" BASIS,
WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
See the License for the specific language governing permissions and
limitations under the License.

Re-exports§

pub use methods::SymplecticMethod;
pub use solve::FirstOrderIvp;
pub use solve::HamiltonianIvp;
pub use solve::Ivp;
pub use solve::JacobianSource;
pub use solve::Method;
pub use solve::SecondOrderIvp;
pub use solve::Solution;

Modules§

dense
Unified dense output interpolation types.
error
Errors for integration methods
ivp
User-supplied dynamical system definitions.
matrix
Matrix types, operations, and utilities.
methods
Numerical method implementations and shared solver utilities.
prelude
Convenient prelude: import the most commonly used traits, types, and functions.
solout
User defined callback hook executed after each accepted step.
solve
High-level solve module: SciPy-like API pieces split into submodules.
status
Status codes for integrators

Macros§

banded_matrix
Create a banded matrix by specifying diagonals. Size and bands are inferred. Usage: banded_matrix!( 0 => [d0…], 1 => [d1…], -1 => [u1…], k => [..], … )
matrix
Create a full dense matrix from rows. Usage: