Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar > Class Template Reference

Time integrator suitable for pseudotransient forward sensitivity analysis. More...

#include <Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp>

Inheritance diagram for Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >:
Tempus::Integrator< Scalar > Teuchos::Describable Teuchos::VerboseObject< Tempus::Integrator< Scalar > > Teuchos::LabeledObject Teuchos::VerboseObjectBase

Public Member Functions

 IntegratorPseudoTransientForwardSensitivity (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar >> &model, const Teuchos::RCP< SensitivityModelEvaluatorBase< Scalar >> &sens_model, const Teuchos::RCP< IntegratorBasic< Scalar >> &fwd_integrator, const Teuchos::RCP< IntegratorBasic< Scalar >> &sens_integrator, const bool reuse_solver, const bool force_W_update)
 Constructor with ParameterList and model, and will be fully initialized. More...
 
 IntegratorPseudoTransientForwardSensitivity ()
 Destructor. More...
 
virtual ~IntegratorPseudoTransientForwardSensitivity ()
 Destructor. More...
 
virtual void initializeSolutionHistory (Scalar t0, Teuchos::RCP< const Thyra::VectorBase< Scalar >> x0, Teuchos::RCP< const Thyra::VectorBase< Scalar >> xdot0=Teuchos::null, Teuchos::RCP< const Thyra::VectorBase< Scalar >> xdotdot0=Teuchos::null, Teuchos::RCP< const Thyra::MultiVectorBase< Scalar >> DxDp0=Teuchos::null, Teuchos::RCP< const Thyra::MultiVectorBase< Scalar >> DxdotDp0=Teuchos::null, Teuchos::RCP< const Thyra::MultiVectorBase< Scalar >> DxdotdotDp0=Teuchos::null)
 Set the initial state from Thyra::VectorBase(s) More...
 
virtual Teuchos::RCP< const
Thyra::VectorBase< Scalar > > 
getX () const
 Get current the solution, x. More...
 
virtual Teuchos::RCP< const
Thyra::MultiVectorBase< Scalar > > 
getDxDp () const
 
virtual Teuchos::RCP< const
Thyra::VectorBase< Scalar > > 
getXDot () const
 Get current the time derivative of the solution, xdot. More...
 
virtual Teuchos::RCP< const
Thyra::MultiVectorBase< Scalar > > 
getDXDotDp () const
 
virtual Teuchos::RCP< const
Thyra::VectorBase< Scalar > > 
getXDotDot () const
 Get current the second time derivative of the solution, xdotdot. More...
 
virtual Teuchos::RCP< const
Thyra::MultiVectorBase< Scalar > > 
getDXDotDotDp () const
 
virtual Teuchos::RCP< const
Thyra::VectorBase< Scalar > > 
getG () const
 Return response function g. More...
 
virtual Teuchos::RCP< const
Thyra::MultiVectorBase< Scalar > > 
getDgDp () const
 Return forward sensitivity stored in Jacobian format. More...
 
SensitivityStepMode getStepMode () const
 What mode the current time integration step is in. More...
 
- Public Member Functions inherited from Tempus::Integrator< Scalar >
- Public Member Functions inherited from Teuchos::Describable
void describe (std::ostream &out, const EVerbosityLevel verbLevel=verbLevel_default) const
 
virtual ~Describable ()
 
 LabeledObject ()
 
virtual ~LabeledObject ()
 
virtual void setObjectLabel (const std::string &objectLabel)
 
virtual std::string getObjectLabel () const
 
DescribableStreamManipulatorState describe (const Describable &describable, const EVerbosityLevel verbLevel=Describable::verbLevel_default)
 
std::ostream & operator<< (std::ostream &os, const DescribableStreamManipulatorState &d)
 
- Public Member Functions inherited from Teuchos::VerboseObject< Tempus::Integrator< Scalar > >
 VerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null)
 
virtual const VerboseObjectsetVerbLevel (const EVerbosityLevel verbLevel) const
 
virtual const VerboseObjectsetOverridingVerbLevel (const EVerbosityLevel verbLevel) const
 
virtual EVerbosityLevel getVerbLevel () const
 
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT
RCP< const ParameterList
getValidVerboseObjectSublist ()
 
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT
void 
setupVerboseObjectSublist (ParameterList *paramList)
 
TEUCHOSPARAMETERLIST_LIB_DLL_EXPORT
void 
readVerboseObjectSublist (ParameterList *paramList, RCP< FancyOStream > *oStream, EVerbosityLevel *verbLevel)
 
void readVerboseObjectSublist (ParameterList *paramList, VerboseObject< ObjectType > *verboseObject)
 
- Public Member Functions inherited from Teuchos::VerboseObjectBase
virtual ~VerboseObjectBase ()
 
 VerboseObjectBase (const RCP< FancyOStream > &oStream=Teuchos::null)
 
virtual const VerboseObjectBasesetOStream (const RCP< FancyOStream > &oStream) const
 
virtual const VerboseObjectBasesetOverridingOStream (const RCP< FancyOStream > &oStream) const
 
virtual VerboseObjectBasesetLinePrefix (const std::string &linePrefix)
 
virtual RCP< FancyOStreamgetOStream () const
 
virtual RCP< FancyOStreamgetOverridingOStream () const
 
virtual std::string getLinePrefix () const
 
virtual OSTab getOSTab (const int tabs=1, const std::string &linePrefix="") const
 

Protected Member Functions

void buildSolutionHistory ()
 
- Protected Member Functions inherited from Teuchos::VerboseObject< Tempus::Integrator< Scalar > >
void initializeVerboseObject (const EVerbosityLevel verbLevel=VERB_DEFAULT, const RCP< FancyOStream > &oStream=Teuchos::null)
 
- Protected Member Functions inherited from Teuchos::VerboseObjectBase
void initializeVerboseObjectBase (const RCP< FancyOStream > &oStream=Teuchos::null)
 
virtual void informUpdatedVerbosityState () const
 

Protected Attributes

Teuchos::RCP
< Thyra::ModelEvaluator
< Scalar > > 
model_
 
Teuchos::RCP
< SensitivityModelEvaluatorBase
< Scalar > > 
sens_model_
 
Teuchos::RCP< IntegratorBasic
< Scalar > > 
state_integrator_
 
Teuchos::RCP< IntegratorBasic
< Scalar > > 
sens_integrator_
 
Teuchos::RCP< SolutionHistory
< Scalar > > 
solutionHistory_
 
bool reuse_solver_
 
bool force_W_update_
 
SensitivityStepMode stepMode_
 

Basic integrator methods

virtual bool advanceTime ()
 Advance the solution to timeMax, and return true if successful. More...
 
virtual bool advanceTime (const Scalar timeFinal) override
 Advance the solution to timeFinal, and return true if successful. More...
 
virtual Scalar getTime () const override
 Get current time. More...
 
virtual int getIndex () const override
 Get current index. More...
 
virtual Status getStatus () const override
 Get Status. More...
 
virtual void setStatus (const Status st) override
 Set Status. More...
 
virtual Teuchos::RCP< Stepper
< Scalar > > 
getStepper () const override
 Get the Stepper. More...
 
Teuchos::RCP< Stepper< Scalar > > getStateStepper () const
 
Teuchos::RCP< Stepper< Scalar > > getSensStepper () const
 
virtual Teuchos::RCP< const
SolutionHistory< Scalar > > 
getSolutionHistory () const override
 Get the SolutionHistory. More...
 
Teuchos::RCP< const
SolutionHistory< Scalar > > 
getStateSolutionHistory () const
 
Teuchos::RCP< const
SolutionHistory< Scalar > > 
getSensSolutionHistory () const
 
virtual Teuchos::RCP
< SolutionHistory< Scalar > > 
getNonConstSolutionHistory () override
 Get the SolutionHistory. More...
 
virtual Teuchos::RCP< const
TimeStepControl< Scalar > > 
getTimeStepControl () const override
 Get the TimeStepControl. More...
 
virtual Teuchos::RCP
< TimeStepControl< Scalar > > 
getNonConstTimeStepControl () override
 
Teuchos::RCP< TimeStepControl
< Scalar > > 
getStateNonConstTimeStepControl ()
 
Teuchos::RCP< TimeStepControl
< Scalar > > 
getSensNonConstTimeStepControl ()
 
virtual Teuchos::RCP
< IntegratorObserver< Scalar > > 
getObserver ()
 Get the Observer. More...
 
virtual void setObserver (Teuchos::RCP< IntegratorObserver< Scalar >> obs=Teuchos::null)
 Set the Observer. More...
 
virtual Teuchos::RCP
< Teuchos::Time
getIntegratorTimer () const override
 Returns the IntegratorTimer_ for this Integrator. More...
 
virtual Teuchos::RCP
< Teuchos::Time
getStepperTimer () const override
 

Overridden from Teuchos::Describable

std::string description () const override
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
 

Additional Inherited Members

- Static Public Member Functions inherited from Teuchos::VerboseObject< Tempus::Integrator< Scalar > >
static void setDefaultVerbLevel (const EVerbosityLevel defaultVerbLevel)
 
static EVerbosityLevel getDefaultVerbLevel ()
 
- Static Public Member Functions inherited from Teuchos::VerboseObjectBase
static void setDefaultOStream (const RCP< FancyOStream > &defaultOStream)
 
static RCP< FancyOStreamgetDefaultOStream ()
 
- Static Public Attributes inherited from Teuchos::Describable
static const EVerbosityLevel verbLevel_default
 

Detailed Description

template<class Scalar>
class Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >

Time integrator suitable for pseudotransient forward sensitivity analysis.

For some problems, time integrators are used to compute steady-state solutions (also known as pseudo-transient solvers). When computing sensitivities, it is not necessary in these cases to propagate sensitivities all the way through the forward time integration. Instead the steady-state is first computed as usual, and then the sensitivities are computed using a similar pseudo-transient time integration applied to the sensitivity equations with the state frozen to the computed steady-state. This integrator specializes the transient sensitivity methods implemented by Tempus::IntegratorForwardSensitivity to this case.

Consider an implicit ODE f(x_dot,x,p) = 0 with a stable steady-state solution x = x^s, x_dot = 0 where f(0,x^s,p) = 0 and all of the eigenvalues of df/dx(0,x^s,p) are in the right half-plane (for an explicit ODE, the eigenvalues must be in the left half-plane). In the pseudo-transient method a time-integrator is applied to f(x_dot,x,p) = 0 until x_dot is sufficiently small. Now consider the forward sensitivity equations: df/dx_dot*z_dot + df/dx*z + df/dp = 0 where z = dx/dp. For pseudo-transient forward sensitivities, the above is integrated from z(0) = 0 until z_dot is sufficiently small, in which case z^s = -(df/dx)^{-1}*(df/dp). Then the final sensitivity of g is dg/dp + dg/dx*z^s. One can see that z^s is the only steady-state solution of the sensitivity equations, since df/dx and df/dp are constant, and must be linearly stable since it has the same Jacobian matrix as the forward equations.

One should use the getX() and getDxDp() methods for extracting the final sultion and its parameter sensitivity as a multi-vector. This data can also be extracted from the solution history, but is stored as a Thyra product vector which requires knowledge of the internal implementation.

Definition at line 57 of file Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp.

Constructor & Destructor Documentation

template<class Scalar >
Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::IntegratorPseudoTransientForwardSensitivity ( const Teuchos::RCP< Thyra::ModelEvaluator< Scalar >> &  model,
const Teuchos::RCP< SensitivityModelEvaluatorBase< Scalar >> &  sens_model,
const Teuchos::RCP< IntegratorBasic< Scalar >> &  fwd_integrator,
const Teuchos::RCP< IntegratorBasic< Scalar >> &  sens_integrator,
const bool  reuse_solver,
const bool  force_W_update 
)

Constructor with ParameterList and model, and will be fully initialized.

In addition to all of the regular integrator options, the supplied parameter list supports the following options contained within a sublist "Sensitivities" from the top-level parameter list:

  • "Reuse State Linear Solver" (default: false) Whether to reuse the model's W matrix, solver, and preconditioner when solving the sensitivity equations. If they can be reused, substantial savings in compute time are possible.
  • "Force W Update" (default: false) When reusing the solver as above whether to force recomputation of W. This can be necessary when the solver overwrites it during the solve-phase (e.g., by a factorization).
  • "Use DfDp as Tangent" (default: false) Reinterpret the df/dp out-arg as the tangent vector (df/dx)(x,p) * dx/dp + df/dp(x,p) as described in the Tempus::CombinedForwardSensitivityModelEvaluator documentation.
  • "Sensitivity Parameter Index" (default: 0) Model evaluator parameter index for which sensitivities will be computed.
  • "Sensitivity X Tangent Index" (default: 1) If "Use DfDp as Tangent" is true, the model evaluator parameter index for passing dx/dp as a Thyra::DefaultMultiVectorProductVector.
  • "Sensitivity X-Dot Tangent Index" (default: 2) If "Use DfDp as Tangent" is true, the model evaluator parameter index for passing dx_dot/dp as a Thyra::DefaultMultiVectorProductVector.
  • "Sensitivity X-Dot-Dot Tangent Index" (default: 3) If "Use DfDp as Tangent" is true, the model evaluator parameter index for passing dx_dot_dot/dp as a Thyra::DefaultMultiVectorProductVector (if the model supports x_dot_dot).

Definition at line 22 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

Destructor.

Constructor that requires a subsequent setStepper, and initialize calls.

Definition at line 40 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
virtual Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::~IntegratorPseudoTransientForwardSensitivity ( )
inlinevirtual

Destructor.

Definition at line 107 of file Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp.

Member Function Documentation

template<class Scalar >
bool Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::advanceTime ( )
virtual

Advance the solution to timeMax, and return true if successful.

Definition at line 50 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
bool Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::advanceTime ( const Scalar  timeFinal)
overridevirtual

Advance the solution to timeFinal, and return true if successful.

Implements Tempus::Integrator< Scalar >.

Definition at line 79 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
Scalar Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getTime ( ) const
overridevirtual

Get current time.

Implements Tempus::Integrator< Scalar >.

Definition at line 125 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
int Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getIndex ( ) const
overridevirtual

Get current index.

Implements Tempus::Integrator< Scalar >.

Definition at line 131 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
Status Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getStatus ( ) const
overridevirtual
template<class Scalar >
void Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::setStatus ( const Status  st)
overridevirtual
template<class Scalar >
Teuchos::RCP< Stepper< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getStepper ( ) const
overridevirtual
template<class Scalar >
Teuchos::RCP< Stepper< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getStateStepper ( ) const
template<class Scalar >
Teuchos::RCP< Stepper< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getSensStepper ( ) const
template<class Scalar >
Teuchos::RCP< const SolutionHistory< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getSolutionHistory ( ) const
overridevirtual
template<class Scalar >
Teuchos::RCP< const SolutionHistory< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getStateSolutionHistory ( ) const
template<class Scalar >
Teuchos::RCP< const SolutionHistory< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getSensSolutionHistory ( ) const
template<class Scalar >
Teuchos::RCP< SolutionHistory< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getNonConstSolutionHistory ( )
overridevirtual
template<class Scalar >
Teuchos::RCP< const TimeStepControl< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getTimeStepControl ( ) const
overridevirtual
template<class Scalar >
Teuchos::RCP< TimeStepControl< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getNonConstTimeStepControl ( )
overridevirtual
template<class Scalar >
Teuchos::RCP< TimeStepControl< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getStateNonConstTimeStepControl ( )
template<class Scalar >
Teuchos::RCP< TimeStepControl< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getSensNonConstTimeStepControl ( )
template<class Scalar >
Teuchos::RCP< IntegratorObserver< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getObserver ( )
virtual

Get the Observer.

Definition at line 239 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
void Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::setObserver ( Teuchos::RCP< IntegratorObserver< Scalar >>  obs = Teuchos::null)
virtual

Set the Observer.

Definition at line 245 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
virtual Teuchos::RCP<Teuchos::Time> Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getIntegratorTimer ( ) const
inlineoverridevirtual

Returns the IntegratorTimer_ for this Integrator.

Implements Tempus::Integrator< Scalar >.

Definition at line 148 of file Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp.

template<class Scalar >
virtual Teuchos::RCP<Teuchos::Time> Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getStepperTimer ( ) const
inlineoverridevirtual
template<class Scalar >
void Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::initializeSolutionHistory ( Scalar  t0,
Teuchos::RCP< const Thyra::VectorBase< Scalar >>  x0,
Teuchos::RCP< const Thyra::VectorBase< Scalar >>  xdot0 = Teuchos::null,
Teuchos::RCP< const Thyra::VectorBase< Scalar >>  xdotdot0 = Teuchos::null,
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar >>  DxDp0 = Teuchos::null,
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar >>  DxdotDp0 = Teuchos::null,
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar >>  DxdotdotDp0 = Teuchos::null 
)
virtual

Set the initial state from Thyra::VectorBase(s)

Definition at line 254 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< const Thyra::VectorBase< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getX ( ) const
virtual

Get current the solution, x.

Definition at line 302 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getDxDp ( ) const
virtual
template<class Scalar >
Teuchos::RCP< const Thyra::VectorBase< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getXDot ( ) const
virtual

Get current the time derivative of the solution, xdot.

Definition at line 321 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getDXDotDp ( ) const
virtual
template<class Scalar >
Teuchos::RCP< const Thyra::VectorBase< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getXDotDot ( ) const
virtual

Get current the second time derivative of the solution, xdotdot.

Definition at line 341 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getDXDotDotDp ( ) const
virtual
template<class Scalar >
Teuchos::RCP< const Thyra::VectorBase< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getG ( ) const
virtual

Return response function g.

Definition at line 361 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getDgDp ( ) const
virtual

Return forward sensitivity stored in Jacobian format.

Definition at line 386 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
std::string Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::description ( ) const
overridevirtual
template<class Scalar >
void Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::describe ( Teuchos::FancyOStream out,
const Teuchos::EVerbosityLevel  verbLevel 
) const
overridevirtual
template<class Scalar >
SensitivityStepMode Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::getStepMode ( ) const

What mode the current time integration step is in.

Definition at line 434 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.

template<class Scalar >
void Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::buildSolutionHistory ( )
protected

Member Data Documentation

template<class Scalar >
Teuchos::RCP<Thyra::ModelEvaluator<Scalar> > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::model_
protected
template<class Scalar >
Teuchos::RCP<SensitivityModelEvaluatorBase<Scalar> > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::sens_model_
protected
template<class Scalar >
Teuchos::RCP<IntegratorBasic<Scalar> > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::state_integrator_
protected
template<class Scalar >
Teuchos::RCP<IntegratorBasic<Scalar> > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::sens_integrator_
protected
template<class Scalar >
Teuchos::RCP<SolutionHistory<Scalar> > Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::solutionHistory_
protected
template<class Scalar >
bool Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::reuse_solver_
protected
template<class Scalar >
bool Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::force_W_update_
protected
template<class Scalar >
SensitivityStepMode Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::stepMode_
protected

The documentation for this class was generated from the following files: