pub struct Sdirk<'a, Eqn, LS, M = <<Eqn as Op>::V as DefaultDenseMatrix>::M, AugmentedEqn = NoAug<Eqn>>where
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
LS: LinearSolver<Eqn::M>,
Eqn: OdeEquationsImplicit,
Eqn::V: DefaultDenseMatrix<T = Eqn::T, C = Eqn::C>,
AugmentedEqn: AugmentedOdeEquations<Eqn>,{ /* private fields */ }
Expand description
A singly diagonally implicit Runge-Kutta method. Can optionally have an explicit first stage for ESDIRK methods.
The particular method is defined by the Tableau used to create the solver.
If the beta
matrix of the Tableau is present this is used for interpolation, otherwise hermite interpolation is used.
Restrictions:
- The upper triangular part of the
a
matrix must be zero (i.e. not fully implicit). - The diagonal of the
a
matrix must be the same non-zero value for all rows (i.e. an SDIRK method), except for the first row which can be zero for ESDIRK methods. - The last row of the
a
matrix must be the same as theb
vector, and the last element of thec
vector must be 1 (i.e. a stiffly accurate method)
Implementations§
Source§impl<'a, M, Eqn, LS, AugmentedEqn> Sdirk<'a, Eqn, LS, M, AugmentedEqn>where
LS: LinearSolver<Eqn::M>,
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
Eqn: OdeEquationsImplicit,
Eqn::V: DefaultDenseMatrix<T = Eqn::T, C = Eqn::C>,
AugmentedEqn: AugmentedOdeEquationsImplicit<Eqn>,
impl<'a, M, Eqn, LS, AugmentedEqn> Sdirk<'a, Eqn, LS, M, AugmentedEqn>where
LS: LinearSolver<Eqn::M>,
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
Eqn: OdeEquationsImplicit,
Eqn::V: DefaultDenseMatrix<T = Eqn::T, C = Eqn::C>,
AugmentedEqn: AugmentedOdeEquationsImplicit<Eqn>,
pub fn new( problem: &'a OdeSolverProblem<Eqn>, state: RkState<Eqn::V>, tableau: Tableau<M>, linear_solver: LS, ) -> Result<Self, DiffsolError>
pub fn new_augmented( problem: &'a OdeSolverProblem<Eqn>, state: RkState<Eqn::V>, tableau: Tableau<M>, linear_solver: LS, augmented_eqn: AugmentedEqn, ) -> Result<Self, DiffsolError>
pub fn get_statistics(&self) -> &BdfStatistics
Trait Implementations§
Source§impl<'a, M, Eqn, LS, AugEqn> AugmentedOdeSolverMethod<'a, Eqn, AugEqn> for Sdirk<'a, Eqn, LS, M, AugEqn>where
Eqn: OdeEquationsImplicit,
AugEqn: AugmentedOdeEquationsImplicit<Eqn>,
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
LS: LinearSolver<Eqn::M>,
Eqn::V: DefaultDenseMatrix<T = Eqn::T, C = Eqn::C>,
for<'b> &'b Eqn::V: VectorRef<Eqn::V>,
for<'b> &'b Eqn::M: MatrixRef<Eqn::M>,
impl<'a, M, Eqn, LS, AugEqn> AugmentedOdeSolverMethod<'a, Eqn, AugEqn> for Sdirk<'a, Eqn, LS, M, AugEqn>where
Eqn: OdeEquationsImplicit,
AugEqn: AugmentedOdeEquationsImplicit<Eqn>,
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
LS: LinearSolver<Eqn::M>,
Eqn::V: DefaultDenseMatrix<T = Eqn::T, C = Eqn::C>,
for<'b> &'b Eqn::V: VectorRef<Eqn::V>,
for<'b> &'b Eqn::M: MatrixRef<Eqn::M>,
fn into_state_and_eqn(self) -> (Self::State, Option<AugEqn>)
fn augmented_eqn(&self) -> Option<&AugEqn>
Source§impl<M, Eqn, LS, AugmentedEqn> Clone for Sdirk<'_, Eqn, LS, M, AugmentedEqn>where
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
LS: LinearSolver<Eqn::M>,
Eqn: OdeEquationsImplicit,
Eqn::V: DefaultDenseMatrix<T = Eqn::T, C = Eqn::C>,
AugmentedEqn: AugmentedOdeEquationsImplicit<Eqn>,
impl<M, Eqn, LS, AugmentedEqn> Clone for Sdirk<'_, Eqn, LS, M, AugmentedEqn>where
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
LS: LinearSolver<Eqn::M>,
Eqn: OdeEquationsImplicit,
Eqn::V: DefaultDenseMatrix<T = Eqn::T, C = Eqn::C>,
AugmentedEqn: AugmentedOdeEquationsImplicit<Eqn>,
Source§impl<'a, M, Eqn, AugmentedEqn, LS> OdeSolverMethod<'a, Eqn> for Sdirk<'a, Eqn, LS, M, AugmentedEqn>where
LS: LinearSolver<Eqn::M>,
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
Eqn: OdeEquationsImplicit,
Eqn::V: DefaultDenseMatrix<T = Eqn::T, C = Eqn::C>,
AugmentedEqn: AugmentedOdeEquationsImplicit<Eqn>,
impl<'a, M, Eqn, AugmentedEqn, LS> OdeSolverMethod<'a, Eqn> for Sdirk<'a, Eqn, LS, M, AugmentedEqn>where
LS: LinearSolver<Eqn::M>,
M: DenseMatrix<T = Eqn::T, V = Eqn::V, C = Eqn::C>,
Eqn: OdeEquationsImplicit,
Eqn::V: DefaultDenseMatrix<T = Eqn::T, C = Eqn::C>,
AugmentedEqn: AugmentedOdeEquationsImplicit<Eqn>,
type State = RkState<<Eqn as Op>::V>
Source§fn problem(&self) -> &'a OdeSolverProblem<Eqn>
fn problem(&self) -> &'a OdeSolverProblem<Eqn>
Get the current problem
Source§fn jacobian(&self) -> Option<Ref<'_, Eqn::M>>
fn jacobian(&self) -> Option<Ref<'_, Eqn::M>>
Returns the current jacobian matrix of the solver, if it has one
Note that this will force a full recalculation of the Jacobian.
Source§fn mass(&self) -> Option<Ref<'_, Eqn::M>>
fn mass(&self) -> Option<Ref<'_, Eqn::M>>
Returns the current mass matrix of the solver, if it has one
Note that this will force a full recalculation of the mass matrix.
Source§fn order(&self) -> usize
fn order(&self) -> usize
Get the current order of accuracy of the solver (e.g. explict euler method is first-order)
Source§fn set_state(&mut self, state: Self::State)
fn set_state(&mut self, state: Self::State)
Replace the current state of the solver with a new state.
Source§fn into_state(self) -> RkState<Eqn::V>
fn into_state(self) -> RkState<Eqn::V>
Take the current state of the solver, if it exists, returning it to the user. This is useful if you want to use this
state in another solver or problem. Note that this will unset the current problem and solver state, so you will need to call
set_problem
again before calling step
or solve
.Source§fn checkpoint(&mut self) -> Self::State
fn checkpoint(&mut self) -> Self::State
Take a checkpoint of the current state of the solver, returning it to the user. This is useful if you want to use this
state in another solver or problem but want to keep this solver active. If you don’t need to use this solver again, you can use
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>
Step the solution forward by one step, altering the internal state of the solver.
The return value is a
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>
Set a stop time for the solver. The solver will stop when the internal time reaches this time.
Once it stops, the stop time is unset. If
tstop
is at or before the current internal time, an error is returned.Source§fn interpolate_sens(
&self,
t: <Eqn as Op>::T,
) -> Result<Vec<<Eqn as Op>::V>, DiffsolError>
fn interpolate_sens( &self, t: <Eqn as Op>::T, ) -> Result<Vec<<Eqn as Op>::V>, DiffsolError>
Interpolate the sensitivity vectors at a given time. This time should be between the current time and the last solver time step
Source§fn interpolate(&self, t: Eqn::T) -> Result<Eqn::V, DiffsolError>
fn interpolate(&self, t: Eqn::T) -> Result<Eqn::V, DiffsolError>
Interpolate the solution at a given time. This time should be between the current time and the last solver time step
Source§fn interpolate_out(&self, t: Eqn::T) -> Result<Eqn::V, DiffsolError>
fn interpolate_out(&self, t: Eqn::T) -> Result<Eqn::V, DiffsolError>
Interpolate the integral of the output function at a given time. This time should be between the current time and the last solver time step
Source§fn state_mut(&mut self) -> StateRefMut<'_, Eqn::V>
fn state_mut(&mut self) -> StateRefMut<'_, Eqn::V>
Get a mutable reference to the current state of the solver
Note that calling this will cause the next call to
step
to perform some reinitialisation to take into
account the mutated state, this could be expensive for multi-step methods.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>
Using the provided state, solve the problem up to time
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>
Using the provided state, solve the problem up to time
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>
Solve the problem and write out the solution at the given timepoints, using checkpointing so that
the solution can be interpolated at any timepoint.
See Self::solve_dense for a similar method that does not use checkpointing.
Auto Trait Implementations§
impl<'a, Eqn, LS, M = <<Eqn as Op>::V as DefaultDenseMatrix>::M, AugmentedEqn = NoAug<Eqn>> !Freeze for Sdirk<'a, Eqn, LS, M, AugmentedEqn>
impl<'a, Eqn, LS, M = <<Eqn as Op>::V as DefaultDenseMatrix>::M, AugmentedEqn = NoAug<Eqn>> !RefUnwindSafe for Sdirk<'a, Eqn, LS, M, AugmentedEqn>
impl<'a, Eqn, LS, M, AugmentedEqn> Send for Sdirk<'a, Eqn, LS, M, AugmentedEqn>
impl<'a, Eqn, LS, M = <<Eqn as Op>::V as DefaultDenseMatrix>::M, AugmentedEqn = NoAug<Eqn>> !Sync for Sdirk<'a, Eqn, LS, M, AugmentedEqn>
impl<'a, Eqn, LS, M, AugmentedEqn> Unpin for Sdirk<'a, Eqn, LS, M, AugmentedEqn>
impl<'a, Eqn, LS, M, AugmentedEqn> UnwindSafe for Sdirk<'a, Eqn, LS, M, AugmentedEqn>where
<Eqn as Op>::V: Sized + UnwindSafe + RefUnwindSafe,
M: UnwindSafe,
LS: UnwindSafe,
<Eqn as Op>::T: UnwindSafe + RefUnwindSafe,
AugmentedEqn: UnwindSafe,
Eqn: RefUnwindSafe,
<<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>
Backwards pass for adjoint sensitivity analysis Read more
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
Mutably borrows from an owned value. Read more
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>
Converts
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>
Converts
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,
Using the provided state, solve the problem up to time
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>
The inverse inclusion map: attempts to construct
self
from the equivalent element of its
superset. Read moreSource§fn is_in_subset(&self) -> bool
fn is_in_subset(&self) -> bool
Checks if
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
Use with care! Same as
self.to_subset
but without any property checks. Always succeeds.Source§fn from_subset(element: &SS) -> SP
fn from_subset(element: &SS) -> SP
The inclusion map: converts
self
to the equivalent element of its superset.