Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus_IntegratorAdjointSensitivity_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_IntegratorAdjointSensitivity_decl_hpp
10 #define Tempus_IntegratorAdjointSensitivity_decl_hpp
11 
12 // Tempus
13 #include "Tempus_IntegratorBasic.hpp"
14 #include "Tempus_AdjointAuxSensitivityModelEvaluator.hpp"
15 
16 namespace Tempus {
17 
18 
19 /** \brief Time integrator suitable for adjoint sensitivity analysis */
20 /**
21  * This integrator implements transient adjoint sensitivities. Given a model
22  * evaluator encapsulating the equations f(x_dot,x,p) = 0, and a response
23  * function g(x,p), these equations are integrated forward in time to some
24  * final time T. Then the adjoint equations
25  * F(y) = d/dt( df/dx_dot^T*y ) - df/dx^T*y = 0
26  * are integrated backward in time starting at t = T where y is the adjoint
27  * variable. Nominally y is a multi-vector belonging to the vector space of f
28  * and the number of columns given by the number of entries in the response g.
29  * The initial conditions for y at t = T are given by
30  * y(T) = df/dx_dot(x_dot(T),x(T),p)^{-T} * dg/dx(x(T),p)^T.
31  * Then the final sensitivity of g is
32  * dg/dp(T) - int_{0}^T(df/dp^T*y)dt + dx/dp(0)^T*df/dx_dot(0)^T*y(0).
33  *
34  * This integrator supports both implicit and explicit steppers, and provides
35  * a method to compute dg/dp as described above after the reverse integration.
36  * The solution history contains solutions stored as product vectors (x,y),
37  * with y further stored as a product vector for each component of g.
38  *
39  * Because of limitations on the steppers, the implementation currently assumes
40  * df/dxdot is a constant matrix.
41  */
42 template<class Scalar>
44  virtual public Tempus::Integrator<Scalar>
45 {
46 public:
47 
48  /** \brief Constructor with ParameterList and model, and will be fully
49  * initialized. */
50  /*!
51  * In addition to all of the regular integrator options, the supplied
52  * parameter list supports the following options contained within a sublist
53  * "Sensitivities" from the top-level parameter list:
54  * <ul>
55  * <li> "Sensitivity Parameter Index", (default: 0) The model evaluator
56  * parameter index for which sensitivities will be computed.
57  * <li> "Response Function Index", (default: 0) The model evaluator
58  * response index for which sensitivities will be computed.
59  * <li> "Response Depends on Parameters", (default: true) Whether the
60  * response function depends on the parameter vector p for which
61  * sensitivities will be computed. If set to false, the dg/dp
62  * term in the sensitivity formula will not be computed.
63  * <li> "Residual Depends on Parameters", (default: true) Whether the
64  * model residual f depends on the parameter vector p for which
65  * sensitivities will be computed. If set to false, its contribution
66  * to the sensitivities will not be computed.
67  * <li> "IC Depends on Parameters", (default: true) Whether the initial
68  * conditions depend on the parameter vector p for which
69  * sensitivities will be computed. If set to false, its contribution
70  * to the sensitivities will not be computed.
71  * <li> "Mass Matrix Is Constant" (default: true) Whether the mass matrix
72  * df/dx_dot is a constant matrix. As describe above, this is
73  * currently required to be true.
74  * <li> "Mass Matrix Is Identity" (default: false) Whether the mass matrix
75  * is the identity matrix, in which some computations can be skipped.
76  * </ul>
77  */
79  Teuchos::RCP<Teuchos::ParameterList> pList,
80  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model);
81 
82  /// Destructor
83  /** \brief Constructor that requires a subsequent setParameterList, setStepper, and initialize calls. */
85 
86  /// Destructor
88 
89  /// \name Basic integrator methods
90  //@{
91 
92  /// Advance the solution to timeMax, and return true if successful.
93  virtual bool advanceTime();
94  /// Advance the solution to timeFinal, and return true if successful.
95  virtual bool advanceTime(const Scalar timeFinal) override;
96  /// Get current time
97  virtual Scalar getTime() const override;
98  /// Get current index
99  virtual int getIndex() const override;
100  /// Get Status
101  virtual Status getStatus() const override;
102  /// Get the Stepper
103  virtual Teuchos::RCP<Stepper<Scalar> > getStepper() const override;
104  /// Return a copy of the Tempus ParameterList
105  virtual Teuchos::RCP<Teuchos::ParameterList> getTempusParameterList() override;
106  virtual void setTempusParameterList(Teuchos::RCP<Teuchos::ParameterList> pl) override;
107  /// Get the SolutionHistory
108  virtual Teuchos::RCP<const SolutionHistory<Scalar> > getSolutionHistory() const override;
109  /// Get the TimeStepControl
110  virtual Teuchos::RCP<const TimeStepControl<Scalar> > getTimeStepControl() const override;
111  virtual Teuchos::RCP<TimeStepControl<Scalar> > getNonConstTimeStepControl() override;
112  /// Returns the IntegratorTimer_ for this Integrator
113  virtual Teuchos::RCP<Teuchos::Time> getIntegratorTimer() const override
114  { return state_integrator_->getIntegratorTimer();}
115  virtual Teuchos::RCP<Teuchos::Time> getStepperTimer() const override
116  { return state_integrator_->getStepperTimer();}
117 
118  //@}
119 
120  /// Set the initial state from Thyra::VectorBase(s)
121  virtual void initializeSolutionHistory(
122  Scalar t0,
123  Teuchos::RCP<const Thyra::VectorBase<Scalar> > x0,
124  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xdot0 = Teuchos::null,
125  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xdotdot0 = Teuchos::null,
126  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > DxDp0 = Teuchos::null,
127  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > DxdotDp0 = Teuchos::null,
128  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > DxdotdotDp0 = Teuchos::null);
129 
130  /// Set the Observer
131  virtual void setObserver(
132  Teuchos::RCP<IntegratorObserver<Scalar> > obs = Teuchos::null);
133  /// Initializes the Integrator after set* function calls
134  virtual void initialize();
135 
136  /// Get current the solution, x
137  virtual Teuchos::RCP<const Thyra::VectorBase<Scalar> > getX() const;
138  /// Get current the time derivative of the solution, xdot
139  virtual Teuchos::RCP<const Thyra::VectorBase<Scalar> > getXdot() const;
140  /// Get current the second time derivative of the solution, xdotdot
141  virtual Teuchos::RCP<const Thyra::VectorBase<Scalar> > getXdotdot() const;
142 
143  /// Return adjoint sensitivity stored in gradient format
144  virtual Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > getDgDp() const;
145 
146  /// \name Overridden from Teuchos::ParameterListAcceptor
147  //@{
148  void setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & pl)
149  override;
150  Teuchos::RCP<Teuchos::ParameterList> getNonconstParameterList() override;
151  Teuchos::RCP<Teuchos::ParameterList> unsetParameterList() override;
152 
153  Teuchos::RCP<const Teuchos::ParameterList> getValidParameters()
154  const override;
155  //@}
156 
157  /// \name Overridden from Teuchos::Describable
158  //@{
159  std::string description() const override;
160  void describe(Teuchos::FancyOStream & out,
161  const Teuchos::EVerbosityLevel verbLevel) const override;
162  //@}
163 
164 protected:
165 
166  // Create sensitivity model evaluator from application model
167  Teuchos::RCP<AdjointAuxSensitivityModelEvaluator<Scalar> >
169  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
170  const Teuchos::RCP<Teuchos::ParameterList>& inputPL);
171 
173  const Teuchos::RCP<const SolutionHistory<Scalar> >& state_solution_history,
174  const Teuchos::RCP<const SolutionHistory<Scalar> >& adjoint_solution_history);
175 
176  Teuchos::RCP<Thyra::ModelEvaluator<Scalar> > model_;
177  Teuchos::RCP<AdjointAuxSensitivityModelEvaluator<Scalar> > adjoint_model_;
178  Teuchos::RCP<IntegratorBasic<Scalar> > state_integrator_;
179  Teuchos::RCP<IntegratorBasic<Scalar> > adjoint_integrator_;
180  Teuchos::RCP<SolutionHistory<Scalar> > solutionHistory_;
181  int p_index_;
182  int g_index_;
187  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > dxdp_init_;
188  Teuchos::RCP<Thyra::MultiVectorBase<Scalar> > dgdp_;
189 };
190 
191 /// Non-member constructor
192 template<class Scalar>
193 Teuchos::RCP<IntegratorAdjointSensitivity<Scalar> >
195  Teuchos::RCP<Teuchos::ParameterList> pList,
196  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model);
197 
198 /// Non-member constructor
199 template<class Scalar>
200 Teuchos::RCP<IntegratorAdjointSensitivity<Scalar> >
202 
203 } // namespace Tempus
204 
205 #endif // Tempus_IntegratorAdjointSensitivity_decl_hpp
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
Teuchos::RCP< Thyra::MultiVectorBase< Scalar > > dgdp_
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList() override
virtual Teuchos::RCP< TimeStepControl< Scalar > > getNonConstTimeStepControl() override
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< AdjointAuxSensitivityModelEvaluator< Scalar > > adjoint_model_
virtual Teuchos::RCP< Teuchos::Time > getStepperTimer() const override
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList() override
Teuchos::RCP< IntegratorAdjointSensitivity< Scalar > > integratorAdjointSensitivity(Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model)
Non-member constructor.
Teuchos::RCP< SolutionHistory< Scalar > > solutionHistory_
void buildSolutionHistory(const Teuchos::RCP< const SolutionHistory< Scalar > > &state_solution_history, const Teuchos::RCP< const SolutionHistory< Scalar > > &adjoint_solution_history)
Teuchos::RCP< IntegratorBasic< Scalar > > adjoint_integrator_
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getX() const
Get current the solution, x.
virtual Teuchos::RCP< Teuchos::Time > getIntegratorTimer() const override
Returns the IntegratorTimer_ for this Integrator.
virtual Teuchos::RCP< const SolutionHistory< Scalar > > getSolutionHistory() const override
Get the SolutionHistory.
Status
Status for the Integrator, the Stepper and the SolutionState.
Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > dxdp_init_
virtual Teuchos::RCP< Teuchos::ParameterList > getTempusParameterList() override
Return a copy of the Tempus ParameterList.
IntegratorObserver class for time integrators.
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 Teuchos::RCP< const TimeStepControl< Scalar > > getTimeStepControl() const override
Get the TimeStepControl.
virtual void setTempusParameterList(Teuchos::RCP< Teuchos::ParameterList > pl) override
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
virtual bool advanceTime()
Advance the solution to timeMax, and return true if successful.
virtual void setObserver(Teuchos::RCP< IntegratorObserver< Scalar > > obs=Teuchos::null)
Set the Observer.
Teuchos::RCP< IntegratorBasic< Scalar > > state_integrator_
Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > model_
Thyra Base interface for time integrators. Time integrators are designed to advance the solution from...
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &pl) override
Time integrator suitable for adjoint sensitivity analysis.
virtual int getIndex() const override
Get current index.
Teuchos::RCP< AdjointAuxSensitivityModelEvaluator< Scalar > > createAdjointModel(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.
virtual Scalar getTime() const override
Get current time.
virtual Status getStatus() const override
Get Status.
virtual void initialize()
Initializes the Integrator after set* function calls.
virtual Teuchos::RCP< Stepper< Scalar > > getStepper() const override
Get the Stepper.
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDgDp() const
Return adjoint sensitivity stored in gradient format.