Tempus
Version of the Day
Time Integration
|
Time integrator suitable for adjoint sensitivity analysis. More...
#include <Tempus_IntegratorAdjointSensitivity_decl.hpp>
Public Member Functions | |
IntegratorAdjointSensitivity (Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model) | |
Constructor with ParameterList and model, and will be fully initialized. More... | |
IntegratorAdjointSensitivity () | |
Destructor. 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 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 current the solution, x. More... | |
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > | getXDot () const |
Get current the time derivative of the solution, xdot. More... | |
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 > > | getDgDp () const |
Return adjoint sensitivity stored in gradient format. 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 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 Member Functions | |
Teuchos::RCP < AdjointAuxSensitivityModelEvaluator < Scalar > > | createAdjointModel (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &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 < AdjointAuxSensitivityModelEvaluator < Scalar > > | adjoint_model_ |
Teuchos::RCP < IntegratorBasicOld< Scalar > > | state_integrator_ |
Teuchos::RCP < IntegratorBasicOld< 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_ |
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 Teuchos::RCP< Stepper < Scalar > > | getStepper () const override |
Get the Stepper. More... | |
virtual Teuchos::RCP < Teuchos::ParameterList > | getTempusParameterList () override |
Return a copy of the Tempus ParameterList. More... | |
virtual void | setTempusParameterList (Teuchos::RCP< Teuchos::ParameterList > pl) override |
virtual Teuchos::RCP< const SolutionHistory< Scalar > > | getSolutionHistory () const 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 |
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 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.
Definition at line 44 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
Tempus::IntegratorAdjointSensitivity< Scalar >::IntegratorAdjointSensitivity | ( | Teuchos::RCP< Teuchos::ParameterList > | pList, |
const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > & | 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:
Definition at line 24 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
Tempus::IntegratorAdjointSensitivity< Scalar >::IntegratorAdjointSensitivity | ( | ) |
Destructor.
Constructor that requires a subsequent setParameterList, setStepper, and initialize calls.
Definition at line 46 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
inlinevirtual |
Destructor.
Definition at line 89 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
|
virtual |
Advance the solution to timeMax, and return true if successful.
Definition at line 55 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Advance the solution to timeFinal, and return true if successful.
Implements Tempus::Integrator< Scalar >.
Definition at line 65 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Get current time.
Implements Tempus::Integrator< Scalar >.
Definition at line 223 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Get current index.
Implements Tempus::Integrator< Scalar >.
Definition at line 231 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Get Status.
Implements Tempus::Integrator< Scalar >.
Definition at line 239 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Get the Stepper.
Implements Tempus::Integrator< Scalar >.
Definition at line 253 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Return a copy of the Tempus ParameterList.
Implements Tempus::Integrator< Scalar >.
Definition at line 261 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Implements Tempus::Integrator< Scalar >.
Definition at line 269 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Get the SolutionHistory.
Implements Tempus::Integrator< Scalar >.
Definition at line 278 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Get the TimeStepControl.
Implements Tempus::Integrator< Scalar >.
Definition at line 286 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Implements Tempus::Integrator< Scalar >.
Definition at line 294 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
inlineoverridevirtual |
Returns the IntegratorTimer_ for this Integrator.
Implements Tempus::Integrator< Scalar >.
Definition at line 115 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
|
inlineoverridevirtual |
Implements Tempus::Integrator< Scalar >.
Definition at line 117 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
|
virtual |
Set the initial state from Thyra::VectorBase(s)
Definition at line 301 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
virtual |
Set the Observer.
Definition at line 315 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
virtual |
Initializes the Integrator after set* function calls.
Definition at line 325 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
virtual |
Get current the solution, x.
Definition at line 334 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
virtual |
Get current the time derivative of the solution, xdot.
Definition at line 342 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
virtual |
Get current the second time derivative of the solution, xdotdot.
Definition at line 350 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
virtual |
Return adjoint sensitivity stored in gradient format.
Definition at line 358 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Implements Teuchos::ParameterListAcceptor.
Definition at line 390 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Implements Teuchos::ParameterListAcceptor.
Definition at line 441 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Implements Teuchos::ParameterListAcceptor.
Definition at line 408 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Reimplemented from Teuchos::ParameterListAcceptor.
Definition at line 417 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Reimplemented from Teuchos::Describable.
Definition at line 366 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
overridevirtual |
Reimplemented from Teuchos::Describable.
Definition at line 375 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
protected |
Definition at line 449 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
protected |
Definition at line 472 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
|
protected |
Definition at line 178 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 179 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 180 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 181 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 182 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 183 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 184 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 185 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 186 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 187 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 188 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 189 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 190 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.