| Tempus
    Version of the Day
    Time Integration | 
Time integrator suitable for pseudotransient adjoint sensitivity analysis. More...
#include <Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp>
 
  
 | Public Member Functions | |
| IntegratorPseudoTransientAdjointSensitivity (Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &adjoint_residual_model, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &adjoint_solve_model) | |
| Constructor with ParameterList and model, and will be fully initialized.  More... | |
| IntegratorPseudoTransientAdjointSensitivity (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &adjoint_residual_model, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &adjoint_solve_model, std::string stepperType) | |
| Constructor with model and "Stepper Type" and is fully initialized with default settings.  More... | |
| IntegratorPseudoTransientAdjointSensitivity (Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &adjoint_model) | |
| Version of the constructor taking a single adjoint model evaluator.  More... | |
| IntegratorPseudoTransientAdjointSensitivity (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &adjoint_model, std::string stepperType) | |
| Version of the constructor taking a single adjoint model evaluator.  More... | |
| IntegratorPseudoTransientAdjointSensitivity (Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model) | |
| Version of the constructor taking a single model evaluator.  More... | |
| IntegratorPseudoTransientAdjointSensitivity (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, std::string stepperType) | |
| Version of the constructor taking a single model evaluator.  More... | |
| IntegratorPseudoTransientAdjointSensitivity () | |
| Destructor.  More... | |
| virtual | ~IntegratorPseudoTransientAdjointSensitivity () | 
| 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 > > y0=Teuchos::null, Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > ydot0=Teuchos::null, Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > ydotdot0=Teuchos::null) | 
| Set the initial state from Thyra::VectorBase(s)  More... | |
| virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > | getX () const | 
| Get the current solution, x.  More... | |
| virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > | getXDot () const | 
| Get the current time derivative of the solution, xdot.  More... | |
| virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > | getXDotDot () const | 
| Get the current second time derivative of the solution, xdotdot.  More... | |
| virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > | getY () const | 
| Get the current adjoint solution, y.  More... | |
| virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > | getYDot () const | 
| Get the current time derivative of the adjoint solution, ydot.  More... | |
| virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > | getYDotDot () const | 
| Get the current second time derivative of the adjoint solution, ydotdot.  More... | |
| virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > | getG () const | 
| Return response function g.  More... | |
| virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > | getDgDp () const | 
| Return adjoint sensitivity stored in gradient format.  More... | |
| SensitivityStepMode | getStepMode () const | 
| What mode the current time integration step is in.  More... | |
| void | setDoForwardIntegration (const bool f) | 
| Set/get whether to do the forward integration.  More... | |
| bool | getDoForwardIntegration () const | 
| void | setDoAdjointIntegration (const bool f) | 
| Set/get whether to do the adjoint integration.  More... | |
| bool | getDoAdjointIntegration () const | 
|  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 VerboseObject & | setVerbLevel (const EVerbosityLevel verbLevel) const | 
| virtual const VerboseObject & | setOverridingVerbLevel (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 VerboseObjectBase & | setOStream (const RCP< FancyOStream > &oStream) const | 
| virtual const VerboseObjectBase & | setOverridingOStream (const RCP< FancyOStream > &oStream) const | 
| virtual VerboseObjectBase & | setLinePrefix (const std::string &linePrefix) | 
| virtual RCP< FancyOStream > | getOStream () const | 
| virtual RCP< FancyOStream > | getOverridingOStream () const | 
| virtual std::string | getLinePrefix () const | 
| virtual OSTab | getOSTab (const int tabs=1, const std::string &linePrefix="") const | 
|  Public Member Functions inherited from Teuchos::ParameterListAcceptor | |
| virtual RCP< const ParameterList > | getParameterList () const | 
| Protected Types | |
| typedef Thyra::DefaultMultiVectorProductVector < Scalar > | DMVPV | 
| Protected Member Functions | |
| Teuchos::RCP < AdjointSensitivityModelEvaluator < Scalar > > | createSensitivityModel (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &adjoint_residual_model, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &adjoint_solve_model, const Teuchos::RCP< Teuchos::ParameterList > &inputPL) | 
| 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 < Thyra::ModelEvaluator < Scalar > > | adjoint_residual_model_ | 
| Teuchos::RCP < Thyra::ModelEvaluator < Scalar > > | adjoint_solve_model_ | 
| Teuchos::RCP < AdjointSensitivityModelEvaluator < Scalar > > | sens_model_ | 
| Teuchos::RCP< IntegratorBasic < Scalar > > | state_integrator_ | 
| Teuchos::RCP< IntegratorBasic < Scalar > > | sens_integrator_ | 
| Teuchos::RCP< SolutionHistory < Scalar > > | solutionHistory_ | 
| Teuchos::RCP < Thyra::VectorBase< Scalar > > | g_ | 
| Teuchos::RCP< DMVPV > | dgdp_ | 
| SensitivityStepMode | stepMode_ | 
| bool | do_forward_integration_ | 
| bool | do_adjoint_integration_ | 
| 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::ParameterListAcceptor | |
| void | setParameterList (const Teuchos::RCP< Teuchos::ParameterList > &pl) override | 
| Teuchos::RCP < Teuchos::ParameterList > | getNonconstParameterList () override | 
| Teuchos::RCP < Teuchos::ParameterList > | unsetParameterList () override | 
| Teuchos::RCP< const Teuchos::ParameterList > | getValidParameters () 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< FancyOStream > | getDefaultOStream () | 
|  Static Public Attributes inherited from Teuchos::Describable | |
| static const EVerbosityLevel | verbLevel_default | 
Time integrator suitable for pseudotransient adjoint 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 adjoint sensitivity equations with the state frozen to the computed steady-state. This integrator specializes the transient sensitivity methods implemented by Tempus::IntegratorAdjointSensitivity 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 adjoint sensitivity equations for some response function g(x,p): df/dx_dot^T*y_dot + df/dx^T*y - dg/dx^T = 0 after the transformation tau = T - t has been applied, where T is the final time. For pseudo-transient adjoint sensitivities, the above is integrated from y(0) = 0 until y_dot is sufficiently small, in which case y^s = (df/dx)^{-T}*(dg/dx)^T. Then the final sensitivity of g is dg/dp^T - df/dp^T*y^s. One can see that y^s is the only steady-state solution of the adjoint equations, since df/dx and dg/dx are constant, and must be linearly stable (since the eigenvalues of df/dx^T are the same as df/dx).
To extract the final solution x(T) and sensitivity dg/dp one should use the getX() and getDgDp() methods, which return these quantities directly. One can also extract this data for all times from the solution history, however the data is stored in Thyra product vectors which requires knowledge of the internal implementation.
Definition at line 59 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | protected | 
Definition at line 261 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::IntegratorPseudoTransientAdjointSensitivity | ( | Teuchos::RCP< Teuchos::ParameterList > | pList, | 
| const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > & | model, | ||
| const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > & | adjoint_residual_model, | ||
| const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > & | adjoint_solve_model | ||
| ) | 
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:
To support use-cases with explicitly computed adjoint operators, the constructor takes two additional model evaluators for computing the adjoint W/W_op. It is assumed the operator returned by these model evaluators is the adjoint, and so will not be transposed. It is also assumed these model evaluators accept the same inArgs as the forward model, however they only requires supporting the adjoint W/W_op outArgs. The first adjoint model evaluator will only be used for computing the adjoint residual, whereas the second will produce the operator used in linear solves. The passed RCP's may point to the same model evalutor object in cases where these operators are the same
| Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::IntegratorPseudoTransientAdjointSensitivity | ( | const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > & | model, | 
| const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > & | adjoint_residual_model, | ||
| const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > & | adjoint_solve_model, | ||
| std::string | stepperType | ||
| ) | 
Constructor with model and "Stepper Type" and is fully initialized with default settings.
| Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::IntegratorPseudoTransientAdjointSensitivity | ( | Teuchos::RCP< Teuchos::ParameterList > | pList, | 
| const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > & | model, | ||
| const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > & | adjoint_model | ||
| ) | 
Version of the constructor taking a single adjoint model evaluator.
| Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::IntegratorPseudoTransientAdjointSensitivity | ( | const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > & | model, | 
| const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > & | adjoint_model, | ||
| std::string | stepperType | ||
| ) | 
Version of the constructor taking a single adjoint model evaluator.
| Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::IntegratorPseudoTransientAdjointSensitivity | ( | Teuchos::RCP< Teuchos::ParameterList > | pList, | 
| const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > & | model | ||
| ) | 
Version of the constructor taking a single model evaluator.
This version takes a single model evaluator for the case when the adjoint is implicitly determined from the forward operator by the (conjugate) transpose.
| Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::IntegratorPseudoTransientAdjointSensitivity | ( | const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > & | model, | 
| std::string | stepperType | ||
| ) | 
Version of the constructor taking a single model evaluator.
This version takes a single model evaluator for the case when the adjoint is implicitly determined from the forward operator by the (conjugate) transpose.
| Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::IntegratorPseudoTransientAdjointSensitivity | ( | ) | 
Destructor.
Constructor that requires a subsequent setParameterList, setStepper, and initialize calls.
Definition at line 106 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | inlinevirtual | 
Destructor.
Definition at line 148 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | virtual | 
Advance the solution to timeMax, and return true if successful.
Definition at line 114 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Advance the solution to timeFinal, and return true if successful.
Implements Tempus::Integrator< Scalar >.
Definition at line 121 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Get current time.
Implements Tempus::Integrator< Scalar >.
Definition at line 193 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Get current index.
Implements Tempus::Integrator< Scalar >.
Definition at line 199 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Get Status.
Implements Tempus::Integrator< Scalar >.
Definition at line 205 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Set Status.
Implements Tempus::Integrator< Scalar >.
Definition at line 215 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Get the Stepper.
Implements Tempus::Integrator< Scalar >.
Definition at line 224 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| Teuchos::RCP< Stepper< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getStateStepper | ( | ) | const | 
Definition at line 231 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| Teuchos::RCP< Stepper< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getSensStepper | ( | ) | const | 
Definition at line 238 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Get the SolutionHistory.
Implements Tempus::Integrator< Scalar >.
Definition at line 245 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| Teuchos::RCP< const SolutionHistory< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getStateSolutionHistory | ( | ) | const | 
Definition at line 252 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| Teuchos::RCP< const SolutionHistory< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getSensSolutionHistory | ( | ) | const | 
Definition at line 260 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Get the SolutionHistory.
Implements Tempus::Integrator< Scalar >.
Definition at line 269 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Get the TimeStepControl.
Implements Tempus::Integrator< Scalar >.
Definition at line 276 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Implements Tempus::Integrator< Scalar >.
Definition at line 284 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| Teuchos::RCP< TimeStepControl< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getStateNonConstTimeStepControl | ( | ) | 
Definition at line 292 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| Teuchos::RCP< TimeStepControl< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getSensNonConstTimeStepControl | ( | ) | 
Definition at line 300 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | virtual | 
Get the Observer.
Definition at line 307 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | virtual | 
Set the Observer.
Definition at line 313 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | inlineoverridevirtual | 
Returns the IntegratorTimer_ for this Integrator.
Implements Tempus::Integrator< Scalar >.
Definition at line 190 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | inlineoverridevirtual | 
Implements Tempus::Integrator< Scalar >.
Definition at line 194 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | virtual | 
Set the initial state from Thyra::VectorBase(s)
Definition at line 322 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | virtual | 
Get the current solution, x.
Definition at line 369 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | virtual | 
Get the current time derivative of the solution, xdot.
Definition at line 376 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | virtual | 
Get the current second time derivative of the solution, xdotdot.
Definition at line 383 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | virtual | 
Get the current adjoint solution, y.
Definition at line 390 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | virtual | 
Get the current time derivative of the adjoint solution, ydot.
Definition at line 401 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | virtual | 
Get the current second time derivative of the adjoint solution, ydotdot.
Definition at line 412 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | virtual | 
Return response function g.
Definition at line 423 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | virtual | 
Return adjoint sensitivity stored in gradient format.
Definition at line 430 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Implements Teuchos::ParameterListAcceptor.
Definition at line 457 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Implements Teuchos::ParameterListAcceptor.
Definition at line 514 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Implements Teuchos::ParameterListAcceptor.
Definition at line 476 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Reimplemented from Teuchos::ParameterListAcceptor.
Definition at line 498 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Reimplemented from Teuchos::Describable.
Definition at line 436 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | overridevirtual | 
Reimplemented from Teuchos::Describable.
Definition at line 444 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| SensitivityStepMode Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getStepMode | ( | ) | const | 
What mode the current time integration step is in.
Definition at line 523 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | inline | 
Set/get whether to do the forward integration.
Definition at line 253 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | inline | 
Definition at line 254 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | inline | 
Set/get whether to do the adjoint integration.
Definition at line 257 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | inline | 
Definition at line 258 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | protected | 
Definition at line 530 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | protected | 
Definition at line 550 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
| 
 | protected | 
Definition at line 274 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | protected | 
Definition at line 275 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | protected | 
Definition at line 276 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | protected | 
Definition at line 277 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | protected | 
Definition at line 278 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | protected | 
Definition at line 279 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | protected | 
Definition at line 280 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | protected | 
Definition at line 281 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | protected | 
Definition at line 282 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | protected | 
Definition at line 283 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | protected | 
Definition at line 284 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
| 
 | protected | 
Definition at line 285 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
 1.8.5
 1.8.5