10#ifndef Tempus_IntegratorForwardSensitivity_impl_hpp
11#define Tempus_IntegratorForwardSensitivity_impl_hpp
14#include "Thyra_DefaultMultiVectorProductVector.hpp"
15#include "Thyra_VectorStdOps.hpp"
16#include "Thyra_MultiVectorStdOps.hpp"
18#include "Tempus_CombinedForwardSensitivityModelEvaluator.hpp"
23template <
class Scalar>
30 bool use_combined_method)
32 integrator_(integrator),
33 sens_model_(sens_model),
34 sens_stepper_(sens_stepper),
35 use_combined_method_(use_combined_method)
40template <
class Scalar>
43 integrator_ = createIntegratorBasic<Scalar>();
44 integrator_->initialize();
47template <
class Scalar>
51 if (use_combined_method_)
52 integrator_->setModel(sens_model_);
54 integrator_->setStepper(sens_stepper_);
57template <
class Scalar>
67 using Teuchos::rcp_dynamic_cast;
69 using Thyra::createMember;
70 using Thyra::VectorSpaceBase;
71 typedef Thyra::DefaultMultiVectorProductVector<Scalar> DMVPV;
76 RCP<const VectorSpaceBase<Scalar>> space;
77 if (use_combined_method_)
78 space = sens_model_->get_x_space();
80 space = sens_stepper_->get_x_space();
81 RCP<DMVPV> X = rcp_dynamic_cast<DMVPV>(createMember(space));
82 RCP<DMVPV> Xdot = rcp_dynamic_cast<DMVPV>(createMember(space));
83 RCP<DMVPV> Xdotdot = rcp_dynamic_cast<DMVPV>(createMember(space));
85 const int num_param = X->getNonconstMultiVector()->domain()->dim() - 1;
86 const Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
87 const Teuchos::Range1D rng(1, num_param);
90 assign(X->getNonconstMultiVector()->col(0).ptr(), *x0);
91 if (DxDp0 == Teuchos::null)
92 assign(X->getNonconstMultiVector()->subView(rng).ptr(), zero);
94 assign(X->getNonconstMultiVector()->subView(rng).ptr(), *DxDp0);
97 if (xdot0 == Teuchos::null)
98 assign(Xdot->getNonconstMultiVector()->col(0).ptr(), zero);
100 assign(Xdot->getNonconstMultiVector()->col(0).ptr(), *xdot0);
101 if (DxdotDp0 == Teuchos::null)
102 assign(Xdot->getNonconstMultiVector()->subView(rng).ptr(), zero);
104 assign(Xdot->getNonconstMultiVector()->subView(rng).ptr(), *DxdotDp0);
107 if (xdotdot0 == Teuchos::null)
108 assign(Xdotdot->getNonconstMultiVector()->col(0).ptr(), zero);
110 assign(Xdotdot->getNonconstMultiVector()->col(0).ptr(), *xdotdot0);
111 if (DxdotDp0 == Teuchos::null)
112 assign(Xdotdot->getNonconstMultiVector()->subView(rng).ptr(), zero);
114 assign(Xdotdot->getNonconstMultiVector()->subView(rng).ptr(), *DxdotdotDp0);
116 integrator_->initializeSolutionHistory(t0, X, Xdot, Xdotdot);
119template <
class Scalar>
120Teuchos::RCP<const Thyra::VectorBase<Scalar>>
124 using Teuchos::rcp_dynamic_cast;
125 typedef Thyra::DefaultMultiVectorProductVector<Scalar> DMVPV;
127 RCP<const DMVPV> X = rcp_dynamic_cast<const DMVPV>(integrator_->getX());
128 return X->getMultiVector()->col(0);
131template <
class Scalar>
132Teuchos::RCP<const Thyra::MultiVectorBase<Scalar>>
136 using Teuchos::rcp_dynamic_cast;
137 typedef Thyra::DefaultMultiVectorProductVector<Scalar> DMVPV;
139 RCP<const DMVPV> X = rcp_dynamic_cast<const DMVPV>(integrator_->getX());
140 const int num_param = X->getMultiVector()->domain()->dim() - 1;
141 const Teuchos::Range1D rng(1, num_param);
142 return X->getMultiVector()->subView(rng);
145template <
class Scalar>
146Teuchos::RCP<const Thyra::VectorBase<Scalar>>
150 using Teuchos::rcp_dynamic_cast;
151 typedef Thyra::DefaultMultiVectorProductVector<Scalar> DMVPV;
153 RCP<const DMVPV> Xdot = rcp_dynamic_cast<const DMVPV>(integrator_->getXDot());
154 return Xdot->getMultiVector()->col(0);
157template <
class Scalar>
158Teuchos::RCP<const Thyra::MultiVectorBase<Scalar>>
162 using Teuchos::rcp_dynamic_cast;
163 typedef Thyra::DefaultMultiVectorProductVector<Scalar> DMVPV;
165 RCP<const DMVPV> Xdot = rcp_dynamic_cast<const DMVPV>(integrator_->getXDot());
166 const int num_param = Xdot->getMultiVector()->domain()->dim() - 1;
167 const Teuchos::Range1D rng(1, num_param);
168 return Xdot->getMultiVector()->subView(rng);
171template <
class Scalar>
172Teuchos::RCP<const Thyra::VectorBase<Scalar>>
176 using Teuchos::rcp_dynamic_cast;
177 typedef Thyra::DefaultMultiVectorProductVector<Scalar> DMVPV;
179 RCP<const DMVPV> Xdotdot =
180 rcp_dynamic_cast<const DMVPV>(integrator_->getXDotDot());
181 return Xdotdot->getMultiVector()->col(0);
184template <
class Scalar>
185Teuchos::RCP<const Thyra::MultiVectorBase<Scalar>>
189 using Teuchos::rcp_dynamic_cast;
190 typedef Thyra::DefaultMultiVectorProductVector<Scalar> DMVPV;
192 RCP<const DMVPV> Xdotdot =
193 rcp_dynamic_cast<const DMVPV>(integrator_->getXDotDot());
194 const int num_param = Xdotdot->getMultiVector()->domain()->dim() - 1;
195 const Teuchos::Range1D rng(1, num_param);
196 return Xdotdot->getMultiVector()->subView(rng);
199template <
class Scalar>
200Teuchos::RCP<const Thyra::VectorBase<Scalar>>
203 typedef Thyra::ModelEvaluatorBase MEB;
207 Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>> smodel;
208 if (use_combined_method_)
209 smodel = sens_model_;
211 smodel = sens_stepper_->getModel();
212 MEB::InArgs<Scalar> inargs = smodel->getNominalValues();
213 MEB::OutArgs<Scalar> outargs = smodel->createOutArgs();
214 inargs.set_t(integrator_->getTime());
215 inargs.set_x(integrator_->getX());
216 if (inargs.supports(MEB::IN_ARG_x_dot))
217 inargs.set_x_dot(integrator_->getXDot());
218 if (inargs.supports(MEB::IN_ARG_x_dot_dot))
219 inargs.set_x_dot_dot(integrator_->getXDotDot());
221 Teuchos::RCP<Thyra::VectorBase<Scalar>> g =
222 Thyra::createMember(smodel->get_g_space(1));
225 smodel->evalModel(inargs, outargs);
229template <
class Scalar>
230Teuchos::RCP<const Thyra::MultiVectorBase<Scalar>>
233 typedef Thyra::ModelEvaluatorBase MEB;
234 typedef Thyra::DefaultMultiVectorProductVector<Scalar> DMVPV;
238 Teuchos::RCP<const Thyra::ModelEvaluator<Scalar>> smodel;
239 if (use_combined_method_)
240 smodel = sens_model_;
242 smodel = sens_stepper_->getModel();
243 MEB::InArgs<Scalar> inargs = smodel->getNominalValues();
244 MEB::OutArgs<Scalar> outargs = smodel->createOutArgs();
245 inargs.set_t(integrator_->getTime());
246 inargs.set_x(integrator_->getX());
247 if (inargs.supports(MEB::IN_ARG_x_dot))
248 inargs.set_x_dot(integrator_->getXDot());
249 if (inargs.supports(MEB::IN_ARG_x_dot_dot))
250 inargs.set_x_dot_dot(integrator_->getXDotDot());
252 Teuchos::RCP<Thyra::VectorBase<Scalar>> G =
253 Thyra::createMember(smodel->get_g_space(0));
254 Teuchos::RCP<DMVPV> dgdp = Teuchos::rcp_dynamic_cast<DMVPV>(G);
257 smodel->evalModel(inargs, outargs);
258 return dgdp->getMultiVector();
261template <
class Scalar>
264 std::string name =
"Tempus::IntegratorForwardSensitivity";
268template <
class Scalar>
270 Teuchos::FancyOStream &out,
const Teuchos::EVerbosityLevel verbLevel)
const
272 auto l_out = Teuchos::fancyOStream(out.getOStream());
273 Teuchos::OSTab ostab(*l_out, 2, this->description());
274 l_out->setOutputToRootOnly(0);
276 *l_out << description() <<
"::describe" << std::endl;
277 integrator_->describe(*l_out, verbLevel);
280template <
class Scalar>
284 return sens_stepper_->getStepMode();
288template <
class Scalar>
289Teuchos::RCP<IntegratorForwardSensitivity<Scalar>>
291 Teuchos::RCP<Teuchos::ParameterList> pList,
296 Teuchos::RCP<SensitivityModelEvaluatorBase<Scalar>> sens_model;
297 Teuchos::RCP<StepperStaggeredForwardSensitivity<Scalar>> sens_stepper;
300 auto fwd_integrator = createIntegratorBasic<Scalar>(pList, model);
303 Teuchos::RCP<Teuchos::ParameterList> pl =
304 Teuchos::rcp(
new Teuchos::ParameterList);
306 Teuchos::RCP<const Teuchos::ParameterList> integrator_pl =
307 fwd_integrator->getValidParameters();
308 Teuchos::RCP<const Teuchos::ParameterList> sensitivity_pl =
310 pl->setParameters(*integrator_pl);
311 Teuchos::ParameterList &spl = pl->sublist(
"Sensitivities");
312 spl.setParameters(*sensitivity_pl);
313 spl.set(
"Sensitivity Method",
"Combined");
314 spl.set(
"Reuse State Linear Solver",
false);
316 pList->setParametersNotAlreadySet(*pl);
318 auto sens_pl = Teuchos::sublist(pList,
"Sensitivities",
false);
319 std::string integratorName =
320 pList->get<std::string>(
"Integrator Name",
"Default Integrator");
321 std::string stepperName =
322 pList->sublist(integratorName).get<std::string>(
"Stepper Name");
323 auto stepper_pl = Teuchos::sublist(pList, stepperName,
true);
324 std::string sensitivity_method =
325 sens_pl->get<std::string>(
"Sensitivity Method");
326 bool use_combined_method = sensitivity_method ==
"Combined";
331 sens_pl->remove(
"Sensitivity Method");
333 if (use_combined_method) {
334 sens_pl->remove(
"Reuse State Linear Solver");
336 sens_solve_model, sens_pl);
337 fwd_integrator = createIntegratorBasic<Scalar>(pList, sens_model);
342 model, sens_residual_model, sens_solve_model, stepper_pl,
344 auto fsa_staggered_me =
345 Teuchos::rcp_const_cast<Thyra::ModelEvaluator<Scalar>>(
346 sens_stepper->getModel());
347 fwd_integrator = createIntegratorBasic<Scalar>(pList, fsa_staggered_me);
348 fwd_integrator->setStepper(sens_stepper);
353 Teuchos::RCP<IntegratorForwardSensitivity<Scalar>> integrator =
355 model, fwd_integrator, sens_model, sens_stepper,
356 use_combined_method));
361template <
class Scalar>
362Teuchos::RCP<IntegratorForwardSensitivity<Scalar>>
365 Teuchos::RCP<IntegratorForwardSensitivity<Scalar>> integrator =
static Teuchos::RCP< const Teuchos::ParameterList > getValidParameters()
Time integrator implementing forward sensitivity analysis.
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDXDotDp() const
virtual void setStepper(Teuchos::RCP< Thyra::ModelEvaluator< Scalar > > model)
Set the Stepper.
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXDotDot() const
Get current the second time derivative of the solution, xdotdot, only. This is the first block only a...
SensitivityStepMode getStepMode() const
What mode the current time integration step is in.
virtual void initializeSolutionHistory(Teuchos::RCP< SolutionState< Scalar > > state=Teuchos::null)
Set the initial state which has the initial conditions.
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDxDp() const
Get the forward sensitivities .
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getX() const
Get the current solution, x, only. If looking for the solution vector and the sensitivities,...
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDXDotDotDp() const
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getG() const
Return response function g.
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDgDp() const
Return forward sensitivity stored in Jacobian format.
Teuchos::RCP< IntegratorBasic< Scalar > > integrator_
std::string description() const override
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getXDot() const
Get current the time derivative of the solution, xdot, only. This is the first block only and not the...
IntegratorForwardSensitivity()
Destructor.
A ModelEvaluator decorator for sensitivity analysis.
A stepper implementing staggered forward sensitivity analysis.
Teuchos::RCP< SensitivityModelEvaluatorBase< Scalar > > wrapCombinedFSAModelEvaluator(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 Teuchos::RCP< const Teuchos::ParameterList > &pList=Teuchos::null)
Teuchos::RCP< IntegratorForwardSensitivity< Scalar > > createIntegratorForwardSensitivity()
Nonmember constructor.