Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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 "Teuchos_VerboseObjectParameterListHelpers.hpp"
13 #include "Thyra_DefaultMultiVectorProductVector.hpp"
14 #include "Thyra_VectorStdOps.hpp"
15 #include "Thyra_MultiVectorStdOps.hpp"
16 #include "NOX_Thyra.H"
17 
18 namespace Tempus {
19 
20 template<class Scalar>
23  Teuchos::RCP<Teuchos::ParameterList> inputPL,
24  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model)
25 {
26  model_ = model;
27  state_integrator_ = integratorBasic<Scalar>(inputPL, model_);
28  sens_model_ = createSensitivityModel(model_, inputPL);
29  sens_integrator_ = integratorBasic<Scalar>(inputPL, sens_model_);
30 }
31 
32 template<class Scalar>
35  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
36  std::string stepperType)
37 {
38  model_ = model;
39  state_integrator_ = integratorBasic<Scalar>(model_, stepperType);
40  sens_model_ = createSensitivityModel(model, Teuchos::null);
41  sens_integrator_ = integratorBasic<Scalar>(sens_model_, stepperType);
42 }
43 
44 template<class Scalar>
47 {
48  state_integrator_ = integratorBasic<Scalar>();
49  sens_integrator_ = integratorBasic<Scalar>();
50 }
51 
52 template<class Scalar>
53 bool
56 {
57  const Scalar tfinal =
58  state_integrator_->getTimeStepControl()->getFinalTime();
59  return advanceTime(tfinal);
60 }
61 
62 template<class Scalar>
63 bool
65 advanceTime(const Scalar timeFinal)
66 {
67  using Teuchos::RCP;
68  using Thyra::VectorBase;
69  typedef Thyra::ModelEvaluatorBase MEB;
70 
71  // Run state integrator and get solution
72  bool state_status = state_integrator_->advanceTime(timeFinal);
73 
74  // Set solution in sensitivity ME
75  sens_model_->setForwardSolutionHistory(
76  state_integrator_->getSolutionHistory());
77 
78  // Run sensitivity integrator
79  bool sens_status = sens_integrator_->advanceTime(timeFinal);
80 
81  // Compute final dg/dp which is computed by response 0 of the adjoint
82  // model evaluator
83  MEB::InArgs<Scalar> inargs = sens_model_->getNominalValues();
84  MEB::OutArgs<Scalar> outargs = sens_model_->createOutArgs();
85  inargs.set_t(sens_integrator_->getTime());
86  inargs.set_x(sens_integrator_->getX());
87  if (inargs.supports(MEB::IN_ARG_x_dot))
88  inargs.set_x_dot(sens_integrator_->getXdot());
89  if (inargs.supports(MEB::IN_ARG_x_dot_dot))
90  inargs.set_x_dot_dot(sens_integrator_->getXdotdot());
91  RCP<VectorBase<Scalar> > G = dgdp_;
92  if (G == Teuchos::null) {
93  G = Thyra::createMember(sens_model_->get_g_space(0));
94  dgdp_ = Teuchos::rcp_dynamic_cast<DMVPV>(G);
95  }
96  outargs.set_g(0, G);
97  sens_model_->evalModel(inargs, outargs);
98 
99  buildSolutionHistory();
100 
101  return state_status && sens_status;
102 }
103 
104 template<class Scalar>
105 Scalar
107 getTime() const
108 {
109  return solutionHistory_->getCurrentTime();
110 }
111 
112 template<class Scalar>
113 Scalar
115 getIndex() const
116 {
117  return solutionHistory_->getCurrentIndex();
118 }
119 
120 template<class Scalar>
121 Status
123 getStatus() const
124 {
125  Status state_status = state_integrator_->getStatus();
126  Status sens_status = sens_integrator_->getStatus();
127  if (state_status == FAILED || sens_status == FAILED)
128  return FAILED;
129  if (state_status == WORKING || sens_status == WORKING)
130  return WORKING;
131  return PASSED;
132 }
133 
134 template<class Scalar>
135 Teuchos::RCP<Stepper<Scalar> >
137 getStepper() const
138 {
139  return state_integrator_->getStepper();
140 }
141 
142 template<class Scalar>
143 Teuchos::RCP<Teuchos::ParameterList>
146 {
147  return state_integrator_->getTempusParameterList();
148 }
149 
150 template<class Scalar>
151 void
153 setTempusParameterList(Teuchos::RCP<Teuchos::ParameterList> pl)
154 {
155  state_integrator_->setTempusParameterList(pl);
156  sens_integrator_->setTempusParameterList(pl);
157 }
158 
159 template<class Scalar>
160 Teuchos::RCP<const SolutionHistory<Scalar> >
163 {
164  return solutionHistory_;
165 }
166 
167 template<class Scalar>
168 Teuchos::RCP<const TimeStepControl<Scalar> >
171 {
172  return state_integrator_->getTimeStepControl();
173 }
174 
175 template<class Scalar>
178  Teuchos::RCP<const Thyra::VectorBase<Scalar> > x0,
179  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xdot0,
180  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xdotdot0,
181  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > y0,
182  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > ydot0,
183  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > ydotdot0)
184 {
185  using Teuchos::RCP;
186  using Teuchos::rcp_dynamic_cast;
187  using Thyra::VectorSpaceBase;
188  using Thyra::assign;
189  using Thyra::createMember;
190 
191  //
192  // Create and initialize product X, Xdot, Xdotdot
193 
194  RCP< const VectorSpaceBase<Scalar> > space = sens_model_->get_x_space();
195  RCP<DMVPV> Y = rcp_dynamic_cast<DMVPV>(createMember(space));
196  RCP<DMVPV> Ydot = rcp_dynamic_cast<DMVPV>(createMember(space));
197  RCP<DMVPV> Ydotdot = rcp_dynamic_cast<DMVPV>(createMember(space));
198  const Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
199 
200  // y
201  if (y0 == Teuchos::null)
202  assign(Y->getNonconstMultiVector().ptr(), zero);
203  else
204  assign(Y->getNonconstMultiVector().ptr(), *y0);
205 
206  // ydot
207  if (ydot0 == Teuchos::null)
208  assign(Ydot->getNonconstMultiVector().ptr(), zero);
209  else
210  assign(Ydot->getNonconstMultiVector().ptr(), *ydot0);
211 
212  // ydotdot
213  if (ydotdot0 == Teuchos::null)
214  assign(Ydotdot->getNonconstMultiVector().ptr(), zero);
215  else
216  assign(Ydotdot->getNonconstMultiVector().ptr(), *ydotdot0);
217 
218  state_integrator_->initializeSolutionHistory(t0, x0, xdot0, xdotdot0);
219  sens_integrator_->initializeSolutionHistory(t0, Y, Ydot, Ydotdot);
220 }
221 
222 template<class Scalar>
223 Teuchos::RCP<const Thyra::VectorBase<Scalar> >
225 getX() const
226 {
227  return state_integrator_->getX();
228 }
229 
230 template<class Scalar>
231 Teuchos::RCP<const Thyra::VectorBase<Scalar> >
233 getXdot() const
234 {
235  return state_integrator_->getXdot();
236 }
237 
238 template<class Scalar>
239 Teuchos::RCP<const Thyra::VectorBase<Scalar> >
241 getXdotdot() const
242 {
243  return state_integrator_->getXdotdot();
244 }
245 
246 template<class Scalar>
247 Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> >
249 getDgDp() const
250 {
251  return dgdp_->getMultiVector();
252 }
253 
254 template<class Scalar>
255 std::string
257 description() const
258 {
259  std::string name = "Tempus::IntegratorPseudoTransientAdjointSensitivity";
260  return(name);
261 }
262 
263 template<class Scalar>
264 void
267  Teuchos::FancyOStream &out,
268  const Teuchos::EVerbosityLevel verbLevel) const
269 {
270  out << description() << "::describe" << std::endl;
271  state_integrator_->describe(out, verbLevel);
272  sens_integrator_->describe(out, verbLevel);
273 }
274 
275 template<class Scalar>
276 void
278 setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & inputPL)
279 {
280  state_integrator_->setParameterList(inputPL);
281  sens_integrator_->setParameterList(inputPL);
282 }
283 
284 template<class Scalar>
285 Teuchos::RCP<Teuchos::ParameterList>
288 {
289  state_integrator_->unsetParameterList();
290  return sens_integrator_->unsetParameterList();
291 }
292 
293 template<class Scalar>
294 Teuchos::RCP<const Teuchos::ParameterList>
297 {
298  Teuchos::RCP<Teuchos::ParameterList> pl =
299  Teuchos::rcp(new Teuchos::ParameterList);
300  Teuchos::RCP<const Teuchos::ParameterList> integrator_pl =
301  state_integrator_->getValidParameters();
302  Teuchos::RCP<const Teuchos::ParameterList> sensitivity_pl =
304  pl->setParameters(*integrator_pl);
305  pl->sublist("Sensitivities").setParameters(*sensitivity_pl);
306 
307  return pl;
308 }
309 
310 template<class Scalar>
311 Teuchos::RCP<Teuchos::ParameterList>
314 {
315  return state_integrator_->getNonconstParameterList();
316 }
317 
318 template <class Scalar>
319 Teuchos::RCP<Tempus::AdjointSensitivityModelEvaluator<Scalar> >
322  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
323  const Teuchos::RCP<Teuchos::ParameterList>& inputPL)
324 {
325  using Teuchos::rcp;
326 
327  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
328  if (inputPL != Teuchos::null) {
329  *pl = inputPL->sublist("Sensitivities");
330  }
331  const Scalar tfinal = state_integrator_->getTimeStepControl()->getFinalTime();
333  model, tfinal, true, pl));
334 }
335 
336 template<class Scalar>
337 void
340 {
341  using Teuchos::RCP;
342  using Teuchos::rcp;
343  using Teuchos::rcp_dynamic_cast;
344  using Teuchos::ParameterList;
345  using Thyra::VectorBase;
346  using Thyra::VectorSpaceBase;
347  using Thyra::assign;
348  using Thyra::defaultProductVector;
349  typedef Thyra::DefaultProductVector<Scalar> DPV;
350  typedef Thyra::DefaultProductVectorSpace<Scalar> DPVS;
351 
352  // Create combined solution histories, first for the states with zero
353  // adjoint and then for the adjoint with frozen states
354  RCP<ParameterList> shPL =
355  Teuchos::sublist(state_integrator_->getIntegratorParameterList(),
356  "Solution History", true);
357  solutionHistory_ = rcp(new SolutionHistory<Scalar>(shPL));
358 
359  RCP<const VectorSpaceBase<Scalar> > x_space =
360  model_->get_x_space();
361  RCP<const VectorSpaceBase<Scalar> > adjoint_space =
362  sens_model_->get_x_space();
363  Teuchos::Array< RCP<const VectorSpaceBase<Scalar> > > spaces(2);
364  spaces[0] = x_space;
365  spaces[1] = adjoint_space;
366  RCP<const DPVS > prod_space = Thyra::productVectorSpace(spaces());
367  const Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
368 
369  RCP<const SolutionHistory<Scalar> > state_solution_history =
370  state_integrator_->getSolutionHistory();
371  int num_states = state_solution_history->getNumStates();
372  for (int i=0; i<num_states; ++i) {
373  RCP<const SolutionState<Scalar> > state = (*state_solution_history)[i];
374 
375  // X
376  RCP<DPV> x = defaultProductVector(prod_space);
377  assign(x->getNonconstVectorBlock(0).ptr(), *(state->getX()));
378  assign(x->getNonconstVectorBlock(1).ptr(), zero);
379  RCP<VectorBase<Scalar> > x_b = x;
380 
381  // X-Dot
382  RCP<VectorBase<Scalar> > x_dot_b;
383  if (state->getXDot() != Teuchos::null) {
384  RCP<DPV> x_dot = defaultProductVector(prod_space);
385  assign(x_dot->getNonconstVectorBlock(0).ptr(), *(state->getXDot()));
386  assign(x_dot->getNonconstVectorBlock(1).ptr(), zero);
387  x_dot_b = x_dot;
388  }
389 
390  // X-Dot-Dot
391  RCP<VectorBase<Scalar> > x_dot_dot_b;
392  if (state->getXDotDot() != Teuchos::null) {
393  RCP<DPV> x_dot_dot = defaultProductVector(prod_space);
394  assign(x_dot_dot->getNonconstVectorBlock(0).ptr(),*(state->getXDotDot()));
395  assign(x_dot_dot->getNonconstVectorBlock(1).ptr(), zero);
396  x_dot_dot_b = x_dot_dot;
397  }
398 
399  RCP<SolutionState<Scalar> > prod_state =
400  rcp(new SolutionState<Scalar>(state->getMetaData()->clone(),
401  x_b, x_dot_b, x_dot_dot_b,
402  state->getStepperState()->clone()));
403  solutionHistory_->addState(prod_state);
404  }
405 
406  RCP<const VectorBase<Scalar> > frozen_x =
407  state_solution_history->getCurrentState()->getX();
408  RCP<const VectorBase<Scalar> > frozen_x_dot =
409  state_solution_history->getCurrentState()->getXDot();
410  RCP<const VectorBase<Scalar> > frozen_x_dot_dot =
411  state_solution_history->getCurrentState()->getXDotDot();
412  RCP<const SolutionHistory<Scalar> > sens_solution_history =
413  sens_integrator_->getSolutionHistory();
414  num_states = sens_solution_history->getNumStates();
415  for (int i=0; i<num_states; ++i) {
416  RCP<const SolutionState<Scalar> > state = (*sens_solution_history)[i];
417 
418  // X
419  RCP<DPV> x = defaultProductVector(prod_space);
420  assign(x->getNonconstVectorBlock(0).ptr(), *frozen_x);
421  assign(x->getNonconstVectorBlock(1).ptr(), *(state->getX()));
422  RCP<VectorBase<Scalar> > x_b = x;
423 
424  // X-Dot
425  RCP<VectorBase<Scalar> > x_dot_b;
426  if (state->getXDot() != Teuchos::null) {
427  RCP<DPV> x_dot = defaultProductVector(prod_space);
428  assign(x_dot->getNonconstVectorBlock(0).ptr(), *frozen_x_dot);
429  assign(x_dot->getNonconstVectorBlock(1).ptr(), *(state->getXDot()));
430  x_dot_b = x_dot;
431  }
432 
433  // X-Dot-Dot
434  RCP<VectorBase<Scalar> > x_dot_dot_b;
435  if (state->getXDotDot() != Teuchos::null) {
436  RCP<DPV> x_dot_dot = defaultProductVector(prod_space);
437  assign(x_dot_dot->getNonconstVectorBlock(0).ptr(), *frozen_x_dot_dot);
438  assign(x_dot_dot->getNonconstVectorBlock(1).ptr(),*(state->getXDotDot()));
439  x_dot_dot_b = x_dot_dot;
440  }
441 
442  RCP<SolutionState<Scalar> > prod_state =
443  rcp(new SolutionState<Scalar>(state->getMetaData()->clone(),
444  x_b, x_dot_b, x_dot_dot_b,
445  state->getStepperState()->clone()));
446  solutionHistory_->addState(prod_state);
447  }
448 }
449 
450 /// Non-member constructor
451 template<class Scalar>
452 Teuchos::RCP<Tempus::IntegratorPseudoTransientAdjointSensitivity<Scalar> >
454  Teuchos::RCP<Teuchos::ParameterList> pList,
455  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model)
456 {
457  Teuchos::RCP<Tempus::IntegratorPseudoTransientAdjointSensitivity<Scalar> > integrator =
458  Teuchos::rcp(new Tempus::IntegratorPseudoTransientAdjointSensitivity<Scalar>(pList, model));
459  return(integrator);
460 }
461 
462 /// Non-member constructor
463 template<class Scalar>
464 Teuchos::RCP<Tempus::IntegratorPseudoTransientAdjointSensitivity<Scalar> >
466  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
467  std::string stepperType)
468 {
469  Teuchos::RCP<Tempus::IntegratorPseudoTransientAdjointSensitivity<Scalar> > integrator =
470  Teuchos::rcp(new Tempus::IntegratorPseudoTransientAdjointSensitivity<Scalar>(model, stepperType));
471  return(integrator);
472 }
473 
474 /// Non-member constructor
475 template<class Scalar>
476 Teuchos::RCP<Tempus::IntegratorPseudoTransientAdjointSensitivity<Scalar> >
478 {
479  Teuchos::RCP<Tempus::IntegratorPseudoTransientAdjointSensitivity<Scalar> > integrator =
481  return(integrator);
482 }
483 
484 } // namespace Tempus
485 #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
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDgDp() const
Return adjoint sensitivity stored in gradient format.
Time integrator suitable for pseudotransient adjoint sensitivity analysis.
Teuchos::RCP< Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar > > integratorPseudoTransientAdjointSensitivity(Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model)
Non-member constructor.
Status
Status for the Integrator, the Stepper and the SolutionState.
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXdot() const
Get current the time derivative of the solution, xdot.
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)
Teuchos::RCP< AdjointSensitivityModelEvaluator< Scalar > > createSensitivityModel(const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< Teuchos::ParameterList > &inputPL)
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXdotdot() const
Get current the second time derivative of the solution, xdotdot.
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.
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...