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"
23 template <
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)
40 template <
class Scalar>
43 integrator_ = createIntegratorBasic<Scalar>();
44 integrator_->initialize();
47 template <
class Scalar>
51 if (use_combined_method_)
52 integrator_->setModel(sens_model_);
54 integrator_->setStepper(sens_stepper_);
57 template <
class Scalar>
67 using Teuchos::rcp_dynamic_cast;
69 using Thyra::createMember;
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));
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);
119 template <
class Scalar>
124 using Teuchos::rcp_dynamic_cast;
127 RCP<const DMVPV> X = rcp_dynamic_cast<
const DMVPV>(integrator_->getX());
131 template <
class Scalar>
136 using Teuchos::rcp_dynamic_cast;
139 RCP<const DMVPV> X = rcp_dynamic_cast<
const DMVPV>(integrator_->getX());
142 return X->getMultiVector()->subView(rng);
145 template <
class Scalar>
150 using Teuchos::rcp_dynamic_cast;
153 RCP<const DMVPV> Xdot = rcp_dynamic_cast<
const DMVPV>(integrator_->getXDot());
157 template <
class Scalar>
162 using Teuchos::rcp_dynamic_cast;
165 RCP<const DMVPV> Xdot = rcp_dynamic_cast<
const DMVPV>(integrator_->getXDot());
166 const int num_param = Xdot->
getMultiVector()->domain()->dim() - 1;
168 return Xdot->getMultiVector()->subView(rng);
171 template <
class Scalar>
176 using Teuchos::rcp_dynamic_cast;
179 RCP<const DMVPV> Xdotdot =
180 rcp_dynamic_cast<
const DMVPV>(integrator_->getXDotDot());
184 template <
class Scalar>
189 using Teuchos::rcp_dynamic_cast;
192 RCP<const DMVPV> Xdotdot =
193 rcp_dynamic_cast<
const DMVPV>(integrator_->getXDotDot());
194 const int num_param = Xdotdot->
getMultiVector()->domain()->dim() - 1;
196 return Xdotdot->getMultiVector()->subView(rng);
199 template <
class Scalar>
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());
222 Thyra::createMember(smodel->get_g_space(1));
225 smodel->evalModel(inargs, outargs);
229 template <
class Scalar>
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());
253 Thyra::createMember(smodel->get_g_space(0));
257 smodel->evalModel(inargs, outargs);
258 return dgdp->getMultiVector();
261 template <
class Scalar>
264 std::string name =
"Tempus::IntegratorForwardSensitivity";
268 template <
class Scalar>
272 auto l_out = Teuchos::fancyOStream(out.
getOStream());
274 l_out->setOutputToRootOnly(0);
276 *l_out << description() <<
"::describe" << std::endl;
277 integrator_->describe(*l_out, verbLevel);
280 template <
class Scalar>
284 return sens_stepper_->getStepMode();
288 template <
class Scalar>
300 auto fwd_integrator = createIntegratorBasic<Scalar>(pList, model);
307 fwd_integrator->getValidParameters();
310 pl->setParameters(*integrator_pl);
313 spl.
set(
"Sensitivity Method",
"Combined");
314 spl.
set(
"Reuse State Linear Solver",
false);
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 =
346 sens_stepper->getModel());
347 fwd_integrator = createIntegratorBasic<Scalar>(pList, fsa_staggered_me);
348 fwd_integrator->setStepper(sens_stepper);
355 model, fwd_integrator, sens_model, sens_stepper,
356 use_combined_method));
361 template <
class Scalar>
371 #endif // Tempus_IntegratorForwardSensitivity_impl_hpp
Teuchos::RCP< IntegratorBasic< Scalar > > integrator_
Teuchos::RCP< IntegratorForwardSensitivity< Scalar > > createIntegratorForwardSensitivity(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.
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...
std::string description() const override
virtual void setStepper(Teuchos::RCP< Thyra::ModelEvaluator< Scalar >> model)
Set the Stepper.
SensitivityStepMode getStepMode() const
What mode the current time integration step is in.
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDXDotDp() const
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDXDotDotDp() const
static Teuchos::RCP< const Teuchos::ParameterList > getValidParameters()
A stepper implementing staggered forward sensitivity analysis.
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDgDp() const
Return forward sensitivity stored in Jacobian format.
IntegratorForwardSensitivity()
Destructor.
virtual Teuchos::RCP< const Thyra::MultiVectorBase< Scalar > > getDxDp() const
Get the forward sensitivities .
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
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)
RCP< const MultiVectorBase< Scalar > > getMultiVector() const
void describe(Teuchos::FancyOStream &out, const Teuchos::EVerbosityLevel verbLevel) const override
ParameterList & setParameters(const ParameterList &source)
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getG() const
Return response function g.
virtual Teuchos::RCP< const Thyra::VectorBase< Scalar > > getX() const
Get the current solution, x, only. If looking for the solution vector and the sensitivities, use SolutionState->getX() which will return a Block MultiVector with the first block containing the current solution, x, and the remaining blocks are the forward sensitivities .
RCP< std::basic_ostream< char_type, traits_type > > getOStream()
ParameterList & setParametersNotAlreadySet(const ParameterList &source)
RCP< MultiVectorBase< Scalar > > getNonconstMultiVector()
ParameterList & sublist(const std::string &name, bool mustAlreadyExist=false, const std::string &docString="")
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...
A ModelEvaluator decorator for sensitivity analysis.
Time integrator implementing forward sensitivity analysis.
virtual void initializeSolutionHistory(Teuchos::RCP< SolutionState< Scalar >> state=Teuchos::null)
Set the initial state which has the initial conditions.