[][src]Trait mathru::analysis::ode::Solver

pub trait Solver<T> where
    T: Field
{ fn solve<F>(
        &self,
        function: F,
        initial_cond: Vector<T>,
        t_start: T,
        t_end: T
    ) -> (Vector<T>, Matrix<T>)
    where
        F: Fn(&T, &Vector<T>) -> Vector<T>
; }

Trait which is implemented by every ode algorithm

Required methods

fn solve<F>(
    &self,
    function: F,
    initial_cond: Vector<T>,
    t_start: T,
    t_end: T
) -> (Vector<T>, Matrix<T>) where
    F: Fn(&T, &Vector<T>) -> Vector<T>, 

Solve the ordinary diffential equation, returning the results of the calculation.

Loading content...

Implementors

impl<T> Solver<T> for Dopri5<T> where
    T: Real
[src]

fn solve<F>(
    &self,
    func: F,
    init: Vector<T>,
    t_start: T,
    t_end: T
) -> (Vector<T>, Matrix<T>) where
    F: Fn(&T, &Vector<T>) -> Vector<T>, 
[src]

Solves func using the 4th order Runge-Kutta-Dormand-Prince algorithm.

Arguments

  • 'func' is an explict oridnary diffential equation
  • 'init' is the initial value at the time 't_start'
  • 't_start' initial time
  • 't_end'

Return

The solver returns a vector and a matrix, containing the times used in each step of the algorithm and the respectful values for that time.

Example

use mathru::*;
use mathru::algebra::linear::{Vector, Matrix};
use mathru::analysis::ode::{Solver, Dopri5};

let h_0: f64 = 0.0001;
let e_max: f64 = 0.000001;
let n_max: u32 = 500;

fn f(_t: &f64, x: &Vector<f64>) -> Vector<f64>
{
    let result  = vector![1.0] + x.clone().apply(&|e: &f64| -> f64 {return e * e;}) ;

    return result;
}

let init: Vector<f64> = vector![0.0];
let solver: Dopri5<f64> = Dopri5::new(h_0, e_max, n_max);

let (t, y): (Vector<f64>, Matrix<f64>) = solver.solve(f, init, 0.0, 1.4);


let (m, _n): (usize, usize) = y.dim();

assert!((1.40_f64 - *t.get(m - 1)).abs() < 0.0001);
assert!((1.4_f64.tan() - *y.get(m - 1, 0)).abs() < 0.0001);

impl<T> Solver<T> for Euler<T> where
    T: Real
[src]

fn solve<F>(
    &self,
    func: F,
    init: Vector<T>,
    t_start: T,
    t_end: T
) -> (Vector<T>, Matrix<T>) where
    F: Fn(&T, &Vector<T>) -> Vector<T>, 
[src]

Solves func using Euler's method.

Arguments

  • 'func' is an explict oridnary diffential equation
  • 'init' is the initial value at the time 't_start'
  • 't_start' initial time
  • 't_end'

Return

The solver returns a vector and a matrix, containing the times used in each step of the algorithm and the respectful values for that time.

Example

use mathru::*;
use mathru::algebra::linear::{vector, Vector, Matrix};
use mathru::analysis::ode::{Solver, Euler};

let f = |t: &f64, _x: &Vector<f64> | -> Vector<f64> { return Vector::new_row(1, vec![1.0]) * (t * &2.0f64); };

let init: Vector<f64> = vector![1.0];
let solver: Euler<f64> = Euler::new(0.01);

let (t, y): (Vector<f64>, Matrix<f64>) = solver.solve(f, init, 0.0, 2.0);

impl<T> Solver<T> for Heun<T> where
    T: Real
[src]

fn solve<F>(
    &self,
    func: F,
    init: Vector<T>,
    t_start: T,
    t_end: T
) -> (Vector<T>, Matrix<T>) where
    F: Fn(&T, &Vector<T>) -> Vector<T>, 
[src]

Solves func using Heun's method.

Arguments

  • 'func' is an explict oridnary diffential equation
  • 'init' is the initial value at the time 't_start'
  • 't_start' initial time
  • 't_end'

Return

The solver returns a vector and a matrix, containing the times used in each step of the algorithm and the respectful values for that time.

Example

use mathru::*;
use mathru::algebra::linear::{Vector, Matrix};
use mathru::analysis::ode::{Solver, Heun};

let f = |t: &f64, _x: &Vector<f64> | -> Vector<f64> { return Vector::new_row(1, vec![1.0]) * (t * &2.0f64); };

let init: Vector<f64> = vector![1.0];
let solver: Heun<f64> = Heun::new(0.01);

let (t, y): (Vector<f64>, Matrix<f64>) = solver.solve(f, init, 0.0, 2.0);

impl<T> Solver<T> for RK4<T> where
    T: Real
[src]

fn solve<F>(
    &self,
    func: F,
    init: Vector<T>,
    t_start: T,
    t_end: T
) -> (Vector<T>, Matrix<T>) where
    F: Fn(&T, &Vector<T>) -> Vector<T>, 
[src]

Solves func using the 4th order Runge-Kutta algorithm.

Arguments

  • 'func' is an explict oridnary diffential equation
  • 'init' is the initial value at the time 't_start'
  • 't_start' initial time
  • 't_end'

Return

The solver returns a vector and a matrix, containing the times used in each step of the algorithm and the respectful values for that time.

Example

use mathru::*;
use mathru::algebra::linear::{Vector, Matrix};
use mathru::analysis::ode::{Solver, RK4};

let f = |t: &f64, _x: &Vector<f64> | -> Vector<f64> { return Vector::new_row(1, vec![1.0]) * (t * &2.0f64); };

let init: Vector<f64> = vector![1.0];
let solver: RK4<f64> = RK4::new(0.01);

let (t, y): (Vector<f64>, Matrix<f64>) = solver.solve(f, init, 0.0, 2.0);

impl<T> Solver<T> for RKF45<T> where
    T: Real
[src]

fn solve<F>(
    &self,
    func: F,
    init: Vector<T>,
    t_start: T,
    t_end: T
) -> (Vector<T>, Matrix<T>) where
    F: Fn(&T, &Vector<T>) -> Vector<T>, 
[src]

Solves func using the 4th order Runge-Kutta-Fehlberg algorithm.

Arguments

  • 'func' is an explict oridnary diffential equation
  • 'init' is the initial value at the time 't_start'
  • 't_start' initial time
  • 't_end'

Return

The solver returns a vector and a matrix, containing the times used in each step of the algorithm and the respectful values for that time.

Example

use mathru::*;
use mathru::algebra::linear::{Vector, Matrix};
use mathru::analysis::ode::{Solver, RKF45};

let f = |t: &f64, _x: &Vector<f64> | -> Vector<f64> { return Vector::new_row(1, vec![1.0]) * (t * &2.0f64); };

let init: Vector<f64> = vector![1.0];
let solver: RKF45<f64> = RKF45::new(0.001, 0.001, 0.01, 0.1, 3000);

let (t, y): (Vector<f64>, Matrix<f64>) = solver.solve(f, init, 0.0, 2.0);
Loading content...