[−][src]Trait bacon_sci::ivp::RungeKuttaSolver
This trait allows a struct to be used in the Runge-Kutta solver.
Things implementing RungeKuttaSolver should have an RKInfo to handle the actual IVP solving. It should also provide the with_* helper functions for convience.
Examples
See struct RK45
for an example of implementing this trait
Required methods
pub fn t_coefficients() -> VectorN<N::RealField, O>
[src]
Returns a vec of coeffecients to multiply the time step by when getting intermediate results. Upper-left portion of Butch Tableaux
pub fn k_coefficients() -> MatrixMN<N::RealField, O, O>
[src]
Returns the coefficients to use on the k_i's when finding another k_i. Upper-right portion of the Butch Tableax. Should be an NxN-1 matrix, where N is the order of the Runge-Kutta Method (Or order+1 for adaptive methods)
pub fn avg_coefficients() -> VectorN<N::RealField, O>
[src]
Coefficients to use when calculating the final step to take. These are the weights of the weighted average of k_i's. Bottom portion of the Butch Tableaux. For adaptive methods, this is the first row of the bottom portion.
pub fn error_coefficients() -> VectorN<N::RealField, O>
[src]
Coefficients to use on the k_i's to find the error between the two orders of Runge-Kutta methods. In the Butch Tableaux, this is the first row of the bottom portion minus the second row.
pub fn solve_ivp<T: Clone, F: FnMut(N::RealField, &[N], &mut T) -> Result<VectorN<N, S>, String>>(
self,
f: F,
params: &mut T
) -> Result<Vec<(N::RealField, VectorN<N, S>)>, String>
[src]
self,
f: F,
params: &mut T
) -> Result<Vec<(N::RealField, VectorN<N, S>)>, String>
Ideally, call RKInfo.solve_ivp
pub fn with_tolerance(self, tol: N::RealField) -> Result<Self, String>
[src]
Set the error tolerance for this solver.
pub fn with_dt_max(self, max: N::RealField) -> Result<Self, String>
[src]
Set the maximum time step for this solver.
pub fn with_dt_min(self, min: N::RealField) -> Result<Self, String>
[src]
Set the minimum time step for this solver.
pub fn with_start(self, t_initial: N::RealField) -> Result<Self, String>
[src]
Set the initial time for this solver.
pub fn with_end(self, t_final: N::RealField) -> Result<Self, String>
[src]
Set the end time for this solver.
pub fn with_initial_conditions(self, start: &[N]) -> Result<Self, String>
[src]
Set the initial conditions for this solver.
pub fn build(self) -> Self
[src]
Build this solver.
Implementors
impl<N: ComplexField, S: DimName> RungeKuttaSolver<N, S, U6> for RK45<N, S> where
DefaultAllocator: Allocator<N, S>,
DefaultAllocator: Allocator<N, U6>,
DefaultAllocator: Allocator<N, U6, U6>,
DefaultAllocator: Allocator<N::RealField, U6>,
DefaultAllocator: Allocator<N::RealField, U6, U6>,
DefaultAllocator: Allocator<N, S, U6>,
[src]
DefaultAllocator: Allocator<N, S>,
DefaultAllocator: Allocator<N, U6>,
DefaultAllocator: Allocator<N, U6, U6>,
DefaultAllocator: Allocator<N::RealField, U6>,
DefaultAllocator: Allocator<N::RealField, U6, U6>,
DefaultAllocator: Allocator<N, S, U6>,
pub fn t_coefficients() -> VectorN<N::RealField, U6>
[src]
pub fn k_coefficients() -> MatrixMN<N::RealField, U6, U6>
[src]
pub fn avg_coefficients() -> VectorN<N::RealField, U6>
[src]
pub fn error_coefficients() -> VectorN<N::RealField, U6>
[src]
pub fn solve_ivp<T: Clone, F: FnMut(N::RealField, &[N], &mut T) -> Result<VectorN<N, S>, String>>(
self,
f: F,
params: &mut T
) -> Result<Vec<(N::RealField, VectorN<N, S>)>, String>
[src]
self,
f: F,
params: &mut T
) -> Result<Vec<(N::RealField, VectorN<N, S>)>, String>