Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_IntegratorPseudoTransientForwardSensitivity_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_IntegratorPseudoTransientForwardSensitivity_impl_hpp
11 #define Tempus_IntegratorPseudoTransientForwardSensitivity_impl_hpp
12 
13 #include "Thyra_DefaultMultiVectorProductVector.hpp"
14 #include "Thyra_VectorStdOps.hpp"
15 #include "Thyra_MultiVectorStdOps.hpp"
16 
18 
19 namespace Tempus {
20 
21 template <class Scalar>
26  const Teuchos::RCP<IntegratorBasic<Scalar>> &fwd_integrator,
27  const Teuchos::RCP<IntegratorBasic<Scalar>> &sens_integrator,
28  const bool reuse_solver, const bool force_W_update)
29  : model_(model),
30  sens_model_(sens_model),
31  state_integrator_(fwd_integrator),
32  sens_integrator_(sens_integrator),
33  reuse_solver_(reuse_solver),
34  force_W_update_(force_W_update),
35  stepMode_(SensitivityStepMode::Forward)
36 {
37 }
38 
39 template <class Scalar>
42  : reuse_solver_(false),
43  force_W_update_(false),
44  stepMode_(SensitivityStepMode::Forward)
45 {
46  state_integrator_ = createIntegratorBasic<Scalar>();
47  sens_integrator_ = createIntegratorBasic<Scalar>();
48 }
49 
50 template <class Scalar>
52 {
53  using Teuchos::RCP;
54  using Thyra::VectorBase;
55 
56  // Run state integrator and get solution
57  stepMode_ = SensitivityStepMode::Forward;
58  bool state_status = state_integrator_->advanceTime();
59 
60  // Set solution in sensitivity ME
61  sens_model_->setForwardSolutionState(state_integrator_->getCurrentState());
62 
63  // Reuse state solver if requested
64  if (reuse_solver_ &&
65  state_integrator_->getStepper()->getSolver() != Teuchos::null) {
66  sens_model_->setSolver(state_integrator_->getStepper()->getSolver(),
67  force_W_update_);
68  }
69 
70  // Run sensitivity integrator
72  bool sens_status = sens_integrator_->advanceTime();
73 
74  buildSolutionHistory();
75 
76  return state_status && sens_status;
77 }
78 
79 template <class Scalar>
81  const Scalar timeFinal)
82 {
83  TEMPUS_FUNC_TIME_MONITOR_DIFF(
84  "Tempus::IntegratorPseudoTransientForwardSensitivity::advanceTime()",
85  TEMPUS_PTFS_AT);
86 
87  using Teuchos::RCP;
88  using Thyra::VectorBase;
89 
90  // Run state integrator and get solution
91  bool state_status = true;
92  {
93  TEMPUS_FUNC_TIME_MONITOR_DIFF(
94  "Tempus::IntegratorPseudoTransientForwardSensitivity::advanceTime::"
95  "state",
96  TEMPUS_PTFS_AT_FWD);
97  state_status = state_integrator_->advanceTime(timeFinal);
98  }
99 
100  // Set solution in sensitivity ME
101  sens_model_->setForwardSolutionState(state_integrator_->getCurrentState());
102 
103  // Reuse state solver if requested
104  if (reuse_solver_ &&
105  state_integrator_->getStepper()->getSolver() != Teuchos::null) {
106  sens_model_->setSolver(state_integrator_->getStepper()->getSolver(),
107  force_W_update_);
108  }
109 
110  // Run sensitivity integrator
111  bool sens_status = true;
112  {
113  TEMPUS_FUNC_TIME_MONITOR_DIFF(
114  "Tempus::IntegratorPseudoTransientForwardSensitivity::advanceTime::"
115  "sensitivity",
116  TEMPUS_PTFS_AT_SEN);
117  sens_status = sens_integrator_->advanceTime(timeFinal);
118  }
119 
120  buildSolutionHistory();
121 
122  return state_status && sens_status;
123 }
124 
125 template <class Scalar>
127 {
128  return solutionHistory_->getCurrentTime();
129 }
130 
131 template <class Scalar>
133 {
134  return solutionHistory_->getCurrentIndex();
135 }
136 
137 template <class Scalar>
139 {
140  Status state_status = state_integrator_->getStatus();
141  Status sens_status = sens_integrator_->getStatus();
142  if (state_status == FAILED || sens_status == FAILED) return FAILED;
143  if (state_status == WORKING || sens_status == WORKING) return WORKING;
144  return PASSED;
145 }
146 
147 template <class Scalar>
149  const Status st)
150 {
151  state_integrator_->setStatus(st);
152  sens_integrator_->setStatus(st);
153 }
154 
155 template <class Scalar>
158 {
159  return state_integrator_->getStepper();
160 }
161 
162 template <class Scalar>
165 {
166  return state_integrator_->getStepper();
167 }
168 
169 template <class Scalar>
172 {
173  return sens_integrator_->getStepper();
174 }
175 
176 template <class Scalar>
179 {
180  return solutionHistory_;
181 }
182 
183 template <class Scalar>
186  const
187 {
188  return state_integrator_->getSolutionHistory();
189 }
190 
191 template <class Scalar>
194  const
195 {
196  return sens_integrator_->getSolutionHistory();
197 }
198 
199 template <class Scalar>
202  Scalar>::getNonConstSolutionHistory()
203 {
204  return solutionHistory_;
205 }
206 
207 template <class Scalar>
210 {
211  return state_integrator_->getTimeStepControl();
212 }
213 
214 template <class Scalar>
217  Scalar>::getNonConstTimeStepControl()
218 {
219  return state_integrator_->getNonConstTimeStepControl();
220 }
221 
222 template <class Scalar>
225  Scalar>::getStateNonConstTimeStepControl()
226 {
227  return state_integrator_->getNonConstTimeStepControl();
228 }
229 
230 template <class Scalar>
233  Scalar>::getSensNonConstTimeStepControl()
234 {
235  return sens_integrator_->getNonConstTimeStepControl();
236 }
237 
238 template <class Scalar>
241 {
242  return state_integrator_->getObserver();
243 }
244 
245 template <class Scalar>
248 {
249  state_integrator_->setObserver(obs);
250  sens_integrator_->setObserver(obs);
251 }
252 
253 template <class Scalar>
256  Scalar t0, Teuchos::RCP<const Thyra::VectorBase<Scalar>> x0,
258  Teuchos::RCP<const Thyra::VectorBase<Scalar>> xdotdot0,
261  Teuchos::RCP<const Thyra::MultiVectorBase<Scalar>> DxdotdotDp0)
262 {
263  using Teuchos::RCP;
264  using Teuchos::rcp_dynamic_cast;
265  using Thyra::assign;
266  using Thyra::createMember;
269 
270  //
271  // Create and initialize product X, Xdot, Xdotdot
272 
273  RCP<const VectorSpaceBase<Scalar>> space = sens_model_->get_x_space();
274  RCP<DMVPV> X = rcp_dynamic_cast<DMVPV>(createMember(space));
275  RCP<DMVPV> Xdot = rcp_dynamic_cast<DMVPV>(createMember(space));
276  RCP<DMVPV> Xdotdot = rcp_dynamic_cast<DMVPV>(createMember(space));
277  const Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
278 
279  // x
280  if (DxDp0 == Teuchos::null)
281  assign(X->getNonconstMultiVector().ptr(), zero);
282  else
283  assign(X->getNonconstMultiVector().ptr(), *DxDp0);
284 
285  // xdot
286  if (DxdotDp0 == Teuchos::null)
287  assign(Xdot->getNonconstMultiVector().ptr(), zero);
288  else
289  assign(Xdot->getNonconstMultiVector().ptr(), *DxdotDp0);
290 
291  // xdotdot
292  if (DxdotDp0 == Teuchos::null)
293  assign(Xdotdot->getNonconstMultiVector().ptr(), zero);
294  else
295  assign(Xdotdot->getNonconstMultiVector().ptr(), *DxdotdotDp0);
296 
297  state_integrator_->initializeSolutionHistory(t0, x0, xdot0, xdotdot0);
298  sens_integrator_->initializeSolutionHistory(t0, X, Xdot, Xdotdot);
299 }
300 
301 template <class Scalar>
304 {
305  return state_integrator_->getX();
306 }
307 
308 template <class Scalar>
311 {
312  using Teuchos::RCP;
313  using Teuchos::rcp_dynamic_cast;
315 
316  RCP<const DMVPV> X = rcp_dynamic_cast<const DMVPV>(sens_integrator_->getX());
317  return X->getMultiVector();
318 }
319 
320 template <class Scalar>
323 {
324  return state_integrator_->getXDot();
325 }
326 
327 template <class Scalar>
330 {
331  using Teuchos::RCP;
332  using Teuchos::rcp_dynamic_cast;
334 
335  RCP<const DMVPV> Xdot =
336  rcp_dynamic_cast<const DMVPV>(sens_integrator_->getXDot());
337  return Xdot->getMultiVector();
338 }
339 
340 template <class Scalar>
343 {
344  return state_integrator_->getXDotDot();
345 }
346 
347 template <class Scalar>
350 {
351  using Teuchos::RCP;
352  using Teuchos::rcp_dynamic_cast;
354 
355  RCP<const DMVPV> Xdotdot =
356  rcp_dynamic_cast<const DMVPV>(sens_integrator_->getXDotDot());
357  return Xdotdot->getMultiVector();
358 }
359 
360 template <class Scalar>
363 {
364  typedef Thyra::ModelEvaluatorBase MEB;
365 
366  // Compute g which is computed by response 1 of the
367  // sensitivity model evaluator
368  MEB::InArgs<Scalar> inargs = sens_model_->getNominalValues();
369  MEB::OutArgs<Scalar> outargs = sens_model_->createOutArgs();
370  inargs.set_t(sens_integrator_->getTime());
371  inargs.set_x(sens_integrator_->getX());
372  if (inargs.supports(MEB::IN_ARG_x_dot))
373  inargs.set_x_dot(sens_integrator_->getXDot());
374  if (inargs.supports(MEB::IN_ARG_x_dot_dot))
375  inargs.set_x_dot_dot(sens_integrator_->getXDotDot());
376 
378  Thyra::createMember(sens_model_->get_g_space(1));
379  outargs.set_g(1, g);
380 
381  sens_model_->evalModel(inargs, outargs);
382  return g;
383 }
384 
385 template <class Scalar>
388 {
389  typedef Thyra::ModelEvaluatorBase MEB;
391 
392  // Compute final dg/dp which is computed by response 0 of the
393  // sensitivity model evaluator
394  MEB::InArgs<Scalar> inargs = sens_model_->getNominalValues();
395  MEB::OutArgs<Scalar> outargs = sens_model_->createOutArgs();
396  inargs.set_t(sens_integrator_->getTime());
397  inargs.set_x(sens_integrator_->getX());
398  if (inargs.supports(MEB::IN_ARG_x_dot))
399  inargs.set_x_dot(sens_integrator_->getXDot());
400  if (inargs.supports(MEB::IN_ARG_x_dot_dot))
401  inargs.set_x_dot_dot(sens_integrator_->getXDotDot());
402 
404  Thyra::createMember(sens_model_->get_g_space(0));
405  Teuchos::RCP<DMVPV> dgdp = Teuchos::rcp_dynamic_cast<DMVPV>(G);
406  outargs.set_g(0, G);
407 
408  sens_model_->evalModel(inargs, outargs);
409  return dgdp->getMultiVector();
410 }
411 
412 template <class Scalar>
414  const
415 {
416  std::string name = "Tempus::IntegratorPseudoTransientForwardSensitivity";
417  return (name);
418 }
419 
420 template <class Scalar>
422  Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const
423 {
424  auto l_out = Teuchos::fancyOStream(out.getOStream());
425  Teuchos::OSTab ostab(*l_out, 2, this->description());
426  l_out->setOutputToRootOnly(0);
427 
428  *l_out << description() << "::describe" << std::endl;
429  state_integrator_->describe(*l_out, verbLevel);
430  sens_integrator_->describe(*l_out, verbLevel);
431 }
432 
433 template <class Scalar>
436 {
437  return stepMode_;
438 }
439 
440 template <class Scalar>
442 {
444  using Teuchos::RCP;
445  using Teuchos::rcp;
446  using Teuchos::rcp_dynamic_cast;
447  using Thyra::assign;
448  using Thyra::createMembers;
450  using Thyra::multiVectorProductVector;
451  using Thyra::VectorBase;
454 
455  // TODO: get the solution history PL or create it
456 
457  // Create combined solution histories, first for the states with zero
458  // sensitivities and then for the sensitivities with frozen states
459  RCP<ParameterList> shPL;
460  // Teuchos::sublist(state_integrator_->getIntegratorParameterList(), "Solution
461  // History", true);
462  solutionHistory_ = createSolutionHistoryPL<Scalar>(shPL);
463 
464  const int num_param = rcp_dynamic_cast<const DMVPV>(sens_integrator_->getX())
465  ->getMultiVector()
466  ->domain()
467  ->dim();
468  RCP<const VectorSpaceBase<Scalar>> x_space = model_->get_x_space();
469  RCP<const Thyra::DefaultMultiVectorProductVectorSpace<Scalar>> prod_space =
470  Thyra::multiVectorProductVectorSpace(x_space, num_param + 1);
471  const Teuchos::Range1D rng(1, num_param);
472  const Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
473 
474  RCP<const SolutionHistory<Scalar>> state_solution_history =
475  state_integrator_->getSolutionHistory();
476  int num_states = state_solution_history->getNumStates();
477  for (int i = 0; i < num_states; ++i) {
478  RCP<const SolutionState<Scalar>> state = (*state_solution_history)[i];
479 
480  // X
481  RCP<MultiVectorBase<Scalar>> x_mv = createMembers(x_space, num_param + 1);
482  assign(x_mv->col(0).ptr(), *(state->getX()));
483  assign(x_mv->subView(rng).ptr(), zero);
484  RCP<VectorBase<Scalar>> x = multiVectorProductVector(prod_space, x_mv);
485 
486  // X-Dot
487  RCP<VectorBase<Scalar>> x_dot;
488  if (state->getXDot() != Teuchos::null) {
489  RCP<MultiVectorBase<Scalar>> x_dot_mv =
490  createMembers(x_space, num_param + 1);
491  assign(x_dot_mv->col(0).ptr(), *(state->getXDot()));
492  assign(x_dot_mv->subView(rng).ptr(), zero);
493  x_dot = multiVectorProductVector(prod_space, x_dot_mv);
494  }
495 
496  // X-Dot-Dot
497  RCP<VectorBase<Scalar>> x_dot_dot;
498  if (state->getXDotDot() != Teuchos::null) {
499  RCP<MultiVectorBase<Scalar>> x_dot_dot_mv =
500  createMembers(x_space, num_param + 1);
501  assign(x_dot_dot_mv->col(0).ptr(), *(state->getXDotDot()));
502  assign(x_dot_dot_mv->subView(rng).ptr(), zero);
503  x_dot_dot = multiVectorProductVector(prod_space, x_dot_dot_mv);
504  }
505 
506  RCP<SolutionState<Scalar>> prod_state = state->clone();
507  prod_state->setX(x);
508  prod_state->setXDot(x_dot);
509  prod_state->setXDotDot(x_dot_dot);
510  solutionHistory_->addState(prod_state);
511  }
512 
513  RCP<const VectorBase<Scalar>> frozen_x =
514  state_solution_history->getCurrentState()->getX();
515  RCP<const VectorBase<Scalar>> frozen_x_dot =
516  state_solution_history->getCurrentState()->getXDot();
517  RCP<const VectorBase<Scalar>> frozen_x_dot_dot =
518  state_solution_history->getCurrentState()->getXDotDot();
519  RCP<const SolutionHistory<Scalar>> sens_solution_history =
520  sens_integrator_->getSolutionHistory();
521  num_states = sens_solution_history->getNumStates();
522  for (int i = 0; i < num_states; ++i) {
523  RCP<const SolutionState<Scalar>> state = (*sens_solution_history)[i];
524 
525  // X
526  RCP<MultiVectorBase<Scalar>> x_mv = createMembers(x_space, num_param + 1);
527  RCP<const MultiVectorBase<Scalar>> dxdp =
528  rcp_dynamic_cast<const DMVPV>(state->getX())->getMultiVector();
529  assign(x_mv->col(0).ptr(), *(frozen_x));
530  assign(x_mv->subView(rng).ptr(), *dxdp);
531  RCP<VectorBase<Scalar>> x = multiVectorProductVector(prod_space, x_mv);
532 
533  // X-Dot
534  RCP<VectorBase<Scalar>> x_dot;
535  if (state->getXDot() != Teuchos::null) {
536  RCP<MultiVectorBase<Scalar>> x_dot_mv =
537  createMembers(x_space, num_param + 1);
538  RCP<const MultiVectorBase<Scalar>> dxdotdp =
539  rcp_dynamic_cast<const DMVPV>(state->getXDot())->getMultiVector();
540  assign(x_dot_mv->col(0).ptr(), *(frozen_x_dot));
541  assign(x_dot_mv->subView(rng).ptr(), *dxdotdp);
542  x_dot = multiVectorProductVector(prod_space, x_dot_mv);
543  }
544 
545  // X-Dot-Dot
546  RCP<VectorBase<Scalar>> x_dot_dot;
547  if (state->getXDotDot() != Teuchos::null) {
548  RCP<MultiVectorBase<Scalar>> x_dot_dot_mv =
549  createMembers(x_space, num_param + 1);
550  RCP<const MultiVectorBase<Scalar>> dxdotdotdp =
551  rcp_dynamic_cast<const DMVPV>(state->getXDotDot())->getMultiVector();
552  assign(x_dot_dot_mv->col(0).ptr(), *(frozen_x_dot_dot));
553  assign(x_dot_dot_mv->subView(rng).ptr(), *dxdotdotdp);
554  x_dot_dot = multiVectorProductVector(prod_space, x_dot_dot_mv);
555  }
556 
557  RCP<SolutionState<Scalar>> prod_state = state->clone();
558  prod_state->setX(x);
559  prod_state->setXDot(x_dot);
560  prod_state->setXDotDot(x_dot_dot);
561  solutionHistory_->addState(prod_state, false);
562  }
563 }
564 
566 template <class Scalar>
571  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar>> &sens_residual_model,
572  const Teuchos::RCP<Thyra::ModelEvaluator<Scalar>> &sens_solve_model)
573 {
574  auto fwd_integrator = createIntegratorBasic<Scalar>(pList, model);
576  Teuchos::RCP<IntegratorBasic<Scalar>> sens_integrator;
577 
578  {
582  fwd_integrator->getValidParameters();
585  pl->setParameters(*integrator_pl);
586  pl->sublist("Sensitivities").setParameters(*sensitivity_pl);
587  pl->sublist("Sensitivities").set("Reuse State Linear Solver", false);
588  pl->sublist("Sensitivities").set("Force W Update", false);
589  pl->sublist("Sensitivities").set("Cache Matrices", false);
590  pList->setParametersNotAlreadySet(*pl);
591  }
592 
593  bool reuse_solver =
594  pList->sublist("Sensitivities").get("Reuse State Linear Solver", false);
595  bool force_W_update =
596  pList->sublist("Sensitivities").get("Force W Update", false);
597  bool cache_matrices =
598  pList->sublist("Sensitivities").get("Cache Matrices", false);
599 
600  {
601  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
602  if (pList != Teuchos::null) {
603  *pl = pList->sublist("Sensitivities");
604  }
605  pl->remove("Reuse State Linear Solver");
606  pl->remove("Force W Update");
607  pl->remove("Cache Matrices");
608  sens_model = wrapStaggeredFSAModelEvaluator(
609  model, sens_residual_model, sens_solve_model, cache_matrices, pl);
610  sens_integrator = createIntegratorBasic<Scalar>(pList, sens_model);
611  }
612 
614  integrator = Teuchos::rcp(
616  model, sens_model, fwd_integrator, sens_integrator, reuse_solver,
617  force_W_update));
618 
619  return (integrator);
620 }
621 
623 template <class Scalar>
626 {
628  integrator = Teuchos::rcp(
630  return (integrator);
631 }
632 
633 } // namespace Tempus
634 #endif // Tempus_IntegratorPseudoTransientForwardSensitivity_impl_hpp
virtual Teuchos::RCP< Stepper< Scalar > > getStepper() const override
Get the Stepper.
Teuchos::RCP< Tempus::IntegratorPseudoTransientForwardSensitivity< Scalar > > createIntegratorPseudoTransientForwardSensitivity(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.
T & get(const std::string &name, T def_value)
Teuchos::RCP< SensitivityModelEvaluatorBase< Scalar > > wrapStaggeredFSAModelEvaluator(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 bool is_pseudotransient, const Teuchos::RCP< const Teuchos::ParameterList > &pList=Teuchos::null)
virtual Teuchos::RCP< const SolutionHistory< Scalar > > getSolutionHistory() const override
Get the SolutionHistory.
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDxDp() const
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 >> 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)
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
virtual Teuchos::RCP< TimeStepControl< Scalar > > getNonConstTimeStepControl() override
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getG() const
Return response function g.
virtual bool advanceTime()
Advance the solution to timeMax, and return true if successful.
bool remove(std::string const &name, bool throwIfNotExists=true)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual void setObserver(Teuchos::RCP< IntegratorObserver< Scalar >> obs=Teuchos::null)
Set the Observer.
SensitivityStepMode getStepMode() const
What mode the current time integration step is in.
Status
Status for the Integrator, the Stepper and the SolutionState.
IntegratorObserver class for time integrators.
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDgDp() const
Return forward sensitivity stored in Jacobian format.
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXDotDot() const
Get current the second time derivative of the solution, xdotdot.
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
RCP< const MultiVectorBase< Scalar > > getMultiVector() const
virtual Teuchos::RCP< IntegratorObserver< Scalar > > getObserver()
Get the Observer.
ParameterList & setParameters(const ParameterList &source)
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDXDotDp() const
RCP< std::basic_ostream< char_type, traits_type > > getOStream()
Time integrator suitable for pseudotransient forward sensitivity analysis.
ParameterList & setParametersNotAlreadySet(const ParameterList &source)
virtual RCP< const VectorSpaceBase< Scalar > > domain() const
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
virtual Teuchos::RCP< const TimeStepControl< Scalar > > getTimeStepControl() const override
Get the TimeStepControl.
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDXDotDotDp() const
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXDot() const
Get current the time derivative of the solution, xdot.
A ModelEvaluator decorator for sensitivity analysis.
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getX() const
Get current the solution, x.