pub struct Bdf<'a, Eqn: OdeEquationsImplicit, Nls: NonLinearSolver<Eqn::M>, M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C> = <<Eqn as Op>::V as DefaultDenseMatrix>::M, AugmentedEqn: AugmentedOdeEquationsImplicit<Eqn> = NoAug<Eqn>>where
Eqn::V: DefaultDenseMatrix,{ /* private fields */ }
Expand description
Implements a Backward Difference formula (BDF) implicit multistep integrator.
The basic algorithm is derived in [1]. This particular implementation follows that implemented in the Matlab routine ode15s described in [2] and the SciPy implementation /[3/], which features the NDF formulas for improved stability with associated differences in the error constants, and calculates the jacobian at J(t_{n+1}, y^0_{n+1}). This implementation was based on that implemented in the SciPy library [3], which also mainly follows [2] but uses the more standard Jacobian update.
§References
[1] Byrne, G. D., & Hindmarsh, A. C. (1975). A polyalgorithm for the numerical solution of ordinary differential equations. ACM Transactions on Mathematical Software (TOMS), 1(1), 71-96. [2] Shampine, L. F., & Reichelt, M. W. (1997). The matlab ode suite. SIAM journal on scientific computing, 18(1), 1-22. [3] Virtanen, P., Gommers, R., Oliphant, T. E., Haberland, M., Reddy, T., Cournapeau, D., … & Van Mulbregt, P. (2020). SciPy 1.0: fundamental algorithms for scientific computing in Python. Nature methods, 17(3), 261-272.
Implementations§
Source§impl<'a, M, Eqn, Nls, AugmentedEqn> Bdf<'a, Eqn, Nls, M, AugmentedEqn>where
AugmentedEqn: AugmentedOdeEquations<Eqn> + OdeEquationsImplicit,
Eqn: OdeEquationsImplicit,
Eqn::V: DefaultDenseMatrix,
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
for<'b> &'b Eqn::V: VectorRef<Eqn::V>,
for<'b> &'b Eqn::M: MatrixRef<Eqn::M>,
Nls: NonLinearSolver<Eqn::M>,
impl<'a, M, Eqn, Nls, AugmentedEqn> Bdf<'a, Eqn, Nls, M, AugmentedEqn>where
AugmentedEqn: AugmentedOdeEquations<Eqn> + OdeEquationsImplicit,
Eqn: OdeEquationsImplicit,
Eqn::V: DefaultDenseMatrix,
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
for<'b> &'b Eqn::V: VectorRef<Eqn::V>,
for<'b> &'b Eqn::M: MatrixRef<Eqn::M>,
Nls: NonLinearSolver<Eqn::M>,
pub fn new( problem: &'a OdeSolverProblem<Eqn>, state: BdfState<Eqn::V, M>, nonlinear_solver: Nls, ) -> Result<Self, DiffsolError>
pub fn new_augmented( state: BdfState<Eqn::V, M>, problem: &'a OdeSolverProblem<Eqn>, augmented_eqn: AugmentedEqn, nonlinear_solver: Nls, ) -> Result<Self, DiffsolError>
pub fn get_statistics(&self) -> &BdfStatistics
Trait Implementations§
Source§impl<'a, M, Eqn, Nls, AugEqn> AugmentedOdeSolverMethod<'a, Eqn, AugEqn> for Bdf<'a, Eqn, Nls, M, AugEqn>where
Eqn: OdeEquationsImplicit,
AugEqn: AugmentedOdeEquationsImplicit<Eqn>,
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
for<'b> &'b Eqn::V: VectorRef<Eqn::V>,
for<'b> &'b Eqn::M: MatrixRef<Eqn::M>,
Nls: NonLinearSolver<Eqn::M>,
Eqn::V: DefaultDenseMatrix<T = Eqn::T>,
impl<'a, M, Eqn, Nls, AugEqn> AugmentedOdeSolverMethod<'a, Eqn, AugEqn> for Bdf<'a, Eqn, Nls, M, AugEqn>where
Eqn: OdeEquationsImplicit,
AugEqn: AugmentedOdeEquationsImplicit<Eqn>,
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
for<'b> &'b Eqn::V: VectorRef<Eqn::V>,
for<'b> &'b Eqn::M: MatrixRef<Eqn::M>,
Nls: NonLinearSolver<Eqn::M>,
Eqn::V: DefaultDenseMatrix<T = Eqn::T>,
fn into_state_and_eqn(self) -> (Self::State, Option<AugEqn>)
fn augmented_eqn(&self) -> Option<&AugEqn>
Source§impl<M, Eqn, Nls, AugmentedEqn> Clone for Bdf<'_, Eqn, Nls, M, AugmentedEqn>where
Eqn: OdeEquationsImplicit,
Nls: NonLinearSolver<Eqn::M>,
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
AugmentedEqn: AugmentedOdeEquationsImplicit<Eqn>,
Eqn::V: DefaultDenseMatrix,
impl<M, Eqn, Nls, AugmentedEqn> Clone for Bdf<'_, Eqn, Nls, M, AugmentedEqn>where
Eqn: OdeEquationsImplicit,
Nls: NonLinearSolver<Eqn::M>,
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
AugmentedEqn: AugmentedOdeEquationsImplicit<Eqn>,
Eqn::V: DefaultDenseMatrix,
Source§impl<'a, M, Eqn, Nls, AugmentedEqn> OdeSolverMethod<'a, Eqn> for Bdf<'a, Eqn, Nls, M, AugmentedEqn>where
Eqn: OdeEquationsImplicit,
AugmentedEqn: AugmentedOdeEquations<Eqn> + OdeEquationsImplicit,
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
Eqn::V: DefaultDenseMatrix,
Nls: NonLinearSolver<Eqn::M>,
for<'b> &'b Eqn::V: VectorRef<Eqn::V>,
for<'b> &'b Eqn::M: MatrixRef<Eqn::M>,
impl<'a, M, Eqn, Nls, AugmentedEqn> OdeSolverMethod<'a, Eqn> for Bdf<'a, Eqn, Nls, M, AugmentedEqn>where
Eqn: OdeEquationsImplicit,
AugmentedEqn: AugmentedOdeEquations<Eqn> + OdeEquationsImplicit,
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
Eqn::V: DefaultDenseMatrix,
Nls: NonLinearSolver<Eqn::M>,
for<'b> &'b Eqn::V: VectorRef<Eqn::V>,
for<'b> &'b Eqn::M: MatrixRef<Eqn::M>,
type State = BdfState<<Eqn as Op>::V, M>
Source§fn order(&self) -> usize
fn order(&self) -> usize
Source§fn jacobian(&self) -> Option<Ref<'_, Eqn::M>>
fn jacobian(&self) -> Option<Ref<'_, Eqn::M>>
Source§fn mass(&self) -> Option<Ref<'_, Eqn::M>>
fn mass(&self) -> Option<Ref<'_, Eqn::M>>
Source§fn set_state(&mut self, state: Self::State)
fn set_state(&mut self, state: Self::State)
Source§fn interpolate(&self, t: Eqn::T) -> Result<Eqn::V, DiffsolError>
fn interpolate(&self, t: Eqn::T) -> Result<Eqn::V, DiffsolError>
Source§fn interpolate_out(&self, t: Eqn::T) -> Result<Eqn::V, DiffsolError>
fn interpolate_out(&self, t: Eqn::T) -> Result<Eqn::V, DiffsolError>
Source§fn interpolate_sens(
&self,
t: <Eqn as Op>::T,
) -> Result<Vec<Eqn::V>, DiffsolError>
fn interpolate_sens( &self, t: <Eqn as Op>::T, ) -> Result<Vec<Eqn::V>, DiffsolError>
Source§fn problem(&self) -> &'a OdeSolverProblem<Eqn>
fn problem(&self) -> &'a OdeSolverProblem<Eqn>
Source§fn into_state(self) -> BdfState<Eqn::V, M>
fn into_state(self) -> BdfState<Eqn::V, M>
set_problem
again before calling step
or solve
.Source§fn state_mut(&mut self) -> StateRefMut<'_, Eqn::V>
fn state_mut(&mut self) -> StateRefMut<'_, Eqn::V>
step
to perform some reinitialisation to take into
account the mutated state, this could be expensive for multi-step methods.Source§fn checkpoint(&mut self) -> Self::State
fn checkpoint(&mut self) -> Self::State
take_state
instead.
Note that this will force a reinitialisation of the internal Jacobian for the solver, if it has one.Source§fn step(&mut self) -> Result<OdeSolverStopReason<Eqn::T>, DiffsolError>
fn step(&mut self) -> Result<OdeSolverStopReason<Eqn::T>, DiffsolError>
Result
containing the reason for stopping the solver, possible reasons are: Read moreSource§fn set_stop_time(&mut self, tstop: <Eqn as Op>::T) -> Result<(), DiffsolError>
fn set_stop_time(&mut self, tstop: <Eqn as Op>::T) -> Result<(), DiffsolError>
tstop
is at or before the current internal time, an error is returned.Source§fn solve(
&mut self,
final_time: Eqn::T,
) -> Result<(<Eqn::V as DefaultDenseMatrix>::M, Vec<Eqn::T>), DiffsolError>
fn solve( &mut self, final_time: Eqn::T, ) -> Result<(<Eqn::V as DefaultDenseMatrix>::M, Vec<Eqn::T>), DiffsolError>
final_time
Returns a Vec of solution values at timepoints chosen by the solver.
After the solver has finished, the internal state of the solver is at time final_time
.Source§fn solve_dense(
&mut self,
t_eval: &[Eqn::T],
) -> Result<<Eqn::V as DefaultDenseMatrix>::M, DiffsolError>
fn solve_dense( &mut self, t_eval: &[Eqn::T], ) -> Result<<Eqn::V as DefaultDenseMatrix>::M, DiffsolError>
t_eval[t_eval.len()-1]
Returns a Vec of solution values at timepoints given by t_eval
.
After the solver has finished, the internal state of the solver is at time t_eval[t_eval.len()-1]
.fn solve_with_checkpointing( &mut self, final_time: Eqn::T, max_steps_between_checkpoints: Option<usize>, ) -> Result<(Checkpointing<'a, Eqn, Self>, <Eqn::V as DefaultDenseMatrix>::M, Vec<Eqn::T>), DiffsolError>
Source§fn solve_dense_with_checkpointing(
&mut self,
t_eval: &[Eqn::T],
max_steps_between_checkpoints: Option<usize>,
) -> Result<(Checkpointing<'a, Eqn, Self>, <Eqn::V as DefaultDenseMatrix>::M), DiffsolError>
fn solve_dense_with_checkpointing( &mut self, t_eval: &[Eqn::T], max_steps_between_checkpoints: Option<usize>, ) -> Result<(Checkpointing<'a, Eqn, Self>, <Eqn::V as DefaultDenseMatrix>::M), DiffsolError>
Auto Trait Implementations§
impl<'a, Eqn, Nls, M = <<Eqn as Op>::V as DefaultDenseMatrix>::M, AugmentedEqn = NoAug<Eqn>> !Freeze for Bdf<'a, Eqn, Nls, M, AugmentedEqn>
impl<'a, Eqn, Nls, M = <<Eqn as Op>::V as DefaultDenseMatrix>::M, AugmentedEqn = NoAug<Eqn>> !RefUnwindSafe for Bdf<'a, Eqn, Nls, M, AugmentedEqn>
impl<'a, Eqn, Nls, M, AugmentedEqn> Send for Bdf<'a, Eqn, Nls, M, AugmentedEqn>
impl<'a, Eqn, Nls, M = <<Eqn as Op>::V as DefaultDenseMatrix>::M, AugmentedEqn = NoAug<Eqn>> !Sync for Bdf<'a, Eqn, Nls, M, AugmentedEqn>
impl<'a, Eqn, Nls, M, AugmentedEqn> Unpin for Bdf<'a, Eqn, Nls, M, AugmentedEqn>
impl<'a, Eqn, Nls, M, AugmentedEqn> UnwindSafe for Bdf<'a, Eqn, Nls, M, AugmentedEqn>where
<Eqn as Op>::V: Sized + UnwindSafe + RefUnwindSafe,
Nls: UnwindSafe,
<Eqn as Op>::T: UnwindSafe + RefUnwindSafe,
M: UnwindSafe,
Eqn: RefUnwindSafe,
AugmentedEqn: UnwindSafe,
<<Eqn as Op>::M as Matrix>::Sparsity: UnwindSafe,
<Eqn as Op>::M: UnwindSafe,
Blanket Implementations§
Source§impl<'a, Eqn, S, Solver> AdjointOdeSolverMethod<'a, Eqn, S> for Solverwhere
Eqn: OdeEquationsImplicitAdjoint + 'a,
S: OdeSolverMethod<'a, Eqn>,
Solver: AugmentedOdeSolverMethod<'a, Eqn, AdjointEquations<'a, Eqn, S>>,
impl<'a, Eqn, S, Solver> AdjointOdeSolverMethod<'a, Eqn, S> for Solverwhere
Eqn: OdeEquationsImplicitAdjoint + 'a,
S: OdeSolverMethod<'a, Eqn>,
Solver: AugmentedOdeSolverMethod<'a, Eqn, AdjointEquations<'a, Eqn, S>>,
Source§fn solve_adjoint_backwards_pass(
self,
t_eval: &[Eqn::T],
dgdu_eval: &[&<Eqn::V as DefaultDenseMatrix>::M],
) -> Result<Self::State, DiffsolError>
fn solve_adjoint_backwards_pass( self, t_eval: &[Eqn::T], dgdu_eval: &[&<Eqn::V as DefaultDenseMatrix>::M], ) -> Result<Self::State, DiffsolError>
Source§impl<T> BorrowMut<T> for Twhere
T: ?Sized,
impl<T> BorrowMut<T> for Twhere
T: ?Sized,
Source§fn borrow_mut(&mut self) -> &mut T
fn borrow_mut(&mut self) -> &mut T
Source§impl<T> CloneToUninit for Twhere
T: Clone,
impl<T> CloneToUninit for Twhere
T: Clone,
Source§impl<T> DistributionExt for Twhere
T: ?Sized,
impl<T> DistributionExt for Twhere
T: ?Sized,
Source§impl<T> IntoEither for T
impl<T> IntoEither for T
Source§fn into_either(self, into_left: bool) -> Either<Self, Self>
fn into_either(self, into_left: bool) -> Either<Self, Self>
self
into a Left
variant of Either<Self, Self>
if into_left
is true
.
Converts self
into a Right
variant of Either<Self, Self>
otherwise. Read moreSource§fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
fn into_either_with<F>(self, into_left: F) -> Either<Self, Self>
self
into a Left
variant of Either<Self, Self>
if into_left(&self)
returns true
.
Converts self
into a Right
variant of Either<Self, Self>
otherwise. Read moreSource§impl<T> Pointable for T
impl<T> Pointable for T
Source§impl<'a, M, Eqn> SensitivitiesOdeSolverMethod<'a, Eqn> for Mwhere
M: AugmentedOdeSolverMethod<'a, Eqn, SensEquations<'a, Eqn>>,
Eqn: OdeEquationsImplicitSens + 'a,
impl<'a, M, Eqn> SensitivitiesOdeSolverMethod<'a, Eqn> for Mwhere
M: AugmentedOdeSolverMethod<'a, Eqn, SensEquations<'a, Eqn>>,
Eqn: OdeEquationsImplicitSens + 'a,
Source§fn solve_dense_sensitivities(
&mut self,
t_eval: &[Eqn::T],
) -> Result<(<Eqn::V as DefaultDenseMatrix>::M, Vec<<Eqn::V as DefaultDenseMatrix>::M>), DiffsolError>where
Eqn: OdeEquationsImplicitSens,
Eqn::M: DefaultSolver,
Eqn::V: DefaultDenseMatrix,
Self: Sized,
fn solve_dense_sensitivities(
&mut self,
t_eval: &[Eqn::T],
) -> Result<(<Eqn::V as DefaultDenseMatrix>::M, Vec<<Eqn::V as DefaultDenseMatrix>::M>), DiffsolError>where
Eqn: OdeEquationsImplicitSens,
Eqn::M: DefaultSolver,
Eqn::V: DefaultDenseMatrix,
Self: Sized,
t_eval[t_eval.len()-1]
Returns a tuple (y, sens)
, where y
is a dense matrix of solution values at timepoints given by t_eval
,
and sens
is a Vec of dense matrices, the ith element of the Vec are the the sensitivities with respect to the ith parameter.
After the solver has finished, the internal state of the solver is at time t_eval[t_eval.len()-1]
.Source§impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
impl<SS, SP> SupersetOf<SS> for SPwhere
SS: SubsetOf<SP>,
Source§fn to_subset(&self) -> Option<SS>
fn to_subset(&self) -> Option<SS>
self
from the equivalent element of its
superset. Read moreSource§fn is_in_subset(&self) -> bool
fn is_in_subset(&self) -> bool
self
is actually part of its subset T
(and can be converted to it).Source§fn to_subset_unchecked(&self) -> SS
fn to_subset_unchecked(&self) -> SS
self.to_subset
but without any property checks. Always succeeds.Source§fn from_subset(element: &SS) -> SP
fn from_subset(element: &SS) -> SP
self
to the equivalent element of its superset.