Tempus  Version of the Day
Time Integration
Tempus_IntegratorPseudoTransientAdjointSensitivity_impl.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_IntegratorPseudoTransientAdjointSensitivity_impl_hpp
10 #define Tempus_IntegratorPseudoTransientAdjointSensitivity_impl_hpp
11 
12 #include "Thyra_DefaultMultiVectorProductVector.hpp"
13 #include "Thyra_VectorStdOps.hpp"
14 #include "Thyra_MultiVectorStdOps.hpp"
16 
17 
18 namespace Tempus {
19 
20 template<class Scalar>
23  Teuchos::RCP<Teuchos::ParameterList> inputPL,
24  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
25  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& adjoint_model)
26 {
27  model_ = model;
28  adjoint_model_ = adjoint_model;
29  state_integrator_ = integratorBasic<Scalar>(inputPL, model_);
30  sens_model_ = createSensitivityModel(model_, adjoint_model_, inputPL);
31  sens_integrator_ = integratorBasic<Scalar>(inputPL, sens_model_);
32 }
33 
34 template<class Scalar>
37  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
38  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& adjoint_model,
39  std::string stepperType)
40 {
41  model_ = model;
42  adjoint_model_ = adjoint_model;
43  state_integrator_ = integratorBasic<Scalar>(model_, stepperType);
44  sens_model_ = createSensitivityModel(model_, adjoint_model_, Teuchos::null);
45  sens_integrator_ = integratorBasic<Scalar>(sens_model_, stepperType);
46 }
47 
48 template<class Scalar>
51  Teuchos::RCP<Teuchos::ParameterList> inputPL,
52  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model) :
54  inputPL, model, Thyra::implicitAdjointModelEvaluator(model))
55 {
56 }
57 
58 template<class Scalar>
61  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
62  std::string stepperType) :
64  model, Thyra::implicitAdjointModelEvaluator(model), stepperType)
65 {
66 }
67 
68 template<class Scalar>
71 {
72  state_integrator_ = integratorBasic<Scalar>();
73  sens_integrator_ = integratorBasic<Scalar>();
74 }
75 
76 template<class Scalar>
77 bool
80 {
81  const Scalar tfinal =
82  state_integrator_->getTimeStepControl()->getFinalTime();
83  return advanceTime(tfinal);
84 }
85 
86 template<class Scalar>
87 bool
89 advanceTime(const Scalar timeFinal)
90 {
91  using Teuchos::RCP;
92  using Thyra::VectorBase;
93  typedef Thyra::ModelEvaluatorBase MEB;
94 
95  // Run state integrator and get solution
96  bool state_status = state_integrator_->advanceTime(timeFinal);
97 
98  // Set solution in sensitivity ME
99  sens_model_->setForwardSolutionHistory(
100  state_integrator_->getSolutionHistory());
101 
102  // Run sensitivity integrator
103  bool sens_status = sens_integrator_->advanceTime(timeFinal);
104 
105  // Compute final dg/dp which is computed by response 0 of the adjoint
106  // model evaluator
107  MEB::InArgs<Scalar> inargs = sens_model_->getNominalValues();
108  MEB::OutArgs<Scalar> outargs = sens_model_->createOutArgs();
109  inargs.set_t(sens_integrator_->getTime());
110  inargs.set_x(sens_integrator_->getX());
111  if (inargs.supports(MEB::IN_ARG_x_dot))
112  inargs.set_x_dot(sens_integrator_->getXDot());
113  if (inargs.supports(MEB::IN_ARG_x_dot_dot))
114  inargs.set_x_dot_dot(sens_integrator_->getXDotDot());
115  RCP<VectorBase<Scalar> > G = dgdp_;
116  if (G == Teuchos::null) {
117  G = Thyra::createMember(sens_model_->get_g_space(0));
118  dgdp_ = Teuchos::rcp_dynamic_cast<DMVPV>(G);
119  }
120  outargs.set_g(0, G);
121  sens_model_->evalModel(inargs, outargs);
122 
123  buildSolutionHistory();
124 
125  return state_status && sens_status;
126 }
127 
128 template<class Scalar>
129 Scalar
131 getTime() const
132 {
133  return solutionHistory_->getCurrentTime();
134 }
135 
136 template<class Scalar>
137 int
139 getIndex() const
140 {
141  return solutionHistory_->getCurrentIndex();
142 }
143 
144 template<class Scalar>
145 Status
147 getStatus() const
148 {
149  Status state_status = state_integrator_->getStatus();
150  Status sens_status = sens_integrator_->getStatus();
151  if (state_status == FAILED || sens_status == FAILED)
152  return FAILED;
153  if (state_status == WORKING || sens_status == WORKING)
154  return WORKING;
155  return PASSED;
156 }
157 
158 template <class Scalar>
160  const Status st) {
161  state_integrator_->setStatus(st);
162  sens_integrator_->setStatus(st);
163 }
164 
165 template<class Scalar>
166 Teuchos::RCP<Stepper<Scalar> >
168 getStepper() const
169 {
170  return state_integrator_->getStepper();
171 }
172 
173 template<class Scalar>
174 Teuchos::RCP<Teuchos::ParameterList>
177 {
178  return state_integrator_->getTempusParameterList();
179 }
180 
181 template<class Scalar>
182 void
184 setTempusParameterList(Teuchos::RCP<Teuchos::ParameterList> pl)
185 {
186  state_integrator_->setTempusParameterList(pl);
187  sens_integrator_->setTempusParameterList(pl);
188 }
189 
190 template<class Scalar>
191 Teuchos::RCP<const SolutionHistory<Scalar> >
194 {
195  return solutionHistory_;
196 }
197 
198 template<class Scalar>
199 Teuchos::RCP<SolutionHistory<Scalar> >
202 {
203  return solutionHistory_;
204 }
205 
206 template<class Scalar>
207 Teuchos::RCP<const TimeStepControl<Scalar> >
210 {
211  return state_integrator_->getTimeStepControl();
212 }
213 
214 template<class Scalar>
215 Teuchos::RCP<TimeStepControl<Scalar> >
218 {
219  return state_integrator_->getNonConstTimeStepControl();
220 }
221 
222 template<class Scalar>
225  Teuchos::RCP<const Thyra::VectorBase<Scalar> > x0,
226  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xdot0,
227  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xdotdot0,
228  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > y0,
229  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > ydot0,
230  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > ydotdot0)
231 {
232  using Teuchos::RCP;
233  using Teuchos::rcp_dynamic_cast;
234  using Thyra::VectorSpaceBase;
235  using Thyra::assign;
236  using Thyra::createMember;
237 
238  //
239  // Create and initialize product X, Xdot, Xdotdot
240 
241  RCP< const VectorSpaceBase<Scalar> > space = sens_model_->get_x_space();
242  RCP<DMVPV> Y = rcp_dynamic_cast<DMVPV>(createMember(space));
243  RCP<DMVPV> Ydot = rcp_dynamic_cast<DMVPV>(createMember(space));
244  RCP<DMVPV> Ydotdot = rcp_dynamic_cast<DMVPV>(createMember(space));
245  const Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
246 
247  // y
248  if (y0 == Teuchos::null)
249  assign(Y->getNonconstMultiVector().ptr(), zero);
250  else
251  assign(Y->getNonconstMultiVector().ptr(), *y0);
252 
253  // ydot
254  if (ydot0 == Teuchos::null)
255  assign(Ydot->getNonconstMultiVector().ptr(), zero);
256  else
257  assign(Ydot->getNonconstMultiVector().ptr(), *ydot0);
258 
259  // ydotdot
260  if (ydotdot0 == Teuchos::null)
261  assign(Ydotdot->getNonconstMultiVector().ptr(), zero);
262  else
263  assign(Ydotdot->getNonconstMultiVector().ptr(), *ydotdot0);
264 
265  state_integrator_->initializeSolutionHistory(t0, x0, xdot0, xdotdot0);
266  sens_integrator_->initializeSolutionHistory(t0, Y, Ydot, Ydotdot);
267 }
268 
269 template<class Scalar>
270 Teuchos::RCP<const Thyra::VectorBase<Scalar> >
272 getX() const
273 {
274  return state_integrator_->getX();
275 }
276 
277 template<class Scalar>
278 Teuchos::RCP<const Thyra::VectorBase<Scalar> >
280 getXDot() const
281 {
282  return state_integrator_->getXDot();
283 }
284 
285 template<class Scalar>
286 Teuchos::RCP<const Thyra::VectorBase<Scalar> >
288 getXDotDot() const
289 {
290  return state_integrator_->getXDotDot();
291 }
292 
293 template<class Scalar>
294 Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> >
296 getDgDp() const
297 {
298  return dgdp_->getMultiVector();
299 }
300 
301 template<class Scalar>
302 std::string
304 description() const
305 {
306  std::string name = "Tempus::IntegratorPseudoTransientAdjointSensitivity";
307  return(name);
308 }
309 
310 template<class Scalar>
311 void
314  Teuchos::FancyOStream &out,
315  const Teuchos::EVerbosityLevel verbLevel) const
316 {
317  auto l_out = Teuchos::fancyOStream( out.getOStream() );
318  Teuchos::OSTab ostab(*l_out, 2, this->description());
319  l_out->setOutputToRootOnly(0);
320 
321  *l_out << description() << "::describe" << std::endl;
322  state_integrator_->describe(*l_out, verbLevel);
323  sens_integrator_->describe(*l_out, verbLevel);
324 }
325 
326 template<class Scalar>
327 void
329 setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & inputPL)
330 {
331  state_integrator_->setParameterList(inputPL);
332  sens_integrator_->setParameterList(inputPL);
333 }
334 
335 template<class Scalar>
336 Teuchos::RCP<Teuchos::ParameterList>
339 {
340  state_integrator_->unsetParameterList();
341  return sens_integrator_->unsetParameterList();
342 }
343 
344 template<class Scalar>
345 Teuchos::RCP<const Teuchos::ParameterList>
348 {
349  Teuchos::RCP<Teuchos::ParameterList> pl =
350  Teuchos::rcp(new Teuchos::ParameterList);
351  Teuchos::RCP<const Teuchos::ParameterList> integrator_pl =
352  state_integrator_->getValidParameters();
353  Teuchos::RCP<const Teuchos::ParameterList> sensitivity_pl =
355  pl->setParameters(*integrator_pl);
356  pl->sublist("Sensitivities").setParameters(*sensitivity_pl);
357 
358  return pl;
359 }
360 
361 template<class Scalar>
362 Teuchos::RCP<Teuchos::ParameterList>
365 {
366  return state_integrator_->getNonconstParameterList();
367 }
368 
369 template <class Scalar>
370 Teuchos::RCP<AdjointSensitivityModelEvaluator<Scalar> >
373  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
374  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& adjoint_model,
375  const Teuchos::RCP<Teuchos::ParameterList>& inputPL)
376 {
377  using Teuchos::rcp;
378 
379  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
380  if (inputPL != Teuchos::null) {
381  *pl = inputPL->sublist("Sensitivities");
382  }
383  const Scalar tfinal = state_integrator_->getTimeStepControl()->getFinalTime();
385  model, adjoint_model, tfinal, true, pl));
386 }
387 
388 template<class Scalar>
389 void
392 {
393  using Teuchos::RCP;
394  using Teuchos::rcp;
395  using Teuchos::rcp_dynamic_cast;
396  using Teuchos::ParameterList;
397  using Thyra::VectorBase;
398  using Thyra::VectorSpaceBase;
399  using Thyra::assign;
400  using Thyra::defaultProductVector;
401  typedef Thyra::DefaultProductVector<Scalar> DPV;
402  typedef Thyra::DefaultProductVectorSpace<Scalar> DPVS;
403 
404  // Create combined solution histories, first for the states with zero
405  // adjoint and then for the adjoint with frozen states
406  RCP<ParameterList> shPL =
407  Teuchos::sublist(state_integrator_->getIntegratorParameterList(),
408  "Solution History", true);
409  solutionHistory_ = createSolutionHistoryPL<Scalar>(shPL);
410 
411  RCP<const VectorSpaceBase<Scalar> > x_space =
412  model_->get_x_space();
413  RCP<const VectorSpaceBase<Scalar> > adjoint_space =
414  sens_model_->get_x_space();
415  Teuchos::Array< RCP<const VectorSpaceBase<Scalar> > > spaces(2);
416  spaces[0] = x_space;
417  spaces[1] = adjoint_space;
418  RCP<const DPVS > prod_space = Thyra::productVectorSpace(spaces());
419  const Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
420 
421  RCP<const SolutionHistory<Scalar> > state_solution_history =
422  state_integrator_->getSolutionHistory();
423  int num_states = state_solution_history->getNumStates();
424  for (int i=0; i<num_states; ++i) {
425  RCP<const SolutionState<Scalar> > state = (*state_solution_history)[i];
426 
427  // X
428  RCP<DPV> x = defaultProductVector(prod_space);
429  assign(x->getNonconstVectorBlock(0).ptr(), *(state->getX()));
430  assign(x->getNonconstVectorBlock(1).ptr(), zero);
431  RCP<VectorBase<Scalar> > x_b = x;
432 
433  // X-Dot
434  RCP<VectorBase<Scalar> > x_dot_b;
435  if (state->getXDot() != Teuchos::null) {
436  RCP<DPV> x_dot = defaultProductVector(prod_space);
437  assign(x_dot->getNonconstVectorBlock(0).ptr(), *(state->getXDot()));
438  assign(x_dot->getNonconstVectorBlock(1).ptr(), zero);
439  x_dot_b = x_dot;
440  }
441 
442  // X-Dot-Dot
443  RCP<VectorBase<Scalar> > x_dot_dot_b;
444  if (state->getXDotDot() != Teuchos::null) {
445  RCP<DPV> x_dot_dot = defaultProductVector(prod_space);
446  assign(x_dot_dot->getNonconstVectorBlock(0).ptr(),*(state->getXDotDot()));
447  assign(x_dot_dot->getNonconstVectorBlock(1).ptr(), zero);
448  x_dot_dot_b = x_dot_dot;
449  }
450 
451  RCP<SolutionState<Scalar> > prod_state = state->clone();
452  prod_state->setX(x_b);
453  prod_state->setXDot(x_dot_b);
454  prod_state->setXDotDot(x_dot_dot_b);
455  solutionHistory_->addState(prod_state);
456  }
457 
458  RCP<const VectorBase<Scalar> > frozen_x =
459  state_solution_history->getCurrentState()->getX();
460  RCP<const VectorBase<Scalar> > frozen_x_dot =
461  state_solution_history->getCurrentState()->getXDot();
462  RCP<const VectorBase<Scalar> > frozen_x_dot_dot =
463  state_solution_history->getCurrentState()->getXDotDot();
464  RCP<const SolutionHistory<Scalar> > sens_solution_history =
465  sens_integrator_->getSolutionHistory();
466  num_states = sens_solution_history->getNumStates();
467  for (int i=0; i<num_states; ++i) {
468  RCP<const SolutionState<Scalar> > state = (*sens_solution_history)[i];
469 
470  // X
471  RCP<DPV> x = defaultProductVector(prod_space);
472  assign(x->getNonconstVectorBlock(0).ptr(), *frozen_x);
473  assign(x->getNonconstVectorBlock(1).ptr(), *(state->getX()));
474  RCP<VectorBase<Scalar> > x_b = x;
475 
476  // X-Dot
477  RCP<VectorBase<Scalar> > x_dot_b;
478  if (state->getXDot() != Teuchos::null) {
479  RCP<DPV> x_dot = defaultProductVector(prod_space);
480  assign(x_dot->getNonconstVectorBlock(0).ptr(), *frozen_x_dot);
481  assign(x_dot->getNonconstVectorBlock(1).ptr(), *(state->getXDot()));
482  x_dot_b = x_dot;
483  }
484 
485  // X-Dot-Dot
486  RCP<VectorBase<Scalar> > x_dot_dot_b;
487  if (state->getXDotDot() != Teuchos::null) {
488  RCP<DPV> x_dot_dot = defaultProductVector(prod_space);
489  assign(x_dot_dot->getNonconstVectorBlock(0).ptr(), *frozen_x_dot_dot);
490  assign(x_dot_dot->getNonconstVectorBlock(1).ptr(),*(state->getXDotDot()));
491  x_dot_dot_b = x_dot_dot;
492  }
493 
494  RCP<SolutionState<Scalar> > prod_state = state->clone();
495  prod_state->setX(x_b);
496  prod_state->setXDot(x_dot_b);
497  prod_state->setXDotDot(x_dot_dot_b);
498  solutionHistory_->addState(prod_state, false);
499  }
500 }
501 
503 template<class Scalar>
504 Teuchos::RCP<IntegratorPseudoTransientAdjointSensitivity<Scalar> >
506  Teuchos::RCP<Teuchos::ParameterList> pList,
507  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model)
508 {
509  Teuchos::RCP<IntegratorPseudoTransientAdjointSensitivity<Scalar> > integrator =
510  Teuchos::rcp(new IntegratorPseudoTransientAdjointSensitivity<Scalar>(pList, model));
511  return(integrator);
512 }
513 
515 template<class Scalar>
516 Teuchos::RCP<IntegratorPseudoTransientAdjointSensitivity<Scalar> >
518  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
519  std::string stepperType)
520 {
521  Teuchos::RCP<IntegratorPseudoTransientAdjointSensitivity<Scalar> > integrator =
522  Teuchos::rcp(new IntegratorPseudoTransientAdjointSensitivity<Scalar>(model, stepperType));
523  return(integrator);
524 }
525 
527 template<class Scalar>
528 Teuchos::RCP<IntegratorPseudoTransientAdjointSensitivity<Scalar> >
530  Teuchos::RCP<Teuchos::ParameterList> pList,
531  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
532  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& adjoint_model)
533 {
534  Teuchos::RCP<IntegratorPseudoTransientAdjointSensitivity<Scalar> > integrator =
535  Teuchos::rcp(new IntegratorPseudoTransientAdjointSensitivity<Scalar>(pList, model, adjoint_model));
536  return(integrator);
537 }
538 
540 template<class Scalar>
541 Teuchos::RCP<IntegratorPseudoTransientAdjointSensitivity<Scalar> >
543  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
544  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& adjoint_model,
545  std::string stepperType)
546 {
547  Teuchos::RCP<IntegratorPseudoTransientAdjointSensitivity<Scalar> > integrator =
548  Teuchos::rcp(new IntegratorPseudoTransientAdjointSensitivity<Scalar>(model, adjoint_model, stepperType));
549  return(integrator);
550 }
551 
553 template<class Scalar>
554 Teuchos::RCP<IntegratorPseudoTransientAdjointSensitivity<Scalar> >
556 {
557  Teuchos::RCP<IntegratorPseudoTransientAdjointSensitivity<Scalar> > integrator =
559  return(integrator);
560 }
561 
562 } // namespace Tempus
563 #endif // Tempus_IntegratorPseudoTransientAdjointSensitivity_impl_hpp
static Teuchos::RCP< const Teuchos::ParameterList > getValidParameters()
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &pl) override
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
virtual void setTempusParameterList(Teuchos::RCP< Teuchos::ParameterList > pl) override
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXDotDot() const
Get current the second time derivative of the solution, xdotdot.
Teuchos::RCP< Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar > > integratorPseudoTransientAdjointSensitivity(Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model)
Nonmember constructor.
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDgDp() const
Return adjoint sensitivity stored in gradient format.
virtual Teuchos::RCP< SolutionHistory< Scalar > > getNonConstSolutionHistory() override
Get the SolutionHistory.
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
Time integrator suitable for pseudotransient adjoint sensitivity analysis.
Teuchos::RCP< AdjointSensitivityModelEvaluator< Scalar > > createSensitivityModel(const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &adjoint_model, const Teuchos::RCP< Teuchos::ParameterList > &inputPL)
virtual Teuchos::RCP< TimeStepControl< Scalar > > getNonConstTimeStepControl() override
Status
Status for the Integrator, the Stepper and the SolutionState.
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 > > y0=Teuchos::null, Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > ydot0=Teuchos::null, Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > ydotdot0=Teuchos::null)
Set the initial state from Thyra::VectorBase(s)
RCP< ImplicitAdjointModelEvaluator< Scalar > > implicitAdjointModelEvaluator(const RCP< const ModelEvaluator< Scalar > > &model)
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXDot() const
Get current the time derivative of the solution, xdot.
virtual Teuchos::RCP< Teuchos::ParameterList > getTempusParameterList() override
Return a copy of the Tempus ParameterList.
virtual Teuchos::RCP< Stepper< Scalar > > getStepper() const override
Get the Stepper.
virtual Teuchos::RCP< const SolutionHistory< Scalar > > getSolutionHistory() const override
Get the SolutionHistory.
virtual Teuchos::RCP< const TimeStepControl< Scalar > > getTimeStepControl() const override
Get the TimeStepControl.
ModelEvaluator for forming adjoint sensitivity equations.
virtual bool advanceTime()
Advance the solution to timeMax, and return true if successful.
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getX() const
Get current the solution, x.