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::IntegratorAdjointSensitivity< Scalar > Class Template Reference

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

#include <Tempus_IntegratorAdjointSensitivity_decl.hpp>

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

Public Member Functions

 IntegratorAdjointSensitivity (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar >> &model, const Teuchos::RCP< IntegratorBasic< Scalar >> &state_integrator, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar >> &adjoint_model, const Teuchos::RCP< AdjointAuxSensitivityModelEvaluator< Scalar >> &adjoint_aux_model, const Teuchos::RCP< IntegratorBasic< Scalar >> &adjoint_integrator, const Teuchos::RCP< SolutionHistory< Scalar >> &solution_history, const int p_index, const int g_index, const bool g_depends_on_p, const bool f_depends_on_p, const bool ic_depends_on_p, const bool mass_matrix_is_identity)
 Full Constructor with model, and will be fully initialized. More...
 
 IntegratorAdjointSensitivity ()
 Constructor that requires a subsequent setParameterList, setStepper, and initialize calls. More...
 
virtual ~IntegratorAdjointSensitivity ()
 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
< IntegratorObserver< Scalar > > 
getObserver ()
 Get the Observer. More...
 
virtual void setObserver (Teuchos::RCP< IntegratorObserver< Scalar >> obs=Teuchos::null)
 Set the Observer. More...
 
virtual void initialize ()
 Initializes the Integrator after set* function calls. 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::MultiVectorBase< Scalar > > 
getDgDp () const
 Return adjoint sensitivity stored in gradient format. More...
 
Teuchos::RCP
< Thyra::ModelEvaluator
< Scalar > > 
getAdjointModel () const
 
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

Teuchos::RCP
< AdjointAuxSensitivityModelEvaluator
< Scalar > > 
createAdjointModel (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar >> &model, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar >> &adjoint_model, const Teuchos::RCP< Teuchos::ParameterList > &inputPL)
 
void buildSolutionHistory (const Teuchos::RCP< const SolutionHistory< Scalar >> &state_solution_history, const Teuchos::RCP< const SolutionHistory< Scalar >> &adjoint_solution_history)
 
- 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< IntegratorBasic
< Scalar > > 
state_integrator_
 
Teuchos::RCP
< Thyra::ModelEvaluator
< Scalar > > 
adjoint_model_
 
Teuchos::RCP
< AdjointAuxSensitivityModelEvaluator
< Scalar > > 
adjoint_aux_model_
 
Teuchos::RCP< IntegratorBasic
< Scalar > > 
adjoint_integrator_
 
Teuchos::RCP< SolutionHistory
< Scalar > > 
solutionHistory_
 
int p_index_
 
int g_index_
 
bool g_depends_on_p_
 
bool f_depends_on_p_
 
bool ic_depends_on_p_
 
bool mass_matrix_is_identity_
 
Teuchos::RCP< const
Thyra::MultiVectorBase< Scalar > > 
dxdp_init_
 
Teuchos::RCP
< Thyra::MultiVectorBase
< Scalar > > 
dgdp_
 
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...
 
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
< 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::IntegratorAdjointSensitivity< Scalar >

Time integrator suitable for adjoint sensitivity analysis.

This integrator implements transient adjoint sensitivities. Given a model evaluator encapsulating the equations f(x_dot,x,p) = 0, and a response function g(x,p), these equations are integrated forward in time to some final time T. Then the adjoint equations F(y) = d/dt( df/dx_dot^T*y ) - df/dx^T*y = 0 are integrated backward in time starting at t = T where y is the adjoint variable. Nominally y is a multi-vector belonging to the vector space of f and the number of columns given by the number of entries in the response g. The initial conditions for y at t = T are given by y(T) = df/dx_dot(x_dot(T),x(T),p)^{-T} * dg/dx(x(T),p)^T. Then the final sensitivity of g is dg/dp(T) - int_{0}^T(df/dp^T*y)dt + dx/dp(0)^T*df/dx_dot(0)^T*y(0).

This integrator supports both implicit and explicit steppers, and provides a method to compute dg/dp as described above after the reverse integration. The solution history contains solutions stored as product vectors (x,y), with y further stored as a product vector for each component of g.

Because of limitations on the steppers, the implementation currently assumes df/dxdot is a constant matrix.

To extract the final solution x(T) and sensitivity dg/dp(T) 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 51 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.

Constructor & Destructor Documentation

template<class Scalar >
Tempus::IntegratorAdjointSensitivity< Scalar >::IntegratorAdjointSensitivity ( const Teuchos::RCP< Thyra::ModelEvaluator< Scalar >> &  model,
const Teuchos::RCP< IntegratorBasic< Scalar >> &  state_integrator,
const Teuchos::RCP< Thyra::ModelEvaluator< Scalar >> &  adjoint_model,
const Teuchos::RCP< AdjointAuxSensitivityModelEvaluator< Scalar >> &  adjoint_aux_model,
const Teuchos::RCP< IntegratorBasic< Scalar >> &  adjoint_integrator,
const Teuchos::RCP< SolutionHistory< Scalar >> &  solution_history,
const int  p_index,
const int  g_index,
const bool  g_depends_on_p,
const bool  f_depends_on_p,
const bool  ic_depends_on_p,
const bool  mass_matrix_is_identity 
)

Full Constructor with model, and will be fully initialized.

Parameters
[in]modelThe forward physics ModelEvaluator
[in]state_integratorForward state Integrator for the forward problem
[in]adjoint_modelModelEvaluator for the adjoint physics/problem
[in]adjoint_aux_modelModelEvaluator for the auxiliary adjoint physics/problem
[in]adjoint_integratorTime integrator for the adjoint problem
[in]solution_historyThe forward state solution history
[in]p_indexSensitivity parameter index
[in]g_indexResponse function index
[in]g_depends_on_pDoes response depends on parameters?
[in]f_depends_on_pDoes residual depends on parameters?
[in]ic_depends_on_pDoes the initial condition depends on parameters?
adjoint_modelModelEvaluator for the adjoint problem. Optional. Default value is null.
[in]mass_matrix_is_identityIs the mass matrix an identity matrix?

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.
  • "Response Depends on Parameters", (default: true) Whether the response function depends on the parameter vector p for which sensitivities will be computed. If set to false, the dg/dp term in the sensitivity formula will not be computed.
  • "Residual Depends on Parameters", (default: true) Whether the model residual f depends on the parameter vector p for which sensitivities will be computed. If set to false, its contribution to the sensitivities will not be computed.
  • "IC Depends on Parameters", (default: true) Whether the initial conditions depend on the parameter vector p for which sensitivities will be computed. If set to false, its contribution to the sensitivities will not 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 an additional model evaluator for computing the adjoint W/W_op. It is assumed the operator returned by this model evaluator is the adjoint, and so will not be transposed. It is also assumed this model evaluator accepts the same inArgs as the forward model, however it only requires supporting the adjoint W/W_op outArgs.

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 25 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

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

Definition at line 61 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

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

Destructor.

Definition at line 129 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.

Member Function Documentation

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

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

Definition at line 69 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

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

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

Implements Tempus::Integrator< Scalar >.

Definition at line 76 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

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

Get current time.

Implements Tempus::Integrator< Scalar >.

Definition at line 307 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

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

Get current index.

Implements Tempus::Integrator< Scalar >.

Definition at line 313 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

template<class Scalar >
Status Tempus::IntegratorAdjointSensitivity< Scalar >::getStatus ( ) const
overridevirtual

Get Status.

Implements Tempus::Integrator< Scalar >.

Definition at line 319 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

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

Set Status.

Implements Tempus::Integrator< Scalar >.

Definition at line 329 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< Stepper< Scalar > > Tempus::IntegratorAdjointSensitivity< Scalar >::getStepper ( ) const
overridevirtual
template<class Scalar >
Teuchos::RCP< const SolutionHistory< Scalar > > Tempus::IntegratorAdjointSensitivity< Scalar >::getSolutionHistory ( ) const
overridevirtual
template<class Scalar >
Teuchos::RCP< const SolutionHistory< Scalar > > Tempus::IntegratorAdjointSensitivity< Scalar >::getStateSolutionHistory ( ) const
template<class Scalar >
Teuchos::RCP< const SolutionHistory< Scalar > > Tempus::IntegratorAdjointSensitivity< Scalar >::getSensSolutionHistory ( ) const
template<class Scalar >
Teuchos::RCP< SolutionHistory< Scalar > > Tempus::IntegratorAdjointSensitivity< Scalar >::getNonConstSolutionHistory ( )
overridevirtual
template<class Scalar >
Teuchos::RCP< const TimeStepControl< Scalar > > Tempus::IntegratorAdjointSensitivity< Scalar >::getTimeStepControl ( ) const
overridevirtual
template<class Scalar >
Teuchos::RCP< TimeStepControl< Scalar > > Tempus::IntegratorAdjointSensitivity< Scalar >::getNonConstTimeStepControl ( )
overridevirtual
template<class Scalar >
Teuchos::RCP< TimeStepControl< Scalar > > Tempus::IntegratorAdjointSensitivity< Scalar >::getStateNonConstTimeStepControl ( )
template<class Scalar >
Teuchos::RCP< TimeStepControl< Scalar > > Tempus::IntegratorAdjointSensitivity< Scalar >::getSensNonConstTimeStepControl ( )
template<class Scalar >
virtual Teuchos::RCP<Teuchos::Time> Tempus::IntegratorAdjointSensitivity< Scalar >::getIntegratorTimer ( ) const
inlineoverridevirtual

Returns the IntegratorTimer_ for this Integrator.

Implements Tempus::Integrator< Scalar >.

Definition at line 164 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.

template<class Scalar >
virtual Teuchos::RCP<Teuchos::Time> Tempus::IntegratorAdjointSensitivity< Scalar >::getStepperTimer ( ) const
inlineoverridevirtual
template<class Scalar >
void Tempus::IntegratorAdjointSensitivity< 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 399 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< IntegratorObserver< Scalar > > Tempus::IntegratorAdjointSensitivity< Scalar >::getObserver ( )
virtual

Get the Observer.

Definition at line 413 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

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

Set the Observer.

Definition at line 419 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

template<class Scalar >
void Tempus::IntegratorAdjointSensitivity< Scalar >::initialize ( )
virtual

Initializes the Integrator after set* function calls.

Definition at line 429 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

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

Get the current solution, x.

Definition at line 437 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

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

Get the current time derivative of the solution, xdot.

Definition at line 444 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

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

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

Definition at line 451 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > Tempus::IntegratorAdjointSensitivity< Scalar >::getY ( ) const
virtual

Get the current adjoint solution, y.

Definition at line 458 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > Tempus::IntegratorAdjointSensitivity< Scalar >::getYDot ( ) const
virtual

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

Definition at line 471 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > Tempus::IntegratorAdjointSensitivity< Scalar >::getYDotDot ( ) const
virtual

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

Definition at line 485 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

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

Return adjoint sensitivity stored in gradient format.

Definition at line 499 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP<Thyra::ModelEvaluator<Scalar> > Tempus::IntegratorAdjointSensitivity< Scalar >::getAdjointModel ( ) const
inline
template<class Scalar >
std::string Tempus::IntegratorAdjointSensitivity< Scalar >::description ( ) const
overridevirtual

Reimplemented from Teuchos::Describable.

Definition at line 505 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

template<class Scalar >
void Tempus::IntegratorAdjointSensitivity< Scalar >::describe ( Teuchos::FancyOStream out,
const Teuchos::EVerbosityLevel  verbLevel 
) const
overridevirtual

Reimplemented from Teuchos::Describable.

Definition at line 512 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

template<class Scalar >
SensitivityStepMode Tempus::IntegratorAdjointSensitivity< Scalar >::getStepMode ( ) const

What mode the current time integration step is in.

Definition at line 525 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.

template<class Scalar >
Teuchos::RCP< AdjointAuxSensitivityModelEvaluator< Scalar > > Tempus::IntegratorAdjointSensitivity< Scalar >::createAdjointModel ( const Teuchos::RCP< Thyra::ModelEvaluator< Scalar >> &  model,
const Teuchos::RCP< Thyra::ModelEvaluator< Scalar >> &  adjoint_model,
const Teuchos::RCP< Teuchos::ParameterList > &  inputPL 
)
protected
template<class Scalar >
void Tempus::IntegratorAdjointSensitivity< Scalar >::buildSolutionHistory ( const Teuchos::RCP< const SolutionHistory< Scalar >> &  state_solution_history,
const Teuchos::RCP< const SolutionHistory< Scalar >> &  adjoint_solution_history 
)
protected

Member Data Documentation

template<class Scalar >
Teuchos::RCP<Thyra::ModelEvaluator<Scalar> > Tempus::IntegratorAdjointSensitivity< Scalar >::model_
protected
template<class Scalar >
Teuchos::RCP<IntegratorBasic<Scalar> > Tempus::IntegratorAdjointSensitivity< Scalar >::state_integrator_
protected
template<class Scalar >
Teuchos::RCP<Thyra::ModelEvaluator<Scalar> > Tempus::IntegratorAdjointSensitivity< Scalar >::adjoint_model_
protected
template<class Scalar >
Teuchos::RCP<AdjointAuxSensitivityModelEvaluator<Scalar> > Tempus::IntegratorAdjointSensitivity< Scalar >::adjoint_aux_model_
protected
template<class Scalar >
Teuchos::RCP<IntegratorBasic<Scalar> > Tempus::IntegratorAdjointSensitivity< Scalar >::adjoint_integrator_
protected
template<class Scalar >
Teuchos::RCP<SolutionHistory<Scalar> > Tempus::IntegratorAdjointSensitivity< Scalar >::solutionHistory_
protected
template<class Scalar >
int Tempus::IntegratorAdjointSensitivity< Scalar >::p_index_
protected
template<class Scalar >
int Tempus::IntegratorAdjointSensitivity< Scalar >::g_index_
protected
template<class Scalar >
bool Tempus::IntegratorAdjointSensitivity< Scalar >::g_depends_on_p_
protected
template<class Scalar >
bool Tempus::IntegratorAdjointSensitivity< Scalar >::f_depends_on_p_
protected
template<class Scalar >
bool Tempus::IntegratorAdjointSensitivity< Scalar >::ic_depends_on_p_
protected
template<class Scalar >
bool Tempus::IntegratorAdjointSensitivity< Scalar >::mass_matrix_is_identity_
protected
template<class Scalar >
Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > Tempus::IntegratorAdjointSensitivity< Scalar >::dxdp_init_
protected
template<class Scalar >
Teuchos::RCP<Thyra::MultiVectorBase<Scalar> > Tempus::IntegratorAdjointSensitivity< Scalar >::dgdp_
protected
template<class Scalar >
SensitivityStepMode Tempus::IntegratorAdjointSensitivity< Scalar >::stepMode_
protected

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