Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus_IntegratorPseudoTransientForwardSensitivity_decl.hpp
Go to the documentation of this file.
1 // @HEADER
2 // ****************************************************************************
3 // Tempus: Copyright (2017) Sandia Corporation
4 //
5 // Distributed under BSD 3-clause license (See accompanying file Copyright.txt)
6 // ****************************************************************************
7 // @HEADER
8 
9 #ifndef Tempus_IntegratorPseudoTransientForwardSensitivity_decl_hpp
10 #define Tempus_IntegratorPseudoTransientForwardSensitivity_decl_hpp
11 
12 // Tempus
13 #include "Tempus_IntegratorBasic.hpp"
15 
16 namespace Tempus {
17 
18 
19 /** \brief Time integrator suitable for pseudotransient forward sensitivity
20  * analysis */
21 /**
22  * For some problems, time integrators are used to compute steady-state
23  * solutions (also known as pseudo-transient solvers). When computing
24  * sensitivities, it is not necessary in these cases to propagate sensitivities
25  * all the way through the forward time integration. Instead the steady-state
26  * is first computed as usual, and then the sensitivities are computed using
27  * a similar pseudo-transient time integration applied to the sensitivity
28  * equations with the state frozen to the computed steady-state. This
29  * integrator specializes the transient sensitivity methods implemented by
30  * Tempus::IntegratorForwardSensitivity to this case.
31  *
32  * Consider an implicit ODE f(x_dot,x,p) = 0 with a stable steady-state solution
33  * x = x^s, x_dot = 0 where f(0,x^s,p) = 0 and all of the eigenvalues of
34  * df/dx(0,x^s,p) are in the right half-plane (for an explicit ODE, the
35  * eigenvalues must be in the left half-plane). In the pseudo-transient method
36  * a time-integrator is applied to f(x_dot,x,p) = 0 until x_dot is sufficiently
37  * small. Now consider the forward sensitivity equations:
38  * df/dx_dot*z_dot + df/dx*z + df/dp = 0
39  * where z = dx/dp. For pseudo-transient forward sensitivities, the above is
40  * integrated from z(0) = 0 until z_dot is sufficiently small, in which case
41  * z^s = -(df/dx)^{-1}*(df/dp).
42  * Then the final sensitivity of g is
43  * dg/dp + dg/dx*z^s.
44  * One can see that z^s is the only steady-state solution of the sensitivity
45  * equations, since df/dx and df/dp are constant, and must be linearly stable
46  * since it has the same Jacobian matrix as the forward equations.
47  */
48 template<class Scalar>
50  virtual public Tempus::Integrator<Scalar>
51 {
52 public:
53 
54  /** \brief Constructor with ParameterList and model, and will be fully
55  * initialized. */
56  /*!
57  * In addition to all of the regular integrator options, the supplied
58  * parameter list supports the following options contained within a sublist
59  * "Sensitivities" from the top-level parameter list:
60  * <ul>
61  * <li> "Reuse State Linear Solver" (default: false) Whether to reuse the
62  * model's W matrix, solver, and preconditioner when solving the
63  * sensitivity equations. If they can be reused, substantial savings
64  * in compute time are possible.
65  * <li> "Force W Update" (default: false) When reusing the solver as above
66  * whether to force recomputation of W. This can be necessary when
67  * the solver overwrites it during the solve-phase (e.g., by a
68  * factorization).
69  * <li> "Use DfDp as Tangent" (default: false) Reinterpret the df/dp
70  * out-arg as the tangent vector (df/dx)(x,p) * dx/dp + df/dp(x,p)
71  * as described in the Tempus::CombinedForwardSensitivityModelEvaluator
72  * documentation.
73  * <li> "Sensitivity Parameter Index" (default: 0) Model evaluator
74  * parameter index for which sensitivities will be computed.
75  * <li> "Sensitivity X Tangent Index" (default: 1) If "Use DfDp as Tangent"
76  * is true, the model evaluator parameter index for passing dx/dp
77  * as a Thyra::DefaultMultiVectorProductVector.
78  * <li> "Sensitivity X-Dot Tangent Index" (default: 2) If
79  * "Use DfDp as Tangent" is true, the model evaluator parameter index
80  * for passing dx_dot/dp as a Thyra::DefaultMultiVectorProductVector.
81  * <li> "Sensitivity X-Dot-Dot Tangent Index" (default: 3) If
82  * "Use DfDp as Tangent" is true, the model evaluator parameter index
83  * for passing dx_dot_dot/dp as a
84  * Thyra::DefaultMultiVectorProductVector (if the model supports
85  * x_dot_dot).
86  * </ul>
87  */
89  Teuchos::RCP<Teuchos::ParameterList> pList,
90  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model);
91 
92  /** \brief Constructor with model and "Stepper Type" and is fully initialized with default settings. */
94  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
95  std::string stepperType);
96 
97  /// Destructor
98  /** \brief Constructor that requires a subsequent setParameterList, setStepper, and initialize calls. */
100 
101  /// Destructor
103 
104  /// \name Basic integrator methods
105  //@{
106 
107  /// Advance the solution to timeMax, and return true if successful.
108  virtual bool advanceTime();
109  /// Advance the solution to timeFinal, and return true if successful.
110  virtual bool advanceTime(const Scalar timeFinal) override;
111  /// Get current time
112  virtual Scalar getTime() const override;
113  /// Get current index
114  virtual Scalar getIndex() const override;
115  /// Get Status
116  virtual Status getStatus() const override;
117  /// Get the Stepper
118  virtual Teuchos::RCP<Stepper<Scalar> > getStepper() const override;
119  /// Return a copy of the Tempus ParameterList
120  virtual Teuchos::RCP<Teuchos::ParameterList> getTempusParameterList() override;
121  virtual void setTempusParameterList(Teuchos::RCP<Teuchos::ParameterList> pl) override;
122  /// Get the SolutionHistory
123  virtual Teuchos::RCP<const SolutionHistory<Scalar> > getSolutionHistory() const override;
124  /// Get the TimeStepControl
125  virtual Teuchos::RCP<const TimeStepControl<Scalar> > getTimeStepControl() const override;
126  virtual Teuchos::RCP<Teuchos::Time> getIntegratorTimer() const override
127  {return state_integrator_->getIntegratorTimer();}
128  virtual Teuchos::RCP<Teuchos::Time> getStepperTimer() const override
129  {return state_integrator_->getStepperTimer();}
130 
131  //@}
132 
133  /// Set the initial state from Thyra::VectorBase(s)
134  virtual void initializeSolutionHistory(
135  Scalar t0,
136  Teuchos::RCP<const Thyra::VectorBase<Scalar> > x0,
137  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xdot0 = Teuchos::null,
138  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xdotdot0 = Teuchos::null,
139  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > DxDp0 = Teuchos::null,
140  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > DxdotDp0 = Teuchos::null,
141  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > DxdotdotDp0 = Teuchos::null);
142 
143  /// Get current the solution, x
144  virtual Teuchos::RCP<const Thyra::VectorBase<Scalar> > getX() const;
145  virtual Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > getDxDp() const;
146  /// Get current the time derivative of the solution, xdot
147  virtual Teuchos::RCP<const Thyra::VectorBase<Scalar> > getXdot() const;
148  virtual Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > getDxdotDp() const;
149  /// Get current the second time derivative of the solution, xdotdot
150  virtual Teuchos::RCP<const Thyra::VectorBase<Scalar> > getXdotdot() const;
151  virtual Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > getDxdotdotDp() const;
152 
153  /// \name Overridden from Teuchos::ParameterListAcceptor
154  //@{
155  void setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & pl)
156  override;
157  Teuchos::RCP<Teuchos::ParameterList> getNonconstParameterList() override;
158  Teuchos::RCP<Teuchos::ParameterList> unsetParameterList() override;
159 
160  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters()
161  const override;
162  //@}
163 
164  /// \name Overridden from Teuchos::Describable
165  //@{
166  std::string description() const override;
167  void describe(Teuchos::FancyOStream & out,
168  const Teuchos::EVerbosityLevel verbLevel) const override;
169  //@}
170 
171 protected:
172 
173  // Create sensitivity model evaluator from application model
174  Teuchos::RCP<SensitivityModelEvaluatorBase<Scalar> >
176  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
177  const Teuchos::RCP<Teuchos::ParameterList>& inputPL);
178 
179  void buildSolutionHistory();
180 
181  Teuchos::RCP<Thyra::ModelEvaluator<Scalar> > model_;
182  Teuchos::RCP<SensitivityModelEvaluatorBase<Scalar> > sens_model_;
183  Teuchos::RCP<IntegratorBasic<Scalar> > state_integrator_;
184  Teuchos::RCP<IntegratorBasic<Scalar> > sens_integrator_;
185  Teuchos::RCP<SolutionHistory<Scalar> > solutionHistory_;
188 };
189 
190 /// Non-member constructor
191 template<class Scalar>
192 Teuchos::RCP<Tempus::IntegratorPseudoTransientForwardSensitivity<Scalar> >
194  Teuchos::RCP<Teuchos::ParameterList> pList,
195  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model);
196 
197 /// Non-member constructor
198 template<class Scalar>
199 Teuchos::RCP<Tempus::IntegratorPseudoTransientForwardSensitivity<Scalar> >
201  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
202  std::string stepperType);
203 
204 /// Non-member constructor
205 template<class Scalar>
206 Teuchos::RCP<Tempus::IntegratorPseudoTransientForwardSensitivity<Scalar> >
208 
209 } // namespace Tempus
210 
211 #endif // Tempus_IntegratorPseudoTransientForwardSensitivity_decl_hpp
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDxdotDp() const
virtual Teuchos::RCP< Stepper< Scalar > > getStepper() const override
Get the Stepper.
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &pl) override
virtual Teuchos::RCP< const SolutionHistory< Scalar > > getSolutionHistory() const override
Get the SolutionHistory.
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDxDp() const
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDxdotdotDp() const
Teuchos::RCP< Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar > > integratorPseudoTransientForwardSensitivity(Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model)
Non-member constructor.
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)
virtual bool advanceTime()
Advance the solution to timeMax, and return true if successful.
Status
Status for the Integrator, the Stepper and the SolutionState.
virtual Teuchos::RCP< Teuchos::Time > getIntegratorTimer() const override
Returns the IntegratorTimer_ for this Integrator.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
virtual Teuchos::RCP< Teuchos::ParameterList > getTempusParameterList() override
Return a copy of the Tempus ParameterList.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXdot() const
Get current the time derivative of the solution, xdot.
Teuchos::RCP< SensitivityModelEvaluatorBase< Scalar > > createSensitivityModel(const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< Teuchos::ParameterList > &inputPL)
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXdotdot() const
Get current the second time derivative of the solution, xdotdot.
Time integrator suitable for pseudotransient forward sensitivity analysis.
Thyra Base interface for time integrators. Time integrators are designed to advance the solution from...
virtual void setTempusParameterList(Teuchos::RCP< Teuchos::ParameterList > pl) override
virtual Teuchos::RCP< const TimeStepControl< Scalar > > getTimeStepControl() const override
Get the TimeStepControl.
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getX() const
Get current the solution, x.