Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Tempus_AdjointSensitivityModelEvaluator_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_AdjointSensitivityModelEvaluator_decl_hpp
10 #define Tempus_AdjointSensitivityModelEvaluator_decl_hpp
11 
12 #include "Thyra_StateFuncModelEvaluatorBase.hpp"
13 #include "Thyra_DefaultMultiVectorProductVectorSpace.hpp"
14 #include "Thyra_DefaultMultiVectorProductVector.hpp"
15 
16 #include "Tempus_SolutionHistory.hpp"
17 
18 namespace Tempus {
19 
20 /** \brief ModelEvaluator for forming adjoint sensitivity equations. */
21 /**
22  * This class wraps a given ModelEvalutor encapsulating f(x_dot,x,p) and creates
23  * a new "residual" for the adjoint sensitivity equations:
24  * F(y) = d/dt( df/dx_dot^T*y ) - df/dx^T*y + dg/dx^T = 0
25  * where y is the adjoint variable. Nominally y is a multi-vector belonging
26  * to the vector space of the residual f with number of columns equal to the
27  * dimension of g. To satisfy the model evaluator interface, y is converted
28  * to a product vector formed by its columns. This is the form of the adjoint
29  * equations suitable for computing pseudotransientsensitivities of a response
30  * function g(x(T),p) or adjoint sensitivities for a time-integrated response
31  * function int_0^T g(x(t),p).
32  *
33  * To compute adjoint sensitivities, the equations f(x_dot,x,p) = 0 are
34  * integrated forward in time to some final time T, with the adjoint equations
35  * above integrated backward in time starting at T. Since the tempus
36  * integrator can only integrate forward in time, we define tau = T - t and
37  * transform the adjoint equations to:
38  * F(y) = d/dtau( df/dx_dot^T*y ) + df/dx^T*y - dg/dx^T = 0.
39  * The initial conditions for y are y(T) = 0. The sensitivity of g(T) is then
40  * int_0^T(dg/dp(T) - df/dp^T*y)dt + dx/dp(0)^T*df/dx_dot(0)^T*y(0)
41  * for transient adjoint sensitivites and is
42  * dg/dp(T) - df/dp^T*y
43  * for pseudotransient.
44  *
45  * This model evaluator supports both implicit and explict forms of the
46  * governing equations. However it assumes df/dx_dot is constant with respect
47  * to x_dot, x, and t so that the adjoint equations become
48  * F(y) = df/dxdot^T*y_dot + df/dx^T*y - dg/dx^T = 0.
49  * Moving beyond this assumption will require modifications to the steppers
50  * in how they generate time-derivative terms.
51  */
52 template <typename Scalar>
54  public Thyra::StateFuncModelEvaluatorBase<Scalar> {
55 public:
56  typedef Thyra::VectorBase<Scalar> Vector;
57  typedef Thyra::MultiVectorBase<Scalar> MultiVector;
58 
59  //! Constructor
60  /*!
61  * t_final is the final integration time used to implement the time
62  * transformation described above. num_adjoint is the number of adjoint
63  * variables defining the number of columns in y, which is determined by
64  * the number of elements in the vector space for the response g.
65  * The optionally supplied parameter list supports the following options:
66  * <ul>
67  * <li> "Sensitivity Parameter Index", (default: 0) The model evaluator
68  * parameter index for which sensitivities will be computed.
69  * <li> "Response Function Index", (default: 0) The model evaluator
70  * response index for which sensitivities will 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  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> > & model,
80  const Scalar& t_final,
81  const bool is_pseudotransient,
82  const Teuchos::RCP<const Teuchos::ParameterList>& pList = Teuchos::null);
83 
84  //! Get the underlying model 'f'
85  Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> > getModel() const
86  { return model_; }
87 
88  //! Set solution history from forward evaluation
90  const Teuchos::RCP<const Tempus::SolutionHistory<Scalar> >& sh);
91 
92  /** \name Public functions overridden from ModelEvaulator. */
93  //@{
94 
95  Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> > get_p_space(int p) const;
96 
97  Teuchos::RCP<const Teuchos::Array<std::string> > get_p_names(int p) const;
98 
99  Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> > get_x_space() const;
100 
101  Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> > get_f_space() const;
102 
103  Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> > get_g_space(int j) const;
104 
105  Teuchos::RCP<Thyra::LinearOpBase<Scalar> > create_W_op() const;
106 
107  Teuchos::RCP<const Thyra::LinearOpWithSolveFactoryBase<Scalar> >
108  get_W_factory() const;
109 
110  Thyra::ModelEvaluatorBase::InArgs<Scalar> createInArgs() const;
111 
112  Thyra::ModelEvaluatorBase::InArgs<Scalar> getNominalValues() const;
113 
114  //@}
115 
116  static Teuchos::RCP<const Teuchos::ParameterList> getValidParameters();
117 
118 private:
119 
120  typedef Thyra::DefaultMultiVectorProductVectorSpace<Scalar> DMVPVS;
121  typedef Thyra::DefaultMultiVectorProductVector<Scalar> DMVPV;
122 
123  Thyra::ModelEvaluatorBase::OutArgs<Scalar> createOutArgsImpl() const;
124 
125  void evalModelImpl(
126  const Thyra::ModelEvaluatorBase::InArgs<Scalar> &inArgs,
127  const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs) const;
128 
129 
130  Thyra::ModelEvaluatorBase::InArgs<Scalar> prototypeInArgs_;
131  Thyra::ModelEvaluatorBase::OutArgs<Scalar> prototypeOutArgs_;
132 
133  Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> > model_;
134 
135  Teuchos::RCP<const DMVPVS> adjoint_space_;
136  Teuchos::RCP<const DMVPVS> residual_space_;
137  Teuchos::RCP<const DMVPVS> response_space_;
138  Teuchos::RCP<const Tempus::SolutionHistory<Scalar> > sh_;
139  Scalar t_final_;
143  int p_index_;
144  int g_index_;
146 
148  mutable Teuchos::RCP<Thyra::VectorBase<Scalar> > my_x_dot_;
149  mutable Teuchos::RCP<Thyra::LinearOpBase<Scalar> > my_dfdx_;
150  mutable Teuchos::RCP<Thyra::LinearOpBase<Scalar> > my_dfdxdot_;
151  mutable Teuchos::RCP<Thyra::LinearOpBase<Scalar> > my_dfdp_op_;
152  mutable Teuchos::RCP<Thyra::MultiVectorBase<Scalar> > my_dfdp_mv_;
153  mutable Teuchos::RCP<Thyra::MultiVectorBase<Scalar> > my_dgdx_mv_;
154  mutable Teuchos::RCP<Tempus::SolutionState<Scalar> > forward_state_;
155  mutable Scalar t_interp_;
156 };
157 
158 } // namespace Tempus
159 
160 #endif
static Teuchos::RCP< const Teuchos::ParameterList > getValidParameters()
AdjointSensitivityModelEvaluator(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model, const Scalar &t_final, const bool is_pseudotransient, const Teuchos::RCP< const Teuchos::ParameterList > &pList=Teuchos::null)
Constructor.
Teuchos::RCP< const Tempus::SolutionHistory< Scalar > > sh_
Thyra::ModelEvaluatorBase::InArgs< Scalar > createInArgs() const
Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > getModel() const
Get the underlying model &#39;f&#39;.
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_x_space() const
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int p) const
Teuchos::RCP< const Thyra::LinearOpWithSolveFactoryBase< Scalar > > get_W_factory() const
void evalModelImpl(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Thyra::ModelEvaluatorBase::OutArgs< Scalar > createOutArgsImpl() const
Thyra::ModelEvaluatorBase::OutArgs< Scalar > prototypeOutArgs_
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_f_space() const
Teuchos::RCP< Thyra::MultiVectorBase< Scalar > > my_dfdp_mv_
Thyra::DefaultMultiVectorProductVector< Scalar > DMVPV
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_g_space(int j) const
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > model_
void setForwardSolutionHistory(const Teuchos::RCP< const Tempus::SolutionHistory< Scalar > > &sh)
Set solution history from forward evaluation.
Thyra::ModelEvaluatorBase::InArgs< Scalar > getNominalValues() const
ModelEvaluator for forming adjoint sensitivity equations.
Thyra::DefaultMultiVectorProductVectorSpace< Scalar > DMVPVS
Teuchos::RCP< Tempus::SolutionState< Scalar > > forward_state_
Teuchos::RCP< Thyra::LinearOpBase< Scalar > > create_W_op() const
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_p_space(int p) const
Teuchos::RCP< Thyra::MultiVectorBase< Scalar > > my_dgdx_mv_