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) | |
Constructor with ParameterList and model, and will be fully initialized. More... | |
IntegratorPseudoTransientAdjointSensitivity (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, std::string stepperType) | |
Constructor with model and "Stepper Type" and is fully initialized with default settings. 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 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... | |
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 Scalar | 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 < 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 |
Basic integrator methods |
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< Teuchos::ParameterList > &inputPL) |
void | buildSolutionHistory () |
Protected Attributes | |
Teuchos::RCP < Thyra::ModelEvaluator < Scalar > > | 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< DMVPV > | dgdp_ |
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).
Definition at line 51 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 159 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::IntegratorPseudoTransientAdjointSensitivity | ( | 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 22 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::IntegratorPseudoTransientAdjointSensitivity | ( | const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > & | model, |
std::string | stepperType | ||
) |
Constructor with model and "Stepper Type" and is fully initialized with default settings.
Definition at line 34 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar >::IntegratorPseudoTransientAdjointSensitivity | ( | ) |
Destructor.
Constructor that requires a subsequent setParameterList, setStepper, and initialize calls.
Definition at line 46 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
inlinevirtual |
Destructor.
Definition at line 88 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
|
virtual |
Advance the solution to timeMax, and return true if successful.
Definition at line 55 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
overridevirtual |
Advance the solution to timeFinal, and return true if successful.
Implements Tempus::Integrator< Scalar >.
Definition at line 65 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
protected |
Definition at line 339 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
protected |
Definition at line 321 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
override |
Definition at line 266 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
override |
Definition at line 257 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
virtual |
Return adjoint sensitivity stored in gradient format.
Definition at line 249 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
overridevirtual |
Get current index.
Implements Tempus::Integrator< Scalar >.
Definition at line 115 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
inlineoverridevirtual |
Returns the IntegratorTimer_ for this Integrator.
Implements Tempus::Integrator< Scalar >.
Definition at line 113 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
|
override |
Definition at line 313 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
overridevirtual |
Get the SolutionHistory.
Implements Tempus::Integrator< Scalar >.
Definition at line 162 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
overridevirtual |
Get Status.
Implements Tempus::Integrator< Scalar >.
Definition at line 123 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
overridevirtual |
Get the Stepper.
Implements Tempus::Integrator< Scalar >.
Definition at line 137 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
inlineoverridevirtual |
Implements Tempus::Integrator< Scalar >.
Definition at line 115 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
|
overridevirtual |
Return a copy of the Tempus ParameterList.
Implements Tempus::Integrator< Scalar >.
Definition at line 145 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
overridevirtual |
Get current time.
Implements Tempus::Integrator< Scalar >.
Definition at line 107 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
overridevirtual |
Get the TimeStepControl.
Implements Tempus::Integrator< Scalar >.
Definition at line 170 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
override |
Definition at line 296 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
virtual |
Get current the solution, x.
Definition at line 225 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
virtual |
Get current the time derivative of the solution, xdot.
Definition at line 233 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
virtual |
Get current the second time derivative of the solution, xdotdot.
Definition at line 241 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
virtual |
Set the initial state from Thyra::VectorBase(s)
Definition at line 177 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
override |
Definition at line 278 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
overridevirtual |
Implements Tempus::Integrator< Scalar >.
Definition at line 153 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
override |
Definition at line 287 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.hpp.
|
protected |
Definition at line 174 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 169 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 172 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 170 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 173 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.
|
protected |
Definition at line 171 of file Tempus_IntegratorPseudoTransientAdjointSensitivity_decl.hpp.