| 
    Tempus
    Version of the Day
    
   Time Integration 
   | 
 
Time integrator suitable for pseudotransient forward sensitivity analysis. More...
#include <Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp>
  
 Public Member Functions | |
| IntegratorPseudoTransientForwardSensitivity (Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model) | |
| Constructor with ParameterList and model, and will be fully initialized.  More... | |
| IntegratorPseudoTransientForwardSensitivity (const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, std::string stepperType) | |
| Constructor with model and "Stepper Type" and is fully initialized with default settings.  More... | |
| IntegratorPseudoTransientForwardSensitivity () | |
| Destructor.  More... | |
| virtual | ~IntegratorPseudoTransientForwardSensitivity () | 
| 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::MultiVectorBase< Scalar > >  | getDxDp () const | 
| virtual Teuchos::RCP< const  Thyra::VectorBase< Scalar > >  | getXdot () const | 
| Get current the time derivative of the solution, xdot.  More... | |
| virtual Teuchos::RCP< const  Thyra::MultiVectorBase< Scalar > >  | getDxdotDp () const | 
| 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 > >  | getDxdotdotDp () const | 
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 < SensitivityModelEvaluatorBase < 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 < SensitivityModelEvaluatorBase < Scalar > >  | sens_model_ | 
| Teuchos::RCP< IntegratorBasic < Scalar > >  | state_integrator_ | 
| Teuchos::RCP< IntegratorBasic < Scalar > >  | sens_integrator_ | 
| Teuchos::RCP< SolutionHistory < Scalar > >  | solutionHistory_ | 
| bool | reuse_solver_ | 
| bool | force_W_update_ | 
Time integrator suitable for pseudotransient forward 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 sensitivity equations with the state frozen to the computed steady-state. This integrator specializes the transient sensitivity methods implemented by Tempus::IntegratorForwardSensitivity 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 forward sensitivity equations: df/dx_dot*z_dot + df/dx*z + df/dp = 0 where z = dx/dp. For pseudo-transient forward sensitivities, the above is integrated from z(0) = 0 until z_dot is sufficiently small, in which case z^s = -(df/dx)^{-1}*(df/dp). Then the final sensitivity of g is dg/dp + dg/dx*z^s. One can see that z^s is the only steady-state solution of the sensitivity equations, since df/dx and df/dp are constant, and must be linearly stable since it has the same Jacobian matrix as the forward equations.
Definition at line 49 of file Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp.
| Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::IntegratorPseudoTransientForwardSensitivity | ( | 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_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
| Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::IntegratorPseudoTransientForwardSensitivity | ( | 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 35 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
| Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar >::IntegratorPseudoTransientForwardSensitivity | ( | ) | 
Destructor.
Constructor that requires a subsequent setParameterList, setStepper, and initialize calls.
Definition at line 49 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  inlinevirtual | 
Destructor.
Definition at line 102 of file Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp.
      
  | 
  virtual | 
Advance the solution to timeMax, and return true if successful.
Definition at line 60 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Advance the solution to timeFinal, and return true if successful.
Implements Tempus::Integrator< Scalar >.
Definition at line 89 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  protected | 
Definition at line 425 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  protected | 
Definition at line 405 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  override | 
Definition at line 344 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  override | 
Definition at line 335 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  virtual | 
Definition at line 319 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  virtual | 
Definition at line 289 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  virtual | 
Definition at line 259 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Get current index.
Implements Tempus::Integrator< Scalar >.
Definition at line 126 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  inlineoverridevirtual | 
Returns the IntegratorTimer_ for this Integrator.
Implements Tempus::Integrator< Scalar >.
Definition at line 127 of file Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp.
      
  | 
  override | 
Definition at line 397 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Implements Tempus::Integrator< Scalar >.
Definition at line 189 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Get the SolutionHistory.
Implements Tempus::Integrator< Scalar >.
Definition at line 173 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Get Status.
Implements Tempus::Integrator< Scalar >.
Definition at line 134 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Get the Stepper.
Implements Tempus::Integrator< Scalar >.
Definition at line 148 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  inlineoverridevirtual | 
Implements Tempus::Integrator< Scalar >.
Definition at line 129 of file Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp.
      
  | 
  overridevirtual | 
Return a copy of the Tempus ParameterList.
Implements Tempus::Integrator< Scalar >.
Definition at line 156 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Get current time.
Implements Tempus::Integrator< Scalar >.
Definition at line 118 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Get the TimeStepControl.
Implements Tempus::Integrator< Scalar >.
Definition at line 181 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  override | 
Definition at line 378 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  virtual | 
Get current the solution, x.
Definition at line 245 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  virtual | 
Get current the time derivative of the solution, xdot.
Definition at line 275 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  virtual | 
Get current the second time derivative of the solution, xdotdot.
Definition at line 305 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  virtual | 
Set the initial state from Thyra::VectorBase(s)
Definition at line 196 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  override | 
Definition at line 356 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  overridevirtual | 
Implements Tempus::Integrator< Scalar >.
Definition at line 164 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  override | 
Definition at line 369 of file Tempus_IntegratorPseudoTransientForwardSensitivity_impl.hpp.
      
  | 
  protected | 
Definition at line 188 of file Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp.
      
  | 
  protected | 
Definition at line 182 of file Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp.
      
  | 
  protected | 
Definition at line 187 of file Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp.
      
  | 
  protected | 
Definition at line 185 of file Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp.
      
  | 
  protected | 
Definition at line 183 of file Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp.
      
  | 
  protected | 
Definition at line 186 of file Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp.
      
  | 
  protected | 
Definition at line 184 of file Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp.