Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups 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 "Thyra_DefaultMultiVectorProductVector.hpp"
13 #include "Thyra_VectorStdOps.hpp"
14 #include "Thyra_MultiVectorStdOps.hpp"
15 
16 
17 namespace Tempus {
18 
19 template<class Scalar>
24 {
25  model_ = model;
26  state_integrator_ = integratorBasic<Scalar>(inputPL, model_);
27  sens_model_ = createSensitivityModel(model_, inputPL);
28  sens_integrator_ = integratorBasic<Scalar>(inputPL, sens_model_);
29 }
30 
31 template<class Scalar>
35  std::string stepperType)
36 {
37  model_ = model;
38  state_integrator_ = integratorBasic<Scalar>(model_, stepperType);
39  sens_model_ = createSensitivityModel(model, Teuchos::null);
40  sens_integrator_ = integratorBasic<Scalar>(sens_model_, stepperType);
41 }
42 
43 template<class Scalar>
46 {
47  state_integrator_ = integratorBasic<Scalar>();
48  sens_integrator_ = integratorBasic<Scalar>();
49 }
50 
51 template<class Scalar>
52 bool
55 {
56  const Scalar tfinal =
57  state_integrator_->getTimeStepControl()->getFinalTime();
58  return advanceTime(tfinal);
59 }
60 
61 template<class Scalar>
62 bool
64 advanceTime(const Scalar timeFinal)
65 {
66  using Teuchos::RCP;
67  using Thyra::VectorBase;
68  typedef Thyra::ModelEvaluatorBase MEB;
69 
70  // Run state integrator and get solution
71  bool state_status = state_integrator_->advanceTime(timeFinal);
72 
73  // Set solution in sensitivity ME
74  sens_model_->setForwardSolutionHistory(
75  state_integrator_->getSolutionHistory());
76 
77  // Run sensitivity integrator
78  bool sens_status = sens_integrator_->advanceTime(timeFinal);
79 
80  // Compute final dg/dp which is computed by response 0 of the adjoint
81  // model evaluator
82  MEB::InArgs<Scalar> inargs = sens_model_->getNominalValues();
83  MEB::OutArgs<Scalar> outargs = sens_model_->createOutArgs();
84  inargs.set_t(sens_integrator_->getTime());
85  inargs.set_x(sens_integrator_->getX());
86  if (inargs.supports(MEB::IN_ARG_x_dot))
87  inargs.set_x_dot(sens_integrator_->getXDot());
88  if (inargs.supports(MEB::IN_ARG_x_dot_dot))
89  inargs.set_x_dot_dot(sens_integrator_->getXDotDot());
90  RCP<VectorBase<Scalar> > G = dgdp_;
91  if (G == Teuchos::null) {
92  G = Thyra::createMember(sens_model_->get_g_space(0));
93  dgdp_ = Teuchos::rcp_dynamic_cast<DMVPV>(G);
94  }
95  outargs.set_g(0, G);
96  sens_model_->evalModel(inargs, outargs);
97 
98  buildSolutionHistory();
99 
100  return state_status && sens_status;
101 }
102 
103 template<class Scalar>
104 Scalar
106 getTime() const
107 {
108  return solutionHistory_->getCurrentTime();
109 }
110 
111 template<class Scalar>
112 int
114 getIndex() const
115 {
116  return solutionHistory_->getCurrentIndex();
117 }
118 
119 template<class Scalar>
120 Status
122 getStatus() const
123 {
124  Status state_status = state_integrator_->getStatus();
125  Status sens_status = sens_integrator_->getStatus();
126  if (state_status == FAILED || sens_status == FAILED)
127  return FAILED;
128  if (state_status == WORKING || sens_status == WORKING)
129  return WORKING;
130  return PASSED;
131 }
132 
133 template<class Scalar>
136 getStepper() const
137 {
138  return state_integrator_->getStepper();
139 }
140 
141 template<class Scalar>
145 {
146  return state_integrator_->getTempusParameterList();
147 }
148 
149 template<class Scalar>
150 void
153 {
154  state_integrator_->setTempusParameterList(pl);
155  sens_integrator_->setTempusParameterList(pl);
156 }
157 
158 template<class Scalar>
162 {
163  return solutionHistory_;
164 }
165 
166 template<class Scalar>
170 {
171  return state_integrator_->getTimeStepControl();
172 }
173 
174 template<class Scalar>
178 {
179  return state_integrator_->getNonConstTimeStepControl();
180 }
181 
182 template<class Scalar>
187  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xdotdot0,
191 {
192  using Teuchos::RCP;
193  using Teuchos::rcp_dynamic_cast;
195  using Thyra::assign;
196  using Thyra::createMember;
197 
198  //
199  // Create and initialize product X, Xdot, Xdotdot
200 
201  RCP< const VectorSpaceBase<Scalar> > space = sens_model_->get_x_space();
202  RCP<DMVPV> Y = rcp_dynamic_cast<DMVPV>(createMember(space));
203  RCP<DMVPV> Ydot = rcp_dynamic_cast<DMVPV>(createMember(space));
204  RCP<DMVPV> Ydotdot = rcp_dynamic_cast<DMVPV>(createMember(space));
205  const Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
206 
207  // y
208  if (y0 == Teuchos::null)
209  assign(Y->getNonconstMultiVector().ptr(), zero);
210  else
211  assign(Y->getNonconstMultiVector().ptr(), *y0);
212 
213  // ydot
214  if (ydot0 == Teuchos::null)
215  assign(Ydot->getNonconstMultiVector().ptr(), zero);
216  else
217  assign(Ydot->getNonconstMultiVector().ptr(), *ydot0);
218 
219  // ydotdot
220  if (ydotdot0 == Teuchos::null)
221  assign(Ydotdot->getNonconstMultiVector().ptr(), zero);
222  else
223  assign(Ydotdot->getNonconstMultiVector().ptr(), *ydotdot0);
224 
225  state_integrator_->initializeSolutionHistory(t0, x0, xdot0, xdotdot0);
226  sens_integrator_->initializeSolutionHistory(t0, Y, Ydot, Ydotdot);
227 }
228 
229 template<class Scalar>
232 getX() const
233 {
234  return state_integrator_->getX();
235 }
236 
237 template<class Scalar>
240 getXDot() const
241 {
242  return state_integrator_->getXDot();
243 }
244 
245 template<class Scalar>
248 getXDotDot() const
249 {
250  return state_integrator_->getXDotDot();
251 }
252 
253 template<class Scalar>
256 getDgDp() const
257 {
258  return dgdp_->getMultiVector();
259 }
260 
261 template<class Scalar>
262 std::string
264 description() const
265 {
266  std::string name = "Tempus::IntegratorPseudoTransientAdjointSensitivity";
267  return(name);
268 }
269 
270 template<class Scalar>
271 void
274  Teuchos::FancyOStream &in_out,
275  const Teuchos::EVerbosityLevel verbLevel) const
276 {
277  auto out = Teuchos::fancyOStream( in_out.getOStream() );
278  *out << description() << "::describe" << std::endl;
279  state_integrator_->describe(*out, verbLevel);
280  sens_integrator_->describe(*out, verbLevel);
281 }
282 
283 template<class Scalar>
284 void
287 {
288  state_integrator_->setParameterList(inputPL);
289  sens_integrator_->setParameterList(inputPL);
290 }
291 
292 template<class Scalar>
296 {
297  state_integrator_->unsetParameterList();
298  return sens_integrator_->unsetParameterList();
299 }
300 
301 template<class Scalar>
305 {
309  state_integrator_->getValidParameters();
312  pl->setParameters(*integrator_pl);
313  pl->sublist("Sensitivities").setParameters(*sensitivity_pl);
314 
315  return pl;
316 }
317 
318 template<class Scalar>
322 {
323  return state_integrator_->getNonconstParameterList();
324 }
325 
326 template <class Scalar>
332 {
333  using Teuchos::rcp;
334 
335  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
336  if (inputPL != Teuchos::null) {
337  *pl = inputPL->sublist("Sensitivities");
338  }
339  const Scalar tfinal = state_integrator_->getTimeStepControl()->getFinalTime();
341  model, tfinal, true, pl));
342 }
343 
344 template<class Scalar>
345 void
348 {
349  using Teuchos::RCP;
350  using Teuchos::rcp;
351  using Teuchos::rcp_dynamic_cast;
353  using Thyra::VectorBase;
355  using Thyra::assign;
356  using Thyra::defaultProductVector;
359 
360  // Create combined solution histories, first for the states with zero
361  // adjoint and then for the adjoint with frozen states
362  RCP<ParameterList> shPL =
363  Teuchos::sublist(state_integrator_->getIntegratorParameterList(),
364  "Solution History", true);
365  solutionHistory_ = createSolutionHistoryPL<Scalar>(shPL);
366 
367  RCP<const VectorSpaceBase<Scalar> > x_space =
368  model_->get_x_space();
369  RCP<const VectorSpaceBase<Scalar> > adjoint_space =
370  sens_model_->get_x_space();
372  spaces[0] = x_space;
373  spaces[1] = adjoint_space;
374  RCP<const DPVS > prod_space = Thyra::productVectorSpace(spaces());
375  const Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
376 
377  RCP<const SolutionHistory<Scalar> > state_solution_history =
378  state_integrator_->getSolutionHistory();
379  int num_states = state_solution_history->getNumStates();
380  for (int i=0; i<num_states; ++i) {
381  RCP<const SolutionState<Scalar> > state = (*state_solution_history)[i];
382 
383  // X
384  RCP<DPV> x = defaultProductVector(prod_space);
385  assign(x->getNonconstVectorBlock(0).ptr(), *(state->getX()));
386  assign(x->getNonconstVectorBlock(1).ptr(), zero);
387  RCP<VectorBase<Scalar> > x_b = x;
388 
389  // X-Dot
390  RCP<VectorBase<Scalar> > x_dot_b;
391  if (state->getXDot() != Teuchos::null) {
392  RCP<DPV> x_dot = defaultProductVector(prod_space);
393  assign(x_dot->getNonconstVectorBlock(0).ptr(), *(state->getXDot()));
394  assign(x_dot->getNonconstVectorBlock(1).ptr(), zero);
395  x_dot_b = x_dot;
396  }
397 
398  // X-Dot-Dot
399  RCP<VectorBase<Scalar> > x_dot_dot_b;
400  if (state->getXDotDot() != Teuchos::null) {
401  RCP<DPV> x_dot_dot = defaultProductVector(prod_space);
402  assign(x_dot_dot->getNonconstVectorBlock(0).ptr(),*(state->getXDotDot()));
403  assign(x_dot_dot->getNonconstVectorBlock(1).ptr(), zero);
404  x_dot_dot_b = x_dot_dot;
405  }
406 
407  RCP<SolutionState<Scalar> > prod_state = state->clone();
408  prod_state->setX(x_b);
409  prod_state->setXDot(x_dot_b);
410  prod_state->setXDotDot(x_dot_dot_b);
411  solutionHistory_->addState(prod_state);
412  }
413 
414  RCP<const VectorBase<Scalar> > frozen_x =
415  state_solution_history->getCurrentState()->getX();
416  RCP<const VectorBase<Scalar> > frozen_x_dot =
417  state_solution_history->getCurrentState()->getXDot();
418  RCP<const VectorBase<Scalar> > frozen_x_dot_dot =
419  state_solution_history->getCurrentState()->getXDotDot();
420  RCP<const SolutionHistory<Scalar> > sens_solution_history =
421  sens_integrator_->getSolutionHistory();
422  num_states = sens_solution_history->getNumStates();
423  for (int i=0; i<num_states; ++i) {
424  RCP<const SolutionState<Scalar> > state = (*sens_solution_history)[i];
425 
426  // X
427  RCP<DPV> x = defaultProductVector(prod_space);
428  assign(x->getNonconstVectorBlock(0).ptr(), *frozen_x);
429  assign(x->getNonconstVectorBlock(1).ptr(), *(state->getX()));
430  RCP<VectorBase<Scalar> > x_b = x;
431 
432  // X-Dot
433  RCP<VectorBase<Scalar> > x_dot_b;
434  if (state->getXDot() != Teuchos::null) {
435  RCP<DPV> x_dot = defaultProductVector(prod_space);
436  assign(x_dot->getNonconstVectorBlock(0).ptr(), *frozen_x_dot);
437  assign(x_dot->getNonconstVectorBlock(1).ptr(), *(state->getXDot()));
438  x_dot_b = x_dot;
439  }
440 
441  // X-Dot-Dot
442  RCP<VectorBase<Scalar> > x_dot_dot_b;
443  if (state->getXDotDot() != Teuchos::null) {
444  RCP<DPV> x_dot_dot = defaultProductVector(prod_space);
445  assign(x_dot_dot->getNonconstVectorBlock(0).ptr(), *frozen_x_dot_dot);
446  assign(x_dot_dot->getNonconstVectorBlock(1).ptr(),*(state->getXDotDot()));
447  x_dot_dot_b = x_dot_dot;
448  }
449 
450  RCP<SolutionState<Scalar> > prod_state = state->clone();
451  prod_state->setX(x_b);
452  prod_state->setXDot(x_dot_b);
453  prod_state->setXDotDot(x_dot_dot_b);
454  solutionHistory_->addState(prod_state, false);
455  }
456 }
457 
459 template<class Scalar>
464 {
467  return(integrator);
468 }
469 
471 template<class Scalar>
475  std::string stepperType)
476 {
479  return(integrator);
480 }
481 
483 template<class Scalar>
486 {
489  return(integrator);
490 }
491 
492 } // namespace Tempus
493 #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 > > getXDot() const
Get current the time derivative of the solution, xdot.
Teuchos::RCP< Tempus::IntegratorPseudoTransientAdjointSensitivity< Scalar > > integratorPseudoTransientAdjointSensitivity(Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model)
Nonmember constructor.
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_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual Teuchos::RCP< TimeStepControl< Scalar > > getNonConstTimeStepControl() override
Status
Status for the Integrator, the Stepper and the SolutionState.
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXDotDot() const
Get current the second time derivative of the solution, xdotdot.
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)
ParameterList & setParameters(const ParameterList &source)
Teuchos::RCP< AdjointSensitivityModelEvaluator< Scalar > > createSensitivityModel(const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< Teuchos::ParameterList > &inputPL)
virtual Teuchos::RCP< Teuchos::ParameterList > getTempusParameterList() override
Return a copy of the Tempus ParameterList.
RCP< std::basic_ostream< char_type, traits_type > > getOStream()
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.
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
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.