Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_IntegratorForwardSensitivity_impl.hpp
Go to the documentation of this file.
1 //@HEADER
2 // *****************************************************************************
3 // Tempus: Time Integration and Sensitivity Analysis Package
4 //
5 // Copyright 2017 NTESS and the Tempus contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 //@HEADER
9 
10 #ifndef Tempus_IntegratorForwardSensitivity_impl_hpp
11 #define Tempus_IntegratorForwardSensitivity_impl_hpp
12 
14 #include "Thyra_DefaultMultiVectorProductVector.hpp"
15 #include "Thyra_VectorStdOps.hpp"
16 #include "Thyra_MultiVectorStdOps.hpp"
17 
18 #include "Tempus_CombinedForwardSensitivityModelEvaluator.hpp"
20 
21 namespace Tempus {
22 
23 template <class Scalar>
26  const Teuchos::RCP<IntegratorBasic<Scalar>> &integrator,
29  &sens_stepper,
30  bool use_combined_method)
31  : model_(model),
32  integrator_(integrator),
33  sens_model_(sens_model),
34  sens_stepper_(sens_stepper),
35  use_combined_method_(use_combined_method)
36 {
37  integrator_->initialize();
38 }
39 
40 template <class Scalar>
42 {
43  integrator_ = createIntegratorBasic<Scalar>();
44  integrator_->initialize();
45 }
46 
47 template <class Scalar>
50 {
51  if (use_combined_method_)
52  integrator_->setModel(sens_model_);
53  else
54  integrator_->setStepper(sens_stepper_);
55 }
56 
57 template <class Scalar>
59  Scalar t0, Teuchos::RCP<const Thyra::VectorBase<Scalar>> x0,
65 {
66  using Teuchos::RCP;
67  using Teuchos::rcp_dynamic_cast;
68  using Thyra::assign;
69  using Thyra::createMember;
72 
73  //
74  // Create and initialize product X, Xdot, Xdotdot
75 
76  RCP<const VectorSpaceBase<Scalar>> space;
77  if (use_combined_method_)
78  space = sens_model_->get_x_space();
79  else
80  space = sens_stepper_->get_x_space();
81  RCP<DMVPV> X = rcp_dynamic_cast<DMVPV>(createMember(space));
82  RCP<DMVPV> Xdot = rcp_dynamic_cast<DMVPV>(createMember(space));
83  RCP<DMVPV> Xdotdot = rcp_dynamic_cast<DMVPV>(createMember(space));
84 
85  const int num_param = X->getNonconstMultiVector()->domain()->dim() - 1;
86  const Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
87  const Teuchos::Range1D rng(1, num_param);
88 
89  // x
90  assign(X->getNonconstMultiVector()->col(0).ptr(), *x0);
91  if (DxDp0 == Teuchos::null)
92  assign(X->getNonconstMultiVector()->subView(rng).ptr(), zero);
93  else
94  assign(X->getNonconstMultiVector()->subView(rng).ptr(), *DxDp0);
95 
96  // xdot
97  if (xdot0 == Teuchos::null)
98  assign(Xdot->getNonconstMultiVector()->col(0).ptr(), zero);
99  else
100  assign(Xdot->getNonconstMultiVector()->col(0).ptr(), *xdot0);
101  if (DxdotDp0 == Teuchos::null)
102  assign(Xdot->getNonconstMultiVector()->subView(rng).ptr(), zero);
103  else
104  assign(Xdot->getNonconstMultiVector()->subView(rng).ptr(), *DxdotDp0);
105 
106  // xdotdot
107  if (xdotdot0 == Teuchos::null)
108  assign(Xdotdot->getNonconstMultiVector()->col(0).ptr(), zero);
109  else
110  assign(Xdotdot->getNonconstMultiVector()->col(0).ptr(), *xdotdot0);
111  if (DxdotDp0 == Teuchos::null)
112  assign(Xdotdot->getNonconstMultiVector()->subView(rng).ptr(), zero);
113  else
114  assign(Xdotdot->getNonconstMultiVector()->subView(rng).ptr(), *DxdotdotDp0);
115 
116  integrator_->initializeSolutionHistory(t0, X, Xdot, Xdotdot);
117 }
118 
119 template <class Scalar>
122 {
123  using Teuchos::RCP;
124  using Teuchos::rcp_dynamic_cast;
126 
127  RCP<const DMVPV> X = rcp_dynamic_cast<const DMVPV>(integrator_->getX());
128  return X->getMultiVector()->col(0);
129 }
130 
131 template <class Scalar>
134 {
135  using Teuchos::RCP;
136  using Teuchos::rcp_dynamic_cast;
138 
139  RCP<const DMVPV> X = rcp_dynamic_cast<const DMVPV>(integrator_->getX());
140  const int num_param = X->getMultiVector()->domain()->dim() - 1;
141  const Teuchos::Range1D rng(1, num_param);
142  return X->getMultiVector()->subView(rng);
143 }
144 
145 template <class Scalar>
148 {
149  using Teuchos::RCP;
150  using Teuchos::rcp_dynamic_cast;
152 
153  RCP<const DMVPV> Xdot = rcp_dynamic_cast<const DMVPV>(integrator_->getXDot());
154  return Xdot->getMultiVector()->col(0);
155 }
156 
157 template <class Scalar>
160 {
161  using Teuchos::RCP;
162  using Teuchos::rcp_dynamic_cast;
164 
165  RCP<const DMVPV> Xdot = rcp_dynamic_cast<const DMVPV>(integrator_->getXDot());
166  const int num_param = Xdot->getMultiVector()->domain()->dim() - 1;
167  const Teuchos::Range1D rng(1, num_param);
168  return Xdot->getMultiVector()->subView(rng);
169 }
170 
171 template <class Scalar>
174 {
175  using Teuchos::RCP;
176  using Teuchos::rcp_dynamic_cast;
178 
179  RCP<const DMVPV> Xdotdot =
180  rcp_dynamic_cast<const DMVPV>(integrator_->getXDotDot());
181  return Xdotdot->getMultiVector()->col(0);
182 }
183 
184 template <class Scalar>
187 {
188  using Teuchos::RCP;
189  using Teuchos::rcp_dynamic_cast;
191 
192  RCP<const DMVPV> Xdotdot =
193  rcp_dynamic_cast<const DMVPV>(integrator_->getXDotDot());
194  const int num_param = Xdotdot->getMultiVector()->domain()->dim() - 1;
195  const Teuchos::Range1D rng(1, num_param);
196  return Xdotdot->getMultiVector()->subView(rng);
197 }
198 
199 template <class Scalar>
202 {
203  typedef Thyra::ModelEvaluatorBase MEB;
204 
205  // Compute g which is computed by response 1 of the
206  // sensitivity model evaluator
208  if (use_combined_method_)
209  smodel = sens_model_;
210  else
211  smodel = sens_stepper_->getModel();
212  MEB::InArgs<Scalar> inargs = smodel->getNominalValues();
213  MEB::OutArgs<Scalar> outargs = smodel->createOutArgs();
214  inargs.set_t(integrator_->getTime());
215  inargs.set_x(integrator_->getX());
216  if (inargs.supports(MEB::IN_ARG_x_dot))
217  inargs.set_x_dot(integrator_->getXDot());
218  if (inargs.supports(MEB::IN_ARG_x_dot_dot))
219  inargs.set_x_dot_dot(integrator_->getXDotDot());
220 
222  Thyra::createMember(smodel->get_g_space(1));
223  outargs.set_g(1, g);
224 
225  smodel->evalModel(inargs, outargs);
226  return g;
227 }
228 
229 template <class Scalar>
232 {
233  typedef Thyra::ModelEvaluatorBase MEB;
235 
236  // Compute final dg/dp which is computed by response 0 of the
237  // sensitivity model evaluator
239  if (use_combined_method_)
240  smodel = sens_model_;
241  else
242  smodel = sens_stepper_->getModel();
243  MEB::InArgs<Scalar> inargs = smodel->getNominalValues();
244  MEB::OutArgs<Scalar> outargs = smodel->createOutArgs();
245  inargs.set_t(integrator_->getTime());
246  inargs.set_x(integrator_->getX());
247  if (inargs.supports(MEB::IN_ARG_x_dot))
248  inargs.set_x_dot(integrator_->getXDot());
249  if (inargs.supports(MEB::IN_ARG_x_dot_dot))
250  inargs.set_x_dot_dot(integrator_->getXDotDot());
251 
253  Thyra::createMember(smodel->get_g_space(0));
254  Teuchos::RCP<DMVPV> dgdp = Teuchos::rcp_dynamic_cast<DMVPV>(G);
255  outargs.set_g(0, G);
256 
257  smodel->evalModel(inargs, outargs);
258  return dgdp->getMultiVector();
259 }
260 
261 template <class Scalar>
263 {
264  std::string name = "Tempus::IntegratorForwardSensitivity";
265  return (name);
266 }
267 
268 template <class Scalar>
270  Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
271 {
272  auto l_out = Teuchos::fancyOStream(out.getOStream());
273  Teuchos::OSTab ostab(*l_out, 2, this->description());
274  l_out->setOutputToRootOnly(0);
275 
276  *l_out << description() << "::describe" << std::endl;
277  integrator_->describe(*l_out, verbLevel);
278 }
279 
280 template <class Scalar>
282 {
283  if (use_combined_method_) return SensitivityStepMode::Combined;
284  return sens_stepper_->getStepMode(); // Staggered case
285 }
286 
288 template <class Scalar>
293  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar>> &sens_residual_model,
294  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar>> &sens_solve_model)
295 {
298 
299  // 1. create integrator
300  auto fwd_integrator = createIntegratorBasic<Scalar>(pList, model);
301 
302  // 2. set parameter list
305  {
307  fwd_integrator->getValidParameters();
310  pl->setParameters(*integrator_pl);
311  Teuchos::ParameterList &spl = pl->sublist("Sensitivities");
312  spl.setParameters(*sensitivity_pl);
313  spl.set("Sensitivity Method", "Combined");
314  spl.set("Reuse State Linear Solver", false);
315  }
316  pList->setParametersNotAlreadySet(*pl);
317 
318  auto sens_pl = Teuchos::sublist(pList, "Sensitivities", false);
319  std::string integratorName =
320  pList->get<std::string>("Integrator Name", "Default Integrator");
321  std::string stepperName =
322  pList->sublist(integratorName).get<std::string>("Stepper Name");
323  auto stepper_pl = Teuchos::sublist(pList, stepperName, true);
324  std::string sensitivity_method =
325  sens_pl->get<std::string>("Sensitivity Method");
326  bool use_combined_method = sensitivity_method == "Combined";
327 
328  // 3. create sensitivity model and stepper
329  // createSensitivityModelAndStepper
330  {
331  sens_pl->remove("Sensitivity Method");
332 
333  if (use_combined_method) {
334  sens_pl->remove("Reuse State Linear Solver");
335  sens_model = wrapCombinedFSAModelEvaluator(model, sens_residual_model,
336  sens_solve_model, sens_pl);
337  fwd_integrator = createIntegratorBasic<Scalar>(pList, sens_model);
338  }
339  else {
340  sens_stepper =
342  model, sens_residual_model, sens_solve_model, stepper_pl,
343  sens_pl));
344  auto fsa_staggered_me =
345  Teuchos::rcp_const_cast<Thyra::ModelEvaluator<Scalar>>(
346  sens_stepper->getModel());
347  fwd_integrator = createIntegratorBasic<Scalar>(pList, fsa_staggered_me);
348  fwd_integrator->setStepper(sens_stepper);
349  }
350  }
351 
352  // 4. initialize propoer integrator
355  model, fwd_integrator, sens_model, sens_stepper,
356  use_combined_method));
357  return (integrator);
358 }
359 
361 template <class Scalar>
364 {
367  return (integrator);
368 }
369 
370 } // namespace Tempus
371 #endif // Tempus_IntegratorForwardSensitivity_impl_hpp
Teuchos::RCP< IntegratorBasic< Scalar > > integrator_
Teuchos::RCP< IntegratorForwardSensitivity< Scalar > > createIntegratorForwardSensitivity(Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar >> &model, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar >> &sens_residual_model, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar >> &sens_solve_model)
Nonmember constructor.
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXDotDot() const
Get current the second time derivative of the solution, xdotdot, only. This is the first block only a...
virtual void setStepper(Teuchos::RCP< Thyra::ModelEvaluator< Scalar >> model)
Set the Stepper.
SensitivityStepMode getStepMode() const
What mode the current time integration step is in.
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDXDotDp() const
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDXDotDotDp() const
A stepper implementing staggered forward sensitivity analysis.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDgDp() const
Return forward sensitivity stored in Jacobian format.
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDxDp() const
Get the forward sensitivities .
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Teuchos::RCP< SensitivityModelEvaluatorBase< Scalar > > wrapCombinedFSAModelEvaluator(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &sens_residual_model, const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &sens_solve_model, const Teuchos::RCP< const Teuchos::ParameterList > &pList=Teuchos::null)
RCP< const MultiVectorBase< Scalar > > getMultiVector() const
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
ParameterList & setParameters(const ParameterList &source)
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getG() const
Return response function g.
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getX() const
Get the current solution, x, only. If looking for the solution vector and the sensitivities, use SolutionState-&gt;getX() which will return a Block MultiVector with the first block containing the current solution, x, and the remaining blocks are the forward sensitivities .
RCP< std::basic_ostream< char_type, traits_type > > getOStream()
ParameterList & setParametersNotAlreadySet(const ParameterList &source)
RCP< MultiVectorBase< Scalar > > getNonconstMultiVector()
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXDot() const
Get current the time derivative of the solution, xdot, only. This is the first block only and not the...
A ModelEvaluator decorator for sensitivity analysis.
Time integrator implementing forward sensitivity analysis.
virtual void initializeSolutionHistory(Teuchos::RCP< SolutionState< Scalar >> state=Teuchos::null)
Set the initial state which has the initial conditions.