| 
    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 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 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 | 
Basic integrator methods  | |
Protected Member Functions | |
| Teuchos::RCP < Tempus::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 Attributes | |
| Teuchos::RCP < Thyra::ModelEvaluator < Scalar > >  | model_ | 
| Teuchos::RCP < Tempus::AdjointAuxSensitivityModelEvaluator < Scalar > >  | adjoint_model_ | 
| Teuchos::RCP< IntegratorBasic < Scalar > >  | state_integrator_ | 
| 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_ | 
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 43 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 25 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
| Tempus::IntegratorAdjointSensitivity< Scalar >::IntegratorAdjointSensitivity | ( | ) | 
Destructor.
Constructor that requires a subsequent setParameterList, setStepper, and initialize calls.
Definition at line 47 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  inlinevirtual | 
Destructor.
Definition at line 87 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
      
  | 
  virtual | 
Advance the solution to timeMax, and return true if successful.
Definition at line 56 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Advance the solution to timeFinal, and return true if successful.
Implements Tempus::Integrator< Scalar >.
Definition at line 66 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  protected | 
Definition at line 452 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  protected | 
Definition at line 429 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  override | 
Definition at line 358 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  override | 
Definition at line 349 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  virtual | 
Return adjoint sensitivity stored in gradient format.
Definition at line 341 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Get current index.
Implements Tempus::Integrator< Scalar >.
Definition at line 232 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  inlineoverridevirtual | 
Returns the IntegratorTimer_ for this Integrator.
Implements Tempus::Integrator< Scalar >.
Definition at line 113 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
      
  | 
  override | 
Definition at line 421 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Implements Tempus::Integrator< Scalar >.
Definition at line 295 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Get the SolutionHistory.
Implements Tempus::Integrator< Scalar >.
Definition at line 279 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Get Status.
Implements Tempus::Integrator< Scalar >.
Definition at line 240 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Get the Stepper.
Implements Tempus::Integrator< Scalar >.
Definition at line 254 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  inlineoverridevirtual | 
Implements Tempus::Integrator< Scalar >.
Definition at line 115 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
      
  | 
  overridevirtual | 
Return a copy of the Tempus ParameterList.
Implements Tempus::Integrator< Scalar >.
Definition at line 262 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Get current time.
Implements Tempus::Integrator< Scalar >.
Definition at line 224 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Get the TimeStepControl.
Implements Tempus::Integrator< Scalar >.
Definition at line 287 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  override | 
Definition at line 397 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  virtual | 
Get current the solution, x.
Definition at line 317 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  virtual | 
Get current the time derivative of the solution, xdot.
Definition at line 325 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  virtual | 
Get current the second time derivative of the solution, xdotdot.
Definition at line 333 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  virtual | 
Set the initial state from Thyra::VectorBase(s)
Definition at line 302 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  override | 
Definition at line 370 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Implements Tempus::Integrator< Scalar >.
Definition at line 270 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  override | 
Definition at line 388 of file Tempus_IntegratorAdjointSensitivity_impl.hpp.
      
  | 
  protected | 
Definition at line 173 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
      
  | 
  protected | 
Definition at line 171 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
      
  | 
  protected | 
Definition at line 182 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
      
  | 
  protected | 
Definition at line 181 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
      
  | 
  protected | 
Definition at line 178 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
      
  | 
  protected | 
Definition at line 177 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
      
  | 
  protected | 
Definition at line 176 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 170 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
      
  | 
  protected | 
Definition at line 175 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
      
  | 
  protected | 
Definition at line 174 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.
      
  | 
  protected | 
Definition at line 172 of file Tempus_IntegratorAdjointSensitivity_decl.hpp.