Rythmos - Transient Integration for Differential Equations
Version of the Day
|
Foward sensitivity stepper concrete subclass. More...
#include <Rythmos_ForwardSensitivityStepper.hpp>
Public Types | |
typedef Teuchos::ScalarTraits < Scalar >::magnitudeType | ScalarMag |
Public Types inherited from Rythmos::InterpolationBufferBase< Scalar > | |
typedef Teuchos::ScalarTraits < Scalar >::magnitudeType | ScalarMag |
Related Functions | |
(Note that these are not member functions.) | |
template<class Scalar > | |
RCP< ForwardSensitivityStepper < Scalar > > | forwardSensitivityStepper () |
Nonmember constructor. More... | |
template<class Scalar > | |
RCP< ForwardSensitivityStepper < Scalar > > | forwardSensitivityStepper (const RCP< const Thyra::ModelEvaluator< Scalar > > &stateModel, const int p_index, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &stateBasePoint, const RCP< StepperBase< Scalar > > &stateStepper, const RCP< Thyra::NonlinearSolverBase< Scalar > > &stateTimeStepSolver, const RCP< StepperBase< Scalar > > &sensStepper=Teuchos::null, const RCP< Thyra::NonlinearSolverBase< Scalar > > &sensTimeStepSolver=Teuchos::null) |
Nonmember constructor. More... | |
template<class Scalar > | |
int | getParameterIndex (const ForwardSensitivityStepper< Scalar > &fwdSensStepper) |
Return the index of the parameter subvector in the underlying state model. More... | |
template<class Scalar > | |
Thyra::ModelEvaluatorBase::InArgs < Scalar > | createStateAndSensInitialCondition (const ForwardSensitivityStepper< Scalar > &fwdSensStepper, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &state_ic, const RCP< const Thyra::MultiVectorBase< Scalar > > S_init=Teuchos::null, const RCP< const Thyra::MultiVectorBase< Scalar > > S_dot_init=Teuchos::null) |
Set up default initial conditions for the state and sensitivity stepper with default zero initial conditions for the sensitivity quantities. More... | |
template<class Scalar > | |
Thyra::ModelEvaluatorBase::InArgs < Scalar > | extractStateInitialCondition (const ForwardSensitivityStepper< Scalar > &fwdSensStepper, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &state_and_sens_ic) |
Extract out the initial condition for just the state model given the initial condition for the state and sensitivity model. More... | |
Related Functions inherited from Rythmos::StepperBase< Scalar > | |
template<class Scalar > | |
bool | isInitialized (const StepperBase< Scalar > &stepper) |
template<class Scalar > | |
bool | isInitialized (const StepperBase< Scalar > &stepper) |
Related Functions inherited from Rythmos::InterpolationBufferBase< Scalar > | |
template<class Scalar > | |
RCP< const Thyra::VectorBase < Scalar > > | get_x (const InterpolationBufferBase< Scalar > &interpBuffer, const Scalar &t) |
Get a single point x(t) from an interpolation buffer. More... | |
template<class Scalar > | |
RCP< const Thyra::VectorBase < Scalar > > | get_xdot (const InterpolationBufferBase< Scalar > &interpBuffer, const Scalar &t) |
Get a single point xdot(t) from an interpolation buffer. More... | |
template<class Scalar > | |
void | get_x_and_x_dot (const InterpolationBufferBase< Scalar > &interpBuffer, const Scalar t, const Ptr< RCP< const Thyra::VectorBase< Scalar > > > &x, const Ptr< RCP< const Thyra::VectorBase< Scalar > > > &x_dot) |
Nonmember helper function to get x and x_dot at t. More... | |
template<class Scalar > | |
void | assertTimePointsAreSorted (const Array< Scalar > &time_vec) |
Assert that a time point vector is sorted. More... | |
template<class Scalar > | |
void | assertNoTimePointsBeforeCurrentTimeRange (const InterpolationBufferBase< Scalar > &interpBuffer, const Array< Scalar > &time_vec, const int &startingTimePointIndex=0) |
Assert that none of the time points fall before the current time range for an interpolation buffer object. More... | |
template<class Scalar > | |
void | assertNoTimePointsInsideCurrentTimeRange (const InterpolationBufferBase< Scalar > &interpBuffer, const Array< Scalar > &time_vec) |
Assert that none of the time points fall inside the current time range for an interpolation buffer object. More... | |
template<class TimeType > | |
void | selectPointsInTimeRange (const Array< TimeType > &points_in, const TimeRange< TimeType > &range, const Ptr< Array< TimeType > > &points_out) |
Select points from an Array that sit in a TimeRange. More... | |
template<class TimeType > | |
void | removePointsInTimeRange (Array< TimeType > *points_in, const TimeRange< TimeType > &range) |
Remove points from an Array that sit in a TimeRange. More... | |
template<class Scalar > | |
bool | getCurrentPoints (const InterpolationBufferBase< Scalar > &interpBuffer, const Array< Scalar > &time_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *x_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *xdot_vec, int *nextTimePointIndex) |
Get time points in the current range of an interpolation buffer object. More... | |
Constructors, Intializers, Misc. | |
ForwardSensitivityStepper () | |
Constructs to uninitialized. More... | |
void | initializeSyncedSteppers (const RCP< const Thyra::ModelEvaluator< Scalar > > &stateModel, const int p_index, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &stateBasePoint, const RCP< StepperBase< Scalar > > &stateStepper, const RCP< Thyra::NonlinearSolverBase< Scalar > > &stateTimeStepSolver, const RCP< StepperBase< Scalar > > &sensStepper=Teuchos::null, const RCP< Thyra::NonlinearSolverBase< Scalar > > &sensTimeStepSolver=Teuchos::null) |
Intialize for synced state and sens steppers. More... | |
void | initializeSyncedSteppersInitCondOnly (const RCP< const Thyra::ModelEvaluator< Scalar > > &stateModel, const RCP< const Thyra::VectorSpaceBase< Scalar > > &p_space, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &stateBasePoint, const RCP< StepperBase< Scalar > > &stateStepper, const RCP< Thyra::NonlinearSolverBase< Scalar > > &stateTimeStepSolver, const RCP< StepperBase< Scalar > > &sensStepper=Teuchos::null, const RCP< Thyra::NonlinearSolverBase< Scalar > > &sensTimeStepSolver=Teuchos::null) |
Intialize for synced state and sens steppers for an initial-condition only parametrized sensitivity problem. More... | |
void | initializeDecoupledSteppers (const RCP< const Thyra::ModelEvaluator< Scalar > > &stateModel, const int p_index, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &stateBasePoint, const RCP< StepperBase< Scalar > > &stateStepper, const RCP< Thyra::NonlinearSolverBase< Scalar > > &stateTimeStepSolver, const RCP< IntegratorBase< Scalar > > &stateIntegrator, const Scalar &finalTime, const RCP< StepperBase< Scalar > > &sensStepper=Teuchos::null, const RCP< Thyra::NonlinearSolverBase< Scalar > > &sensTimeStepSolver=Teuchos::null) |
Intialize for decoupled state and sens steppers. More... | |
bool | stateModelIsConst () const |
Return if the state model is const-only or not. More... | |
RCP< const Thyra::ModelEvaluator< Scalar > > | getStateModel () const |
Return the state model that was passed into an initialize function. More... | |
RCP< StepperBase< Scalar > > | getNonconstStateStepper () |
Return the state stepper that was passed into an initialize function. More... | |
RCP< const ForwardSensitivityModelEvaluatorBase < Scalar > > | getFwdSensModel () const |
Return the forward sensitivity model evaluator object that got created internally when the initialize function was called. More... | |
RCP< const StateAndForwardSensitivityModelEvaluator < Scalar > > | getStateAndFwdSensModel () const |
Return the state and forward sensitivity model evaluator object that got created internally when the nitialize function was called. More... | |
Overridden from Teuchos::ParameterListAcceptor | |
void | setParameterList (RCP< Teuchos::ParameterList > const ¶mList) |
RCP< const Teuchos::ParameterList > | getValidParameters () const |
Overridden from StepperBase | |
bool | acceptsModel () const |
Returns false. More... | |
void | setModel (const RCP< const Thyra::ModelEvaluator< Scalar > > &model) |
Throws exception. More... | |
void | setNonconstModel (const RCP< Thyra::ModelEvaluator< Scalar > > &model) |
Throws exception. More... | |
RCP< const Thyra::ModelEvaluator< Scalar > > | getModel () const |
Returns getStateAndFwdSensModel() . More... | |
RCP< Thyra::ModelEvaluator < Scalar > > | getNonconstModel () |
void | setInitialCondition (const Thyra::ModelEvaluatorBase::InArgs< Scalar > &state_and_sens_ic) |
Sets the full initial condition for x_bar = [ x; s_bar] and x_bar_dot = [ x_dot; s_bar_dot ] as well as the initial time and the parameter values. More... | |
Thyra::ModelEvaluatorBase::InArgs < Scalar > | getInitialCondition () const |
Scalar | takeStep (Scalar dt, StepSizeType stepType) |
const StepStatus< Scalar > | getStepStatus () const |
Overridden from InterpolationBufferBase | |
RCP< const Thyra::VectorSpaceBase< Scalar > > | get_x_space () const |
Returns the space for x_bar and x_bar_dot . More... | |
void | addPoints (const Array< Scalar > &time_vec, const Array< Teuchos::RCP< const Thyra::VectorBase< Scalar > > > &x_vec, const Array< Teuchos::RCP< const Thyra::VectorBase< Scalar > > > &xdot_vec) |
TimeRange< Scalar > | getTimeRange () const |
void | getPoints (const Array< Scalar > &time_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *x_vec, Array< RCP< const Thyra::VectorBase< Scalar > > > *xdot_vec, Array< ScalarMag > *accuracy_vec) const |
void | getNodes (Array< Scalar > *time_vec) const |
void | removeNodes (Array< Scalar > &time_vec) |
int | getOrder () const |
Deprecated. | |
void | initialize (const RCP< const Thyra::ModelEvaluator< Scalar > > &stateModel, const int p_index, const Thyra::ModelEvaluatorBase::InArgs< Scalar > &stateBasePoint, const RCP< StepperBase< Scalar > > &stateStepper, const RCP< Thyra::NonlinearSolverBase< Scalar > > &stateTimeStepSolver, const RCP< StepperBase< Scalar > > &sensStepper=Teuchos::null, const RCP< Thyra::NonlinearSolverBase< Scalar > > &sensTimeStepSolver=Teuchos::null) |
Deprecated. More... | |
Additional Inherited Members | |
Public Member Functions inherited from Rythmos::StepperBase< Scalar > | |
virtual bool | supportsCloning () const |
Return if this stepper supports cloning or not. More... | |
virtual RCP< StepperBase < Scalar > > | cloneStepperAlgorithm () const |
Clone the stepper object if supported. More... | |
virtual bool | isImplicit () const |
Return if this stepper is an implicit stepper. More... | |
virtual bool | modelIsConst () const |
Return of the model is only const or can be returned as a non-const object. More... | |
virtual void | setStepControlData (const StepperBase &stepper) |
Set step control data from another stepper. More... | |
Public Member Functions inherited from Rythmos::InterpolationBufferBase< Scalar > | |
virtual void | addPoints (const Array< Scalar > &time_vec, const Array< RCP< const Thyra::VectorBase< Scalar > > > &x_vec, const Array< RCP< const Thyra::VectorBase< Scalar > > > &xdot_vec)=0 |
Add points to the buffer. More... | |
Foward sensitivity stepper concrete subclass.
This class provides a very general implemenation of a forward sensitivity stepper.
The most general form of the parameterized state equation is:
f(x_dot(t),x(t),p) = 0, over t = [t0,tf] x(t0) = x_init(p)
As shown above, the parameters p
are assumed to be steady state parameters and can enter through the intial condition and/or through the DAE equation itself. This class also supports a form of the problem where parameters p
are only assumed to bin the initial condition x_init(p)
and not in the state equation. In this case, you can just drop out all of the terms d(f)/d(p)
shown in the equations below because they are zero.
The forward sensitivity equations that are solved along with the state equation, written in multi-vector form, are:
d(f)/d(x_dot)*S_dot + d(f)/d(x)*S + d(f)/d(p) = 0, over t = [t0,tf] S(t0) = d(x_init)/d(p)
where S
is a multi-vector with np
columns where each column S(:,j) = d(x)/d(p_j)
is the sensitivity of x(t)
with respect to the p_j
parameter.
The forward sensitivity equations are a DAE and must be solved using a time integrator. Conceptually, the state plus forward sensitivity system can be throught of as a big composite DAE model of the form:
f_bar(x_bar_dot(t),x_bar(t)) = 0, over t = [t0,tf] x_bar(t0) = x_bar_init
where
x_bar = [ x; s_bar ] s_bar = [ S(:,0); S(:,0); ...; S(:,np-1) ]
and f_bar(...)
is the obvious concatenated state and sensitivity system. See the class StateAndForwardSensitivityModelEvaluatorBase
for a description of how to get at the components of x
, s_bar
, and S
given x_bar
.
The InterpolationBufferBase
interface implemented by this class is defined with respect to the full composite solution vector x_bar
which is returned as a product vector with two components. The first component is x
. The second component is another product vector for the implicit concatenation of the columns of the sensitivities shown as s_bar
above. The s_bar
product vector is really implemented directly as a multi-vector and represents an efficient way to represent the forward sensitivities. Therefore, through the interpolation buffer interface function getPoints()
a client can access the state and the sensitivities at any point in the range of the current timestep.
Note that this class does not implement the function setModel()
since any arbitrary combined state + sensitivity model can not be supported.
There are a variety of ways that one can go about implementing a state plus forward sensitivity stepper. Three ways for doing this are described in the report "Design of New DASPK for Sensitivity Analysis" by Shengtai Li and Linda Petzold. The three ways are the simultaneous corrector, the staggered direct and the staggered corrector methods.
The simultaneous corrector method would be equivalent to forming one big ModelEvaluator for the state and sensitivities where the "state" variables would be the x_bar
variables described above and then this big system would be solved with a single stepper object and as single nonlinear solver. The advantage of this approach is that it makes great reuse of all of the timestepping software. Also, by being able to specialize the nonlinear solver (which you can't do in most software) you could set up the nonlinear solver to first solve the nonlinear state timestep equation, and then solve the linear sensitivity equations. The problem with this approach would be that it would be very wasteful if the timestep had to be cut back in order to reduce the local truncation error of the state solution. This would result in the sensitivity solution being thrown away for each cut-back iteration. Because of this fundamental problem, we have not implemented the simultaneous corrector method. Actually, we are not really sure why anyone implements ths method.
The staggered direct and staggered corrector methods are similar in several ways. In each method, the state timestep is first fully solved (including any stepsize reduction iterations that are needed) and then the sensitivities are solved, taking advantage of the pre-computed state timestep Jacobian. One difference between the two methods is that in the staggered direct method, the sensitivities are solved for directly. This can result in numerical accuracy problems and does not allow the reuse of an inexact solve when a direct factorization is used. In the staggered corrector method, an implicit corrector solve is used to compute the change in the sensitivity variables using a Newton method. This results in better numerical stability and allows the reuse of an old Jacobian and factorization in the Newton method. However, if an exact Jacobian is used for the solve, then the Newton method will converge in one iteration (assuming the linear solver tolerance is made tight enough) with no harm done.
Because of the advantages of the staggered corrector method over the other methods, the staggered corrector method is what is implemented in this stepper class. However, the term "corrector" is not really appropriate to be used in the context of this class since this class does not have to assume anything about how timesteps are computed and does not care if a predictor/corrector method is used or not.
While this class does provide a full ModelEvaluator for the full state plus forward sensitivity DAE f_bar(x_bar_hat,x_bar)
it is not solved all at one time as described above. Instead, the step is solved first for the state equation and then a ModelEvaluator for just the linear forward sensitivity equations is formed and is solved over the same time step as the forward solve.
Currently, timestep control is not performed for the forward sensitivity variables. In the future, however, it would not be too difficult to allow for the timestep to be reduced for the sensitivity variables but this would require a "undoStep()" operation be implimented for the state stepper object and this is not currently supported by the StepperBase
interface.
2007/15/21: rabart: ToDo: This class only works for implicit models and steppers right now but it would be easy to get this to work for explicit steppers and models later with a little work. With an explicit method and model, we don't need to reuse W_tilde so this is easier in a way!
ToDo: Finish documentation!
Definition at line 209 of file Rythmos_ForwardSensitivityStepper.hpp.
typedef Teuchos::ScalarTraits<Scalar>::magnitudeType Rythmos::ForwardSensitivityStepper< Scalar >::ScalarMag |
Definition at line 216 of file Rythmos_ForwardSensitivityStepper.hpp.
Rythmos::ForwardSensitivityStepper< Scalar >::ForwardSensitivityStepper | ( | ) |
Constructs to uninitialized.
Definition at line 780 of file Rythmos_ForwardSensitivityStepper.hpp.
void Rythmos::ForwardSensitivityStepper< Scalar >::initializeSyncedSteppers | ( | const RCP< const Thyra::ModelEvaluator< Scalar > > & | stateModel, |
const int | p_index, | ||
const Thyra::ModelEvaluatorBase::InArgs< Scalar > & | stateBasePoint, | ||
const RCP< StepperBase< Scalar > > & | stateStepper, | ||
const RCP< Thyra::NonlinearSolverBase< Scalar > > & | stateTimeStepSolver, | ||
const RCP< StepperBase< Scalar > > & | sensStepper = Teuchos::null , |
||
const RCP< Thyra::NonlinearSolverBase< Scalar > > & | sensTimeStepSolver = Teuchos::null |
||
) |
Intialize for synced state and sens steppers.
stateModel | [in,persisting] The ModelEvaluator that defines the parameterized state model f(x_dot,x,p) . |
p_index | [in] The index of the parameter subvector in stateModel for which sensitivities will be computed for. |
baseStatePoint | [in] Whatever input arguments are needed to define the state of the model including the parameters except x, x_dot, and t! |
stateStepper | [in,persisting] The stepper object that will be used to advance the state solution x(t) . This stepper need not be setup with a model or a nonlinear timestep solver. All this stepper object needs is to be given its parameters to determine exactly what timestepping algorithm will be employed. The model and the timestep solver objects will be set internally. |
stateTimeStepSolver | [in,persisting] The nonlinear solver object that is used to solve for the state timestep equation. This is needed to extract the Jacobian matrix that is used in the sensitivity model. If the stepper is not an implicit stepper and does not use an implicit time step solver, then this argument should be left null. |
sensStepper | [in,persisting] The stepper object that will be used to advance the sensitivity solution S(t) . This stepper need not be setup with a model or a nonlinear timestep solver. All this stepper object needs is to be given its parameters to determine exactly what timestepping algorithm will be employed. The model and the timestep solver objects will be set internally. If this argument is null, then the stateStepper object will be cloned to generate this stepper object. The most common use cases should just pass in Teuchos::null and just use the identical stepper as the state stepper. However, this argument allows a client to specialize exactly what the sensitivity stepper does and therefore this hook is allowed. |
sensTimeStepSolver | [in,persisting] The nonlinear solver object that is used to solve for the (linear) sensitivity timestep equation. If the stepper is not an implicit stepper and does not use an implicit timestep solver, then this argument can be left null. If the stepper is implicit, and this argument is left null, then a Thyra::LinearNonlinearSolver object will be created and used. The most common use cases should just pass in Teuchos::null and just use the simple linear nonlinear solver to will perform just a single linear solve. However, this argument allows a client to specialize exactly what the nonlinear solver in the sensitivity stepper does and therefore this hook is exposed to clients. |
Here *this
is set up to synchronize the state and sensitivity solvers. Currently, error control is only done by the state stepper and not the sensitivity stepper but the overall implementation has a high degree of resuse and will therefore compute sensitivities quite fast.
Definition at line 787 of file Rythmos_ForwardSensitivityStepper.hpp.
void Rythmos::ForwardSensitivityStepper< Scalar >::initializeSyncedSteppersInitCondOnly | ( | const RCP< const Thyra::ModelEvaluator< Scalar > > & | stateModel, |
const RCP< const Thyra::VectorSpaceBase< Scalar > > & | p_space, | ||
const Thyra::ModelEvaluatorBase::InArgs< Scalar > & | stateBasePoint, | ||
const RCP< StepperBase< Scalar > > & | stateStepper, | ||
const RCP< Thyra::NonlinearSolverBase< Scalar > > & | stateTimeStepSolver, | ||
const RCP< StepperBase< Scalar > > & | sensStepper = Teuchos::null , |
||
const RCP< Thyra::NonlinearSolverBase< Scalar > > & | sensTimeStepSolver = Teuchos::null |
||
) |
Intialize for synced state and sens steppers for an initial-condition only parametrized sensitivity problem.
stateModel | [in,persisting] See initializeSyncedSteppers(). |
p_space | [in] The vector space for the parameterized initial condition parameters. |
baseStatePoint | [in] See initializeSyncedSteppers(). |
stateStepper | [in,persisting] See initializeSyncedSteppers(). |
stateTimeStepSolver | [in,persisting] See initializeSyncedSteppers(). |
sensStepper | [in,persisting] See initializeSyncedSteppers(). |
sensTimeStepSolver | [in,persisting] See initializeSyncedSteppers(). |
Here *this
is set up to synchronize the state and sensitivity solvers for an initial-condition only forward sensitivity problem. Currently, error control is only done by the state stepper and not the sensitivity stepper but the overall implementation has a high degree of resuse and will therefore compute sensitivities quite fast.
Definition at line 804 of file Rythmos_ForwardSensitivityStepper.hpp.
void Rythmos::ForwardSensitivityStepper< Scalar >::initializeDecoupledSteppers | ( | const RCP< const Thyra::ModelEvaluator< Scalar > > & | stateModel, |
const int | p_index, | ||
const Thyra::ModelEvaluatorBase::InArgs< Scalar > & | stateBasePoint, | ||
const RCP< StepperBase< Scalar > > & | stateStepper, | ||
const RCP< Thyra::NonlinearSolverBase< Scalar > > & | stateTimeStepSolver, | ||
const RCP< IntegratorBase< Scalar > > & | stateIntegrator, | ||
const Scalar & | finalTime, | ||
const RCP< StepperBase< Scalar > > & | sensStepper = Teuchos::null , |
||
const RCP< Thyra::NonlinearSolverBase< Scalar > > & | sensTimeStepSolver = Teuchos::null |
||
) |
Intialize for decoupled state and sens steppers.
stateModel | [in,persisting] See initializeSyncedSteppers() . |
p_index | [in] See initializeSyncedSteppers() . |
baseStatePoint | [in] See initializeSyncedSteppers() . |
stateStepper | [in,persisting] See initializeSyncedSteppers() . |
stateTimeStepSolver | [in,persisting] See initializeSyncedSteppers() . |
stateIntegrator | [in,persisting] The intergrator that will be used to integrate the state given stateStepper . This integrator must be set up with a trailing interpolation buffer in order to be able to allow for complete flexibility in how the time steps for the sens equations are solved. |
finalTime | [in] The final time for the state integrator. This is needed to initialize stateIntegrator with stateStepper . |
sensStepper | [in,persisting] See initializeSyncedSteppers() . |
sensTimeStepSolver | [in,persisting] See initializeSyncedSteppers() . |
Here *this
is set up to run the state and sens steppers completely independently; each with the their own error control strategies. The state stepper in driven through the state integrator which in turn is driven by the ForwardSensitivityModelEvaluatorBase that is driven by the sens stepper.
Definition at line 820 of file Rythmos_ForwardSensitivityStepper.hpp.
bool Rythmos::ForwardSensitivityStepper< Scalar >::stateModelIsConst | ( | ) | const |
Return if the state model is const-only or not.
Definition at line 841 of file Rythmos_ForwardSensitivityStepper.hpp.
RCP< const Thyra::ModelEvaluator< Scalar > > Rythmos::ForwardSensitivityStepper< Scalar >::getStateModel | ( | ) | const |
Return the state model that was passed into an initialize function.
Definition at line 849 of file Rythmos_ForwardSensitivityStepper.hpp.
RCP< StepperBase< Scalar > > Rythmos::ForwardSensitivityStepper< Scalar >::getNonconstStateStepper | ( | ) |
Return the state stepper that was passed into an initialize function.
Definition at line 857 of file Rythmos_ForwardSensitivityStepper.hpp.
RCP< const ForwardSensitivityModelEvaluatorBase< Scalar > > Rythmos::ForwardSensitivityStepper< Scalar >::getFwdSensModel | ( | ) | const |
Return the forward sensitivity model evaluator object that got created internally when the initialize function was called.
Definition at line 865 of file Rythmos_ForwardSensitivityStepper.hpp.
RCP< const StateAndForwardSensitivityModelEvaluator< Scalar > > Rythmos::ForwardSensitivityStepper< Scalar >::getStateAndFwdSensModel | ( | ) | const |
Return the state and forward sensitivity model evaluator object that got created internally when the nitialize function was called.
This is also the same model returned by the function getModel()
, except through it's concrete subclass type.
Definition at line 873 of file Rythmos_ForwardSensitivityStepper.hpp.
void Rythmos::ForwardSensitivityStepper< Scalar >::setParameterList | ( | RCP< Teuchos::ParameterList > const & | paramList | ) |
Definition at line 883 of file Rythmos_ForwardSensitivityStepper.hpp.
RCP< const Teuchos::ParameterList > Rythmos::ForwardSensitivityStepper< Scalar >::getValidParameters | ( | ) | const |
Definition at line 897 of file Rythmos_ForwardSensitivityStepper.hpp.
|
virtual |
Returns false.
Reimplemented from Rythmos::StepperBase< Scalar >.
Definition at line 919 of file Rythmos_ForwardSensitivityStepper.hpp.
|
virtual |
Throws exception.
Implements Rythmos::StepperBase< Scalar >.
Definition at line 925 of file Rythmos_ForwardSensitivityStepper.hpp.
|
virtual |
Throws exception.
Implements Rythmos::StepperBase< Scalar >.
Definition at line 936 of file Rythmos_ForwardSensitivityStepper.hpp.
|
virtual |
Returns getStateAndFwdSensModel()
.
Warning, currently the returned model does not implement evalModel(...) or define a W object. It is just used for getting the spaces and for creating an InArgs object for setting the initial condition.
Implements Rythmos::StepperBase< Scalar >.
Definition at line 948 of file Rythmos_ForwardSensitivityStepper.hpp.
|
virtual |
Implements Rythmos::StepperBase< Scalar >.
Definition at line 956 of file Rythmos_ForwardSensitivityStepper.hpp.
|
virtual |
Sets the full initial condition for x_bar = [ x; s_bar]
and x_bar_dot = [ x_dot; s_bar_dot ]
as well as the initial time and the parameter values.
The InArgs object must be created using this->getModel()->createInArgs()
and then populated with the initial values. The product vectors for x_bar
and x_bar_dot
can be created using this->getStateAndFwdSensModel()->create_x_bar_vec(...)
. All of the input objects in state_and_sens_ic
will be cloned and therefore no memory of the objects in state_and_sens_ic
will be retained after calling this function.
Implements Rythmos::StepperBase< Scalar >.
Definition at line 963 of file Rythmos_ForwardSensitivityStepper.hpp.
|
virtual |
Implements Rythmos::StepperBase< Scalar >.
Definition at line 1038 of file Rythmos_ForwardSensitivityStepper.hpp.
|
virtual |
Implements Rythmos::StepperBase< Scalar >.
Definition at line 1046 of file Rythmos_ForwardSensitivityStepper.hpp.
|
virtual |
Implements Rythmos::StepperBase< Scalar >.
Definition at line 1064 of file Rythmos_ForwardSensitivityStepper.hpp.
|
virtual |
Returns the space for x_bar
and x_bar_dot
.
This space is a nested product vector space as described above. Dynamic casting is required to get at the ProductVectorSapceBase
and ProductVectorBase
intefaces.
Implements Rythmos::InterpolationBufferBase< Scalar >.
Definition at line 1102 of file Rythmos_ForwardSensitivityStepper.hpp.
void Rythmos::ForwardSensitivityStepper< Scalar >::addPoints | ( | const Array< Scalar > & | time_vec, |
const Array< Teuchos::RCP< const Thyra::VectorBase< Scalar > > > & | x_vec, | ||
const Array< Teuchos::RCP< const Thyra::VectorBase< Scalar > > > & | xdot_vec | ||
) |
Definition at line 1109 of file Rythmos_ForwardSensitivityStepper.hpp.
|
virtual |
Implements Rythmos::InterpolationBufferBase< Scalar >.
Definition at line 1121 of file Rythmos_ForwardSensitivityStepper.hpp.
|
virtual |
Implements Rythmos::InterpolationBufferBase< Scalar >.
Definition at line 1128 of file Rythmos_ForwardSensitivityStepper.hpp.
|
virtual |
Implements Rythmos::InterpolationBufferBase< Scalar >.
Definition at line 1200 of file Rythmos_ForwardSensitivityStepper.hpp.
|
virtual |
Implements Rythmos::InterpolationBufferBase< Scalar >.
Definition at line 1219 of file Rythmos_ForwardSensitivityStepper.hpp.
|
virtual |
Implements Rythmos::InterpolationBufferBase< Scalar >.
Definition at line 1228 of file Rythmos_ForwardSensitivityStepper.hpp.
|
inline |
Deprecated.
Definition at line 518 of file Rythmos_ForwardSensitivityStepper.hpp.
|
related |
Nonmember constructor.
Definition at line 617 of file Rythmos_ForwardSensitivityStepper.hpp.
|
related |
Nonmember constructor.
Definition at line 630 of file Rythmos_ForwardSensitivityStepper.hpp.
|
related |
Return the index of the parameter subvector in the underlying state model.
Definition at line 654 of file Rythmos_ForwardSensitivityStepper.hpp.
|
related |
Set up default initial conditions for the state and sensitivity stepper with default zero initial conditions for the sensitivity quantities.
Definition at line 671 of file Rythmos_ForwardSensitivityStepper.hpp.
|
related |
Extract out the initial condition for just the state model given the initial condition for the state and sensitivity model.
Definition at line 735 of file Rythmos_ForwardSensitivityStepper.hpp.