Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Public Member Functions | List of all members
Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar > Class Template Reference

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

#include <Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp>

Inheritance diagram for Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >:
Tempus::Integrator< Scalar >

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.
 
 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.
 
 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.
 
 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.
 
 IntegratorPseudoTransientAdjointSensitivity (Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model)
 Version of the constructor taking a single model evaluator.
 
 IntegratorPseudoTransientAdjointSensitivity (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, std::string stepperType)
 Version of the constructor taking a single model evaluator.
 
 IntegratorPseudoTransientAdjointSensitivity ()
 Destructor.
 
virtual ~IntegratorPseudoTransientAdjointSensitivity ()
 Destructor.
 
- Public Member Functions inherited from Tempus::Integrator< Scalar >

Overridden from Teuchos::Describable

typedef Thyra::DefaultMultiVectorProductVector< Scalar > DMVPV
 
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< DMVPVdgdp_
 
SensitivityStepMode stepMode_
 
bool do_forward_integration_
 
bool do_adjoint_integration_
 
std::string description () const override
 
void describe (Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
 
SensitivityStepMode getStepMode () const
 What mode the current time integration step is in.
 
void setDoForwardIntegration (const bool f)
 Set/get whether to do the forward integration.
 
bool getDoForwardIntegration () const
 
void setDoAdjointIntegration (const bool f)
 Set/get whether to do the adjoint integration.
 
bool getDoAdjointIntegration () const
 
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 ()
 

Basic integrator methods

virtual bool advanceTime ()
 Advance the solution to timeMax, and return true if successful.
 
virtual bool advanceTime (const Scalar timeFinal) override
 Advance the solution to timeFinal, and return true if successful.
 
virtual Scalar getTime () const override
 Get current time.
 
virtual int getIndex () const override
 Get current index.
 
virtual Status getStatus () const override
 Get Status.
 
virtual void setStatus (const Status st) override
 Set Status.
 
virtual Teuchos::RCP< Stepper< Scalar > > getStepper () const override
 Get the Stepper.
 
Teuchos::RCP< Stepper< Scalar > > getStateStepper () const
 
Teuchos::RCP< Stepper< Scalar > > getSensStepper () const
 
virtual Teuchos::RCP< const SolutionHistory< Scalar > > getSolutionHistory () const override
 Get the SolutionHistory.
 
Teuchos::RCP< const SolutionHistory< Scalar > > getStateSolutionHistory () const
 
Teuchos::RCP< const SolutionHistory< Scalar > > getSensSolutionHistory () const
 
virtual Teuchos::RCP< SolutionHistory< Scalar > > getNonConstSolutionHistory () override
 Get the SolutionHistory.
 
virtual Teuchos::RCP< const TimeStepControl< Scalar > > getTimeStepControl () const override
 Get the TimeStepControl.
 
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.
 
virtual void setObserver (Teuchos::RCP< IntegratorObserver< Scalar > > obs=Teuchos::null)
 Set the Observer.
 
virtual Teuchos::RCP< Teuchos::Time > getIntegratorTimer () const override
 Returns the IntegratorTimer_ for this Integrator.
 
virtual Teuchos::RCP< Teuchos::Time > getStepperTimer () const override
 
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)
 
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getX () const
 Get the current solution, x.
 
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXDot () const
 Get the current time derivative of the solution, xdot.
 
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXDotDot () const
 Get the current second time derivative of the solution, xdotdot.
 
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getY () const
 Get the current adjoint solution, y.
 
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getYDot () const
 Get the current time derivative of the adjoint solution, ydot.
 
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getYDotDot () const
 Get the current second time derivative of the adjoint solution, ydotdot.
 
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getG () const
 Return response function g.
 
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDgDp () const
 Return adjoint sensitivity stored in gradient format.
 

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
 

Detailed Description

template<class Scalar>
class Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >

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 60 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.

Member Typedef Documentation

◆ DMVPV

template<class Scalar >
typedef Thyra::DefaultMultiVectorProductVector<Scalar> Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::DMVPV
protected

Constructor & Destructor Documentation

◆ IntegratorPseudoTransientAdjointSensitivity() [1/7]

template<class Scalar >
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:

  • "Sensitivity Parameter Index", (default: 0) The model evaluator parameter index for which sensitivities will be computed.
  • "Response Function Index", (default: 0) The model evaluator response index for which sensitivities will be computed.
  • "Mass Matrix Is Constant" (default: true) Whether the mass matrix df/dx_dot is a constant matrix. As describe above, this is currently required to be true.
  • "Mass Matrix Is Identity" (default: false) Whether the mass matrix is the identity matrix, in which some computations can be skipped.

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

Definition at line 21 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ IntegratorPseudoTransientAdjointSensitivity() [2/7]

template<class Scalar >
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.

Definition at line 41 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ IntegratorPseudoTransientAdjointSensitivity() [3/7]

template<class Scalar >
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.

Definition at line 61 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ IntegratorPseudoTransientAdjointSensitivity() [4/7]

template<class Scalar >
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.

Definition at line 72 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ IntegratorPseudoTransientAdjointSensitivity() [5/7]

template<class Scalar >
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.

Definition at line 83 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ IntegratorPseudoTransientAdjointSensitivity() [6/7]

template<class Scalar >
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.

Definition at line 93 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ IntegratorPseudoTransientAdjointSensitivity() [7/7]

template<class Scalar >
Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::IntegratorPseudoTransientAdjointSensitivity

Destructor.

Constructor that requires a subsequent setParameterList, setStepper, and initialize calls.

Definition at line 103 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ ~IntegratorPseudoTransientAdjointSensitivity()

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

Destructor.

Definition at line 147 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.

Member Function Documentation

◆ advanceTime() [1/2]

template<class Scalar >
bool Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::advanceTime
virtual

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

Definition at line 113 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ advanceTime() [2/2]

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

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

Implements Tempus::Integrator< Scalar >.

Definition at line 123 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ getTime()

template<class Scalar >
Scalar Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getTime
overridevirtual

Get current time.

Implements Tempus::Integrator< Scalar >.

Definition at line 188 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ getIndex()

template<class Scalar >
int Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getIndex
overridevirtual

Get current index.

Implements Tempus::Integrator< Scalar >.

Definition at line 196 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ getStatus()

template<class Scalar >
Status Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getStatus
overridevirtual

◆ setStatus()

template<class Scalar >
void Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::setStatus ( const Status  st)
overridevirtual

◆ getStepper()

template<class Scalar >
Teuchos::RCP< Stepper< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getStepper
overridevirtual

◆ getStateStepper()

template<class Scalar >
Teuchos::RCP< Stepper< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getStateStepper

◆ getSensStepper()

template<class Scalar >
Teuchos::RCP< Stepper< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getSensStepper

◆ getSolutionHistory()

template<class Scalar >
Teuchos::RCP< const SolutionHistory< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getSolutionHistory
overridevirtual

◆ getStateSolutionHistory()

template<class Scalar >
Teuchos::RCP< const SolutionHistory< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getStateSolutionHistory

◆ getSensSolutionHistory()

template<class Scalar >
Teuchos::RCP< const SolutionHistory< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getSensSolutionHistory

◆ getNonConstSolutionHistory()

template<class Scalar >
Teuchos::RCP< SolutionHistory< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getNonConstSolutionHistory
overridevirtual

◆ getTimeStepControl()

template<class Scalar >
Teuchos::RCP< const TimeStepControl< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getTimeStepControl
overridevirtual

◆ getNonConstTimeStepControl()

template<class Scalar >
Teuchos::RCP< TimeStepControl< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getNonConstTimeStepControl
overridevirtual

◆ getStateNonConstTimeStepControl()

template<class Scalar >
Teuchos::RCP< TimeStepControl< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getStateNonConstTimeStepControl

◆ getSensNonConstTimeStepControl()

template<class Scalar >
Teuchos::RCP< TimeStepControl< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getSensNonConstTimeStepControl

◆ getObserver()

template<class Scalar >
Teuchos::RCP< IntegratorObserver< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getObserver
virtual

Get the Observer.

Definition at line 313 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ setObserver()

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

Set the Observer.

Definition at line 321 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ getIntegratorTimer()

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

Returns the IntegratorTimer_ for this Integrator.

Implements Tempus::Integrator< Scalar >.

Definition at line 185 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.

◆ getStepperTimer()

template<class Scalar >
virtual Teuchos::RCP< Teuchos::Time > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getStepperTimer ( ) const
inlineoverridevirtual

◆ initializeSolutionHistory()

template<class Scalar >
void Tempus::IntegratorPseudoTransientAdjointSensitivity< 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 > >  y0 = Teuchos::null,
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > >  ydot0 = Teuchos::null,
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > >  ydotdot0 = Teuchos::null 
)
virtual

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

Definition at line 329 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ getX()

template<class Scalar >
Teuchos::RCP< const Thyra::VectorBase< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getX
virtual

Get the current solution, x.

Definition at line 377 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ getXDot()

template<class Scalar >
Teuchos::RCP< const Thyra::VectorBase< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getXDot
virtual

Get the current time derivative of the solution, xdot.

Definition at line 385 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ getXDotDot()

template<class Scalar >
Teuchos::RCP< const Thyra::VectorBase< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getXDotDot
virtual

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

Definition at line 393 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ getY()

template<class Scalar >
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getY
virtual

Get the current adjoint solution, y.

Definition at line 401 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ getYDot()

template<class Scalar >
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getYDot
virtual

Get the current time derivative of the adjoint solution, ydot.

Definition at line 413 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ getYDotDot()

template<class Scalar >
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getYDotDot
virtual

Get the current second time derivative of the adjoint solution, ydotdot.

Definition at line 425 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ getG()

template<class Scalar >
Teuchos::RCP< const Thyra::VectorBase< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getG
virtual

Return response function g.

Definition at line 437 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ getDgDp()

template<class Scalar >
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getDgDp
virtual

Return adjoint sensitivity stored in gradient format.

Definition at line 445 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ setParameterList()

template<class Scalar >
void Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::setParameterList ( const Teuchos::RCP< Teuchos::ParameterList > &  pl)
override

◆ getNonconstParameterList()

template<class Scalar >
Teuchos::RCP< Teuchos::ParameterList > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getNonconstParameterList
override

◆ unsetParameterList()

template<class Scalar >
Teuchos::RCP< Teuchos::ParameterList > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::unsetParameterList
override

◆ getValidParameters()

template<class Scalar >
Teuchos::RCP< const Teuchos::ParameterList > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getValidParameters
override

◆ description()

template<class Scalar >
std::string Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::description
override

◆ describe()

template<class Scalar >
void Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::describe ( Teuchos::FancyOStream &  out,
const Teuchos::EVerbosityLevel  verbLevel 
) const
override

◆ getStepMode()

template<class Scalar >
SensitivityStepMode Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getStepMode

What mode the current time integration step is in.

Definition at line 547 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.

◆ setDoForwardIntegration()

template<class Scalar >
void Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::setDoForwardIntegration ( const bool  f)
inline

Set/get whether to do the forward integration.

Definition at line 243 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.

◆ getDoForwardIntegration()

template<class Scalar >
bool Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getDoForwardIntegration ( ) const
inline

◆ setDoAdjointIntegration()

template<class Scalar >
void Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::setDoAdjointIntegration ( const bool  f)
inline

Set/get whether to do the adjoint integration.

Definition at line 247 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.

◆ getDoAdjointIntegration()

template<class Scalar >
bool Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::getDoAdjointIntegration ( ) const
inline

◆ createSensitivityModel()

template<class Scalar >
Teuchos::RCP< AdjointSensitivityModelEvaluator< Scalar > > Tempus::IntegratorPseudoTransientAdjointSensitivity< 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 
)
protected

◆ buildSolutionHistory()

template<class Scalar >
void Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::buildSolutionHistory
protected

Member Data Documentation

◆ model_

template<class Scalar >
Teuchos::RCP<Thyra::ModelEvaluator<Scalar> > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::model_
protected

◆ adjoint_residual_model_

template<class Scalar >
Teuchos::RCP<Thyra::ModelEvaluator<Scalar> > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::adjoint_residual_model_
protected

◆ adjoint_solve_model_

template<class Scalar >
Teuchos::RCP<Thyra::ModelEvaluator<Scalar> > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::adjoint_solve_model_
protected

◆ sens_model_

template<class Scalar >
Teuchos::RCP<AdjointSensitivityModelEvaluator<Scalar> > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::sens_model_
protected

◆ state_integrator_

template<class Scalar >
Teuchos::RCP<IntegratorBasic<Scalar> > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::state_integrator_
protected

◆ sens_integrator_

template<class Scalar >
Teuchos::RCP<IntegratorBasic<Scalar> > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::sens_integrator_
protected

◆ solutionHistory_

template<class Scalar >
Teuchos::RCP<SolutionHistory<Scalar> > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::solutionHistory_
protected

◆ g_

template<class Scalar >
Teuchos::RCP<Thyra::VectorBase<Scalar> > Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::g_
protected

◆ dgdp_

template<class Scalar >
Teuchos::RCP<DMVPV> Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::dgdp_
protected

◆ stepMode_

template<class Scalar >
SensitivityStepMode Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::stepMode_
protected

◆ do_forward_integration_

template<class Scalar >
bool Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::do_forward_integration_
protected

◆ do_adjoint_integration_

template<class Scalar >
bool Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::do_adjoint_integration_
protected

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