Tempus  Version of the Day
Time Integration
Tempus_IntegratorAdjointSensitivity_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_IntegratorAdjointSensitivity_impl_hpp
10 #define Tempus_IntegratorAdjointSensitivity_impl_hpp
11 
12 #include "Teuchos_VerboseObjectParameterListHelpers.hpp"
13 #include "Thyra_DefaultMultiVectorProductVectorSpace.hpp"
14 #include "Thyra_DefaultMultiVectorProductVector.hpp"
15 #include "Thyra_DefaultProductVectorSpace.hpp"
16 #include "Thyra_DefaultProductVector.hpp"
17 #include "Thyra_VectorStdOps.hpp"
18 #include "Thyra_MultiVectorStdOps.hpp"
19 #include "NOX_Thyra.H"
20 
21 namespace Tempus {
22 
23 template<class Scalar>
26  Teuchos::RCP<Teuchos::ParameterList> inputPL,
27  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model)
28 {
29  model_ = model;
30  setParameterList(inputPL);
31  state_integrator_ = integratorBasic<Scalar>(inputPL, model_);
32  adjoint_model_ = createAdjointModel(model_, inputPL);
33  adjoint_integrator_ = integratorBasic<Scalar>(inputPL, adjoint_model_);
34 }
35 
36 template<class Scalar>
39 {
40  state_integrator_ = integratorBasic<Scalar>();
41  adjoint_integrator_ = integratorBasic<Scalar>();
42 }
43 
44 template<class Scalar>
45 bool
48 {
49  const Scalar tfinal =
50  state_integrator_->getTimeStepControl()->getFinalTime();
51  return advanceTime(tfinal);
52 }
53 
54 template<class Scalar>
55 bool
57 advanceTime(const Scalar timeFinal)
58 {
59  using Teuchos::RCP;
60  using Teuchos::rcp_dynamic_cast;
61  using Thyra::VectorBase;
62  using Thyra::VectorSpaceBase;
63  using Thyra::MultiVectorBase;
64  using Thyra::LinearOpBase;
65  using Thyra::LinearOpWithSolveBase;
66  using Thyra::createMember;
67  using Thyra::createMembers;
68  using Thyra::assign;
69  typedef Thyra::ModelEvaluatorBase MEB;
70  typedef Thyra::DefaultMultiVectorProductVector<Scalar> DMVPV;
71  typedef Thyra::DefaultProductVector<Scalar> DPV;
72 
73  // Get initial state for later
74  RCP<const SolutionHistory<Scalar> > state_solution_history =
75  state_integrator_->getSolutionHistory();
76  RCP<const SolutionState<Scalar> > initial_state =
77  (*state_solution_history)[0];
78 
79  // Run state integrator and get solution
80  bool state_status = state_integrator_->advanceTime(timeFinal);
81 
82  // Set solution history in adjoint stepper
83  adjoint_model_->setForwardSolutionHistory(state_solution_history);
84 
85  // Compute dg/dx
86  RCP<const VectorSpaceBase<Scalar> > g_space = model_->get_g_space(g_index_);
87  RCP<const VectorSpaceBase<Scalar> > x_space = model_->get_x_space();
88  const int num_g = g_space->dim();
89  RCP<MultiVectorBase<Scalar> > dgdx = createMembers(x_space, num_g);
90  MEB::InArgs<Scalar> inargs = model_->getNominalValues();
91  RCP<const SolutionState<Scalar> > state =
92  state_solution_history->getCurrentState();
93  inargs.set_t(state->getTime());
94  inargs.set_x(state->getX());
95  inargs.set_x_dot(state->getXDot());
96  MEB::OutArgs<Scalar> outargs = model_->createOutArgs();
97  outargs.set_DgDx(g_index_,
98  MEB::Derivative<Scalar>(dgdx, MEB::DERIV_MV_GRADIENT_FORM));
99  model_->evalModel(inargs, outargs);
100  outargs.set_DgDx(g_index_, MEB::Derivative<Scalar>());
101 
102  // Compute ICs == [ (df/dx_dot)^{-T} (dg/dx)^T; 0 ]
103  // For explicit form, we are relying on the user to inform us the
104  // the mass matrix is the identity. It would be nice to be able to determine
105  // somehow automatically that we are using an explicit stepper.
106  RCP<DPV> adjoint_init =
107  rcp_dynamic_cast<DPV>(Thyra::createMember(adjoint_model_->get_x_space()));
108  RCP<MultiVectorBase<Scalar> > adjoint_init_mv =
109  rcp_dynamic_cast<DMVPV>(adjoint_init->getNonconstVectorBlock(0))->getNonconstMultiVector();
110  assign(adjoint_init->getNonconstVectorBlock(1).ptr(),
111  Teuchos::ScalarTraits<Scalar>::zero());
112  if (mass_matrix_is_identity_)
113  assign(adjoint_init_mv.ptr(), *dgdx);
114  else {
115  inargs.set_alpha(1.0);
116  inargs.set_beta(0.0);
117  RCP<LinearOpWithSolveBase<Scalar> > W = model_->create_W();
118  outargs.set_W(W);
119  model_->evalModel(inargs, outargs);
120  W->solve(Thyra::CONJTRANS, *dgdx, adjoint_init_mv.ptr());
121  outargs.set_W(Teuchos::null);
122  }
123 
124  // Run sensitivity integrator and get solution
125  adjoint_integrator_->setInitialState(Scalar(0.0), adjoint_init);
126  bool sens_status = adjoint_integrator_->advanceTime(timeFinal);
127  RCP<const SolutionHistory<Scalar> > adjoint_solution_history =
128  adjoint_integrator_->getSolutionHistory();
129 
130  // Compute dg/dp at final time T
131  RCP<const VectorSpaceBase<Scalar> > p_space = model_->get_p_space(p_index_);
132  dgdp_ = createMembers(p_space, num_g);
133  if (g_depends_on_p_) {
134  MEB::DerivativeSupport dgdp_support =
135  outargs.supports(MEB::OUT_ARG_DgDp, g_index_, p_index_);
136  if (dgdp_support.supports(MEB::DERIV_MV_GRADIENT_FORM)) {
137  outargs.set_DgDp(g_index_, p_index_,
138  MEB::Derivative<Scalar>(dgdp_,
139  MEB::DERIV_MV_GRADIENT_FORM));
140  model_->evalModel(inargs, outargs);
141  }
142  else if (dgdp_support.supports(MEB::DERIV_MV_JACOBIAN_FORM)) {
143  const int num_p = p_space->dim();
144  RCP<MultiVectorBase<Scalar> > dgdp_trans =
145  createMembers(g_space, num_p);
146  outargs.set_DgDp(g_index_, p_index_,
147  MEB::Derivative<Scalar>(dgdp_trans,
148  MEB::DERIV_MV_JACOBIAN_FORM));
149  model_->evalModel(inargs, outargs);
150  Thyra::DetachedMultiVectorView<Scalar> dgdp_view(*dgdp_);
151  Thyra::DetachedMultiVectorView<Scalar> dgdp_trans_view(*dgdp_trans);
152  for (int i=0; i<num_p; ++i)
153  for (int j=0; j<num_g; ++j)
154  dgdp_view(i,j) = dgdp_trans_view(j,i);
155  }
156  else
157  TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
158  "Invalid dg/dp support");
159  outargs.set_DgDp(g_index_, p_index_, MEB::Derivative<Scalar>());
160  }
161  else
162  assign(dgdp_.ptr(), Scalar(0.0));
163 
164  // Add in initial condition term = (dx/dp^T(0))*(df/dx_dot^T(0))*y(0)
165  // If dxdp_init_ is null, assume it is zero
166  if (ic_depends_on_p_ && dxdp_init_ != Teuchos::null) {
167  RCP<const SolutionState<Scalar> > adjoint_state =
168  adjoint_solution_history->getCurrentState();
169  RCP<const VectorBase<Scalar> > adjoint_x =
170  rcp_dynamic_cast<const DPV>(adjoint_state->getX())->getVectorBlock(0);
171  RCP<const MultiVectorBase<Scalar> > adjoint_mv =
172  rcp_dynamic_cast<const DMVPV>(adjoint_x)->getMultiVector();
173  if (mass_matrix_is_identity_)
174  dxdp_init_->apply(Thyra::CONJTRANS, *adjoint_mv, dgdp_.ptr(), Scalar(1.0),
175  Scalar(1.0));
176  else {
177  inargs.set_t(initial_state->getTime());
178  inargs.set_x(initial_state->getX());
179  inargs.set_x_dot(initial_state->getXDot());
180  inargs.set_alpha(1.0);
181  inargs.set_beta(0.0);
182  RCP<LinearOpBase<Scalar> > W_op = model_->create_W_op();
183  outargs.set_W_op(W_op);
184  model_->evalModel(inargs, outargs);
185  outargs.set_W_op(Teuchos::null);
186  RCP<MultiVectorBase<Scalar> > tmp = createMembers(x_space, num_g);
187  W_op->apply(Thyra::CONJTRANS, *adjoint_mv, tmp.ptr(), Scalar(1.0),
188  Scalar(0.0));
189  dxdp_init_->apply(Thyra::CONJTRANS, *tmp, dgdp_.ptr(), Scalar(1.0),
190  Scalar(1.0));
191  }
192  }
193 
194  // Add in model parameter term = \int_0^T( (df/dp^T(t)*y(t) )dt which
195  // is computed during the adjoint integration as an auxiliary integral
196  // (2nd block of the solution vector)
197  if (f_depends_on_p_) {
198  RCP<const SolutionState<Scalar> > adjoint_state =
199  adjoint_solution_history->getCurrentState();
200  RCP<const VectorBase<Scalar> > z =
201  rcp_dynamic_cast<const DPV>(adjoint_state->getX())->getVectorBlock(1);
202  RCP<const MultiVectorBase<Scalar> > z_mv =
203  rcp_dynamic_cast<const DMVPV>(z)->getMultiVector();
204  Thyra::V_VmV(dgdp_.ptr(), *dgdp_, *z_mv);
205  }
206 
207  buildSolutionHistory(state_solution_history, adjoint_solution_history);
208 
209  return state_status && sens_status;
210 }
211 
212 template<class Scalar>
213 Scalar
215 getTime() const
216 {
217  return solutionHistory_->getCurrentTime();
218 }
219 
220 template<class Scalar>
221 Scalar
223 getIndex() const
224 {
225  return solutionHistory_->getCurrentIndex();
226 }
227 
228 template<class Scalar>
229 Status
231 getStatus() const
232 {
233  Status state_status = state_integrator_->getStatus();
234  Status sens_status = adjoint_integrator_->getStatus();
235  if (state_status == FAILED || sens_status == FAILED)
236  return FAILED;
237  if (state_status == WORKING || sens_status == WORKING)
238  return WORKING;
239  return PASSED;
240 }
241 
242 template<class Scalar>
243 Teuchos::RCP<Stepper<Scalar> >
245 getStepper() const
246 {
247  return state_integrator_->getStepper();
248 }
249 
250 template<class Scalar>
251 Teuchos::RCP<Teuchos::ParameterList>
254 {
255  return state_integrator_->getTempusParameterList();
256 }
257 
258 template<class Scalar>
259 void
261 setTempusParameterList(Teuchos::RCP<Teuchos::ParameterList> pl)
262 {
263  state_integrator_->setTempusParameterList(pl);
264  adjoint_integrator_->setTempusParameterList(pl);
265 }
266 
267 template<class Scalar>
268 Teuchos::RCP<const SolutionHistory<Scalar> >
271 {
272  return solutionHistory_;
273 }
274 
275 template<class Scalar>
276 Teuchos::RCP<const TimeStepControl<Scalar> >
279 {
280  return state_integrator_->getTimeStepControl();
281 }
282 
283 template<class Scalar>
285 setInitialState(Scalar t0,
286  Teuchos::RCP<const Thyra::VectorBase<Scalar> > x0,
287  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xdot0,
288  Teuchos::RCP<const Thyra::VectorBase<Scalar> > xdotdot0,
289  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > DxDp0,
290  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > DxdotDp0,
291  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> > DxdotdotDp0)
292 {
293  state_integrator_->setInitialState(t0, x0, xdot0, xdotdot0);
294  dxdp_init_ = DxDp0;
295 }
296 
297 template<class Scalar>
298 Teuchos::RCP<const Thyra::VectorBase<Scalar> >
300 getX() const
301 {
302  return state_integrator_->getX();
303 }
304 
305 template<class Scalar>
306 Teuchos::RCP<const Thyra::VectorBase<Scalar> >
308 getXdot() const
309 {
310  return state_integrator_->getXdot();
311 }
312 
313 template<class Scalar>
314 Teuchos::RCP<const Thyra::VectorBase<Scalar> >
316 getXdotdot() const
317 {
318  return state_integrator_->getXdotdot();
319 }
320 
321 template<class Scalar>
322 Teuchos::RCP<const Thyra::MultiVectorBase<Scalar> >
324 getDgDp() const
325 {
326  return dgdp_;
327 }
328 
329 template<class Scalar>
330 std::string
332 description() const
333 {
334  std::string name = "Tempus::IntegratorAdjointSensitivity";
335  return(name);
336 }
337 
338 template<class Scalar>
339 void
342  Teuchos::FancyOStream &out,
343  const Teuchos::EVerbosityLevel verbLevel) const
344 {
345  out << description() << "::describe" << std::endl;
346  state_integrator_->describe(out, verbLevel);
347  adjoint_integrator_->describe(out, verbLevel);
348 }
349 
350 template<class Scalar>
351 void
353 setParameterList(const Teuchos::RCP<Teuchos::ParameterList> & inputPL)
354 {
355  if (state_integrator_ != Teuchos::null)
356  state_integrator_->setParameterList(inputPL);
357  if (adjoint_integrator_ != Teuchos::null)
358  adjoint_integrator_->setParameterList(inputPL);
359  Teuchos::ParameterList& spl = inputPL->sublist("Sensitivities");
360  p_index_ = spl.get<int>("Sensitivity Parameter Index", 0);
361  g_index_ = spl.get<int>("Response Function Index", 0);
362  g_depends_on_p_ = spl.get<bool>("Response Depends on Parameters", true);
363  f_depends_on_p_ = spl.get<bool>("Residual Depends on Parameters", true);
364  ic_depends_on_p_ = spl.get<bool>("IC Depends on Parameters", true);
365  mass_matrix_is_identity_ = spl.get<bool>("Mass Matrix Is Identity", false);
366 }
367 
368 template<class Scalar>
369 Teuchos::RCP<Teuchos::ParameterList>
372 {
373  state_integrator_->unsetParameterList();
374  return adjoint_integrator_->unsetParameterList();
375 }
376 
377 template<class Scalar>
378 Teuchos::RCP<const Teuchos::ParameterList>
381 {
382  using Teuchos::RCP;
383  using Teuchos::ParameterList;
384  using Teuchos::parameterList;
385 
386  RCP<ParameterList> pl = parameterList();
387  RCP<const ParameterList> integrator_pl =
388  state_integrator_->getValidParameters();
389  RCP<const ParameterList> sensitivity_pl =
391  pl->setParameters(*integrator_pl);
392  ParameterList& spl = pl->sublist("Sensitivities");
393  spl.setParameters(*sensitivity_pl);
394  spl.set<bool>("Response Depends on Parameters", true);
395  spl.set<bool>("Residual Depends on Parameters", true);
396  spl.set<bool>("IC Depends on Parameters", true);
397 
398  return pl;
399 }
400 
401 template<class Scalar>
402 Teuchos::RCP<Teuchos::ParameterList>
405 {
406  return state_integrator_->getNonconstParameterList();
407 }
408 
409 template <class Scalar>
410 Teuchos::RCP<Tempus::AdjointAuxSensitivityModelEvaluator<Scalar> >
413  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model,
414  const Teuchos::RCP<Teuchos::ParameterList>& inputPL)
415 {
416  using Teuchos::RCP;
417  using Teuchos::rcp;
418  using Teuchos::ParameterList;
419 
420  RCP<ParameterList> spl = Teuchos::parameterList();
421  if (inputPL != Teuchos::null) {
422  *spl = inputPL->sublist("Sensitivities");
423  }
424  spl->remove("Response Depends on Parameters");
425  spl->remove("Residual Depends on Parameters");
426  spl->remove("IC Depends on Parameters");
427  const Scalar tfinal = state_integrator_->getTimeStepControl()->getFinalTime();
429  model, tfinal, spl));
430 }
431 
432 template<class Scalar>
433 void
436  const Teuchos::RCP<const SolutionHistory<Scalar> >& state_solution_history,
437  const Teuchos::RCP<const SolutionHistory<Scalar> >& adjoint_solution_history)
438 {
439  using Teuchos::RCP;
440  using Teuchos::rcp;
441  using Teuchos::rcp_dynamic_cast;
442  using Teuchos::ParameterList;
443  using Thyra::VectorBase;
444  using Thyra::MultiVectorBase;
445  using Thyra::VectorSpaceBase;
446  using Thyra::createMembers;
447  using Thyra::multiVectorProductVector;
448  using Thyra::assign;
449  typedef Thyra::DefaultProductVectorSpace<Scalar> DPVS;
450  typedef Thyra::DefaultProductVector<Scalar> DPV;
451 
452  // Create combined solution histories combining the forward and adjoint
453  // solutions. We do not include the auxiliary part from the adjoint solution.
454  RCP<ParameterList> shPL =
455  Teuchos::sublist(state_integrator_->getIntegratorParameterList(),
456  "Solution History", true);
457  solutionHistory_ = rcp(new SolutionHistory<Scalar>(shPL));
458 
459  RCP<const VectorSpaceBase<Scalar> > x_space = model_->get_x_space();
460  RCP<const VectorSpaceBase<Scalar> > adjoint_space =
461  rcp_dynamic_cast<const DPVS>(adjoint_model_->get_x_space())->getBlock(0);
462  Teuchos::Array< RCP<const VectorSpaceBase<Scalar> > > spaces(2);
463  spaces[0] = x_space;
464  spaces[1] = adjoint_space;
465  RCP<const DPVS > prod_space = Thyra::productVectorSpace(spaces());
466 
467  int num_states = state_solution_history->getNumStates();
468  const Scalar t_final = state_integrator_->getTime();
469  for (int i=0; i<num_states; ++i) {
470  RCP<const SolutionState<Scalar> > forward_state =
471  (*state_solution_history)[i];
472  RCP<const SolutionState<Scalar> > adjoint_state =
473  adjoint_solution_history->findState(t_final-forward_state->getTime());
474 
475  // X
476  RCP<DPV> x = Thyra::defaultProductVector(prod_space);
477  RCP<const VectorBase<Scalar> > adjoint_x =
478  rcp_dynamic_cast<const DPV>(adjoint_state->getX())->getVectorBlock(0);
479  assign(x->getNonconstVectorBlock(0).ptr(), *(forward_state->getX()));
480  assign(x->getNonconstVectorBlock(1).ptr(), *(adjoint_x));
481  RCP<VectorBase<Scalar> > x_b = x;
482 
483  // X-Dot
484  RCP<DPV> x_dot = Thyra::defaultProductVector(prod_space);
485  RCP<const VectorBase<Scalar> > adjoint_x_dot =
486  rcp_dynamic_cast<const DPV>(adjoint_state->getXDot())->getVectorBlock(0);
487  assign(x_dot->getNonconstVectorBlock(0).ptr(), *(forward_state->getXDot()));
488  assign(x_dot->getNonconstVectorBlock(1).ptr(), *(adjoint_x_dot));
489  RCP<VectorBase<Scalar> > x_dot_b = x_dot;
490 
491  // X-Dot-Dot
492  RCP<DPV> x_dot_dot;
493  if (forward_state->getXDotDot() != Teuchos::null) {
494  x_dot_dot = Thyra::defaultProductVector(prod_space);
495  RCP<const VectorBase<Scalar> > adjoint_x_dot_dot =
496  rcp_dynamic_cast<const DPV>(
497  adjoint_state->getXDotDot())->getVectorBlock(0);
498  assign(x_dot_dot->getNonconstVectorBlock(0).ptr(),
499  *(forward_state->getXDotDot()));
500  assign(x_dot_dot->getNonconstVectorBlock(1).ptr(),
501  *(adjoint_x_dot_dot));
502  }
503  RCP<VectorBase<Scalar> > x_dot_dot_b = x_dot_dot;
504 
505  RCP<SolutionState<Scalar> > prod_state =
506  rcp(new SolutionState<Scalar>(forward_state->getMetaData()->clone(),
507  x_b, x_dot_b, x_dot_dot_b,
508  forward_state->getStepperState()->clone(),
509  Teuchos::null));
510  solutionHistory_->addState(prod_state);
511  }
512 }
513 
514 /// Non-member constructor
515 template<class Scalar>
516 Teuchos::RCP<Tempus::IntegratorAdjointSensitivity<Scalar> >
518  Teuchos::RCP<Teuchos::ParameterList> pList,
519  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar> >& model)
520 {
521  Teuchos::RCP<Tempus::IntegratorAdjointSensitivity<Scalar> > integrator =
522  Teuchos::rcp(new Tempus::IntegratorAdjointSensitivity<Scalar>(pList, model));
523  return(integrator);
524 }
525 
526 /// Non-member constructor
527 template<class Scalar>
528 Teuchos::RCP<Tempus::IntegratorAdjointSensitivity<Scalar> >
530 {
531  Teuchos::RCP<Tempus::IntegratorAdjointSensitivity<Scalar> > integrator =
533  return(integrator);
534 }
535 
536 } // namespace Tempus
537 #endif // Tempus_IntegratorAdjointSensitivity_impl_hpp
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
Teuchos::RCP< Tempus::IntegratorAdjointSensitivity< Scalar > > integratorAdjointSensitivity(Teuchos::RCP< Teuchos::ParameterList > pList, const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model)
Non-member constructor.
Teuchos::RCP< Tempus::AdjointAuxSensitivityModelEvaluator< Scalar > > createAdjointModel(const Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< Teuchos::ParameterList > &inputPL)
Teuchos::RCP< Teuchos::ParameterList > getNonconstParameterList() override
Teuchos::RCP< const Teuchos::ParameterList > getValidParameters() const override
static Teuchos::RCP< const Teuchos::ParameterList > getValidParameters()
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDgDp() const
Return adjoint sensitivity stored in gradient format.
ModelEvaluator for forming adjoint sensitivity equations.
Teuchos::RCP< Teuchos::ParameterList > unsetParameterList() override
virtual Scalar getIndex() const override
Get current index.
void buildSolutionHistory(const Teuchos::RCP< const SolutionHistory< Scalar > > &state_solution_history, const Teuchos::RCP< const SolutionHistory< Scalar > > &adjoint_solution_history)
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getX() const
Get current the solution, x.
virtual Teuchos::RCP< const SolutionHistory< Scalar > > getSolutionHistory() const override
Get the SolutionHistory.
Status
Status for the Integrator, the Stepper and the SolutionState.
virtual Teuchos::RCP< Teuchos::ParameterList > getTempusParameterList() override
Return a copy of the Tempus ParameterList.
virtual Teuchos::RCP< const TimeStepControl< Scalar > > getTimeStepControl() const override
Get the TimeStepControl.
virtual void setTempusParameterList(Teuchos::RCP< Teuchos::ParameterList > pl) override
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
virtual bool advanceTime()
Advance the solution to timeMax, and return true if successful.
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXdotdot() const
Get current the second time derivative of the solution, xdotdot.
void setParameterList(const Teuchos::RCP< Teuchos::ParameterList > &pl) override
Time integrator suitable for adjoint sensitivity analysis.
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXdot() const
Get current the time derivative of the solution, xdot.
virtual Scalar getTime() const override
Get current time.
virtual Status getStatus() const override
Get Status.
Solution state for integrators and steppers. SolutionState contains the metadata for solutions and th...
virtual Teuchos::RCP< Stepper< Scalar > > getStepper() const override
Get the Stepper.
virtual void setInitialState(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 > > DxDp0=Teuchos::null, Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > DxdotDp0=Teuchos::null, Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > DxdotdotDp0=Teuchos::null)
Set the initial state from Thyra::VectorBase(s)