Tempus  Version of the Day
Time Integration
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Tempus_StaggeredForwardSensitivityModelEvaluator_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_StaggeredForwardSensitivityModelEvaluator_impl_hpp
10 #define Tempus_StaggeredForwardSensitivityModelEvaluator_impl_hpp
11 
12 #include "Thyra_DefaultMultiVectorProductVector.hpp"
16 #include "Thyra_VectorStdOps.hpp"
17 #include "Thyra_MultiVectorStdOps.hpp"
18 
19 namespace Tempus {
20 
21 template <typename Scalar>
24  const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& model,
26  sens_residual_model,
28  sens_solve_model,
29  const bool is_pseudotransient,
31  const Teuchos::RCP<MultiVector>& dxdp_init,
32  const Teuchos::RCP<MultiVector>& dx_dotdp_init,
33  const Teuchos::RCP<MultiVector>& dx_dotdotdp_init)
34  : model_(model),
35  sens_residual_model_(sens_residual_model),
36  sens_solve_model_(sens_solve_model),
37  dxdp_init_(dxdp_init),
38  dx_dotdp_init_(dx_dotdp_init),
39  dx_dotdotdp_init_(dx_dotdotdp_init),
40  p_index_(0),
41  g_index_(-1),
42  x_tangent_index_(1),
43  xdot_tangent_index_(2),
44  xdotdot_tangent_index_(3),
45  use_dfdp_as_tangent_(false),
46  use_dgdp_as_tangent_(false),
47  num_param_(0),
48  num_response_(0),
49  g_offset_(0),
50  is_pseudotransient_(is_pseudotransient),
51  mass_matrix_is_computed_(false),
52  jacobian_matrix_is_computed_(false),
53  acceleration_matrix_is_computed_(false),
54  residual_sensitivity_is_computed_(false),
55  t_interp_(Teuchos::ScalarTraits<Scalar>::rmax())
56 {
57  typedef Thyra::ModelEvaluatorBase MEB;
58 
59  // Set parameters
62  if (pList != Teuchos::null) *pl = *pList;
64  use_dfdp_as_tangent_ = pl->get<bool>("Use DfDp as Tangent");
65  use_dgdp_as_tangent_ = pl->get<bool>("Use DgDp as Tangent");
66  p_index_ = pl->get<int>("Sensitivity Parameter Index");
67  g_index_ = pl->get<int>("Response Function Index", -1);
68  x_tangent_index_ = pl->get<int>("Sensitivity X Tangent Index");
69  xdot_tangent_index_ = pl->get<int>("Sensitivity X-Dot Tangent Index");
70  xdotdot_tangent_index_ = pl->get<int>("Sensitivity X-Dot-Dot Tangent Index");
71 
72  num_param_ = model_->get_p_space(p_index_)->dim();
73  if (g_index_ >= 0) {
74  num_response_ = model_->get_g_space(g_index_)->dim();
75  g_offset_ = 2;
76  }
77  dxdp_space_ =
78  Thyra::multiVectorProductVectorSpace(model_->get_x_space(), num_param_);
79  dfdp_space_ =
80  Thyra::multiVectorProductVectorSpace(model_->get_f_space(), num_param_);
81  if (g_index_ >= 0)
82  dgdp_space_ = Thyra::multiVectorProductVectorSpace(
83  model_->get_g_space(g_index_), num_param_);
84 
85  // forward and sensitivity models must support same InArgs
86  MEB::InArgs<Scalar> me_inArgs = model_->createInArgs();
87  me_inArgs.assertSameSupport(sens_residual_model_->createInArgs());
88  me_inArgs.assertSameSupport(sens_solve_model_->createInArgs());
89 
90  MEB::InArgsSetup<Scalar> inArgs;
91  inArgs.setModelEvalDescription(this->description());
92  inArgs.setSupports(MEB::IN_ARG_x);
93  if (me_inArgs.supports(MEB::IN_ARG_x_dot))
94  inArgs.setSupports(MEB::IN_ARG_x_dot);
95  if (me_inArgs.supports(MEB::IN_ARG_t)) inArgs.setSupports(MEB::IN_ARG_t);
96  if (me_inArgs.supports(MEB::IN_ARG_alpha))
97  inArgs.setSupports(MEB::IN_ARG_alpha);
98  if (me_inArgs.supports(MEB::IN_ARG_beta))
99  inArgs.setSupports(MEB::IN_ARG_beta);
100  if (me_inArgs.supports(MEB::IN_ARG_W_x_dot_dot_coeff))
101  inArgs.setSupports(MEB::IN_ARG_W_x_dot_dot_coeff);
102 
103  // Support additional parameters for x and xdot
104  inArgs.set_Np(me_inArgs.Np());
105  prototypeInArgs_ = inArgs;
106 
107  MEB::OutArgs<Scalar> me_outArgs = model_->createOutArgs();
108  MEB::OutArgs<Scalar> sens_mer_outArgs = sens_residual_model_->createOutArgs();
109  MEB::OutArgs<Scalar> sens_mes_outArgs = sens_solve_model_->createOutArgs();
110  MEB::OutArgsSetup<Scalar> outArgs;
111  outArgs.setModelEvalDescription(this->description());
112  outArgs.set_Np_Ng(me_inArgs.Np(), me_outArgs.Ng() + g_offset_);
113  outArgs.setSupports(MEB::OUT_ARG_f);
114  if (sens_mer_outArgs.supports(MEB::OUT_ARG_W_op) ||
115  sens_mes_outArgs.supports(MEB::OUT_ARG_W_op))
116  outArgs.setSupports(MEB::OUT_ARG_W_op);
117 
118  // Response 0 is the reduced dg/dp (no sensitivities supported)
119  // Response 1 is the reduced g (no sensitivities supported)
120  for (int j = 0; j < me_outArgs.Ng(); ++j) {
121  outArgs.setSupports(MEB::OUT_ARG_DgDx_dot, j + g_offset_,
122  me_outArgs.supports(MEB::OUT_ARG_DgDx_dot, j));
123  outArgs.setSupports(MEB::OUT_ARG_DgDx, j + g_offset_,
124  me_outArgs.supports(MEB::OUT_ARG_DgDx, j));
125  for (int l = 0; l < me_outArgs.Np(); ++l) {
126  outArgs.setSupports(MEB::OUT_ARG_DgDp, j + g_offset_, l,
127  me_outArgs.supports(MEB::OUT_ARG_DgDp, j, l));
128  }
129  }
130  prototypeOutArgs_ = outArgs;
131 
132  // Sensitivity residual ME must support W_op to define adjoint ODE/DAE.
133  // Must support alpha, beta if it suports x_dot
134  TEUCHOS_ASSERT(me_inArgs.supports(MEB::IN_ARG_x));
136  TEUCHOS_ASSERT(sens_mer_outArgs.supports(MEB::OUT_ARG_W_op));
137  if (me_inArgs.supports(MEB::IN_ARG_x_dot)) {
138  TEUCHOS_ASSERT(me_inArgs.supports(MEB::IN_ARG_alpha));
139  TEUCHOS_ASSERT(me_inArgs.supports(MEB::IN_ARG_beta));
140  }
141  TEUCHOS_ASSERT(me_outArgs.supports(MEB::OUT_ARG_DfDp, p_index_)
142  .supports(MEB::DERIV_MV_JACOBIAN_FORM));
143 }
144 
145 template <typename Scalar>
149 {
150  sh_ = sh;
152 }
153 
154 template <typename Scalar>
157 {
158  sh_ = Teuchos::null;
159  forward_state_ = s;
160 }
161 
162 template <typename Scalar>
165 {
166  TEUCHOS_ASSERT(p < model_->Np());
167  return model_->get_p_space(p);
168 }
169 
170 template <typename Scalar>
173 {
174  TEUCHOS_ASSERT(p < model_->Np());
175  return model_->get_p_names(p);
176 }
177 
178 template <typename Scalar>
181 {
182  return dxdp_space_;
183 }
184 
185 template <typename Scalar>
188 {
189  return dfdp_space_;
190 }
191 
192 template <typename Scalar>
195 {
196  if (g_index_ >= 0) {
197  if (j == 0)
198  return dgdp_space_;
199  else if (j == 1)
200  return model_->get_g_space(g_index_);
201  }
202  return model_->get_g_space(j - g_offset_);
203 }
204 
205 template <typename Scalar>
208 {
209  if (g_index_ >= 0) {
210  if (j == 0) {
211  Teuchos::Array<std::string> names = model_->get_g_names(g_index_);
212  for (int i = 0; i < names.size(); ++i)
213  names[i] = names[i] + "_reduced sensitivity";
214  return names();
215  }
216  else if (j == 1)
217  return model_->get_g_names(g_index_);
218  }
219  return model_->get_g_names(j - g_offset_);
220 }
221 
222 template <typename Scalar>
225 {
227  if (lo_ != Teuchos::null)
228  op = lo_;
229  else
230  op = sens_solve_model_->create_W_op();
231  return Thyra::nonconstMultiVectorLinearOp(op, num_param_);
232 }
233 
234 template <typename Scalar>
237  int j) const
238 {
239  return model_->create_DgDx_dot_op(j - g_offset_);
240 }
241 
242 template <typename Scalar>
245 {
246  return model_->create_DgDx_op(j - g_offset_);
247 }
248 
249 template <typename Scalar>
252  int l) const
253 {
254  return model_->create_DgDp_op(j - g_offset_, l);
255 }
256 
257 template <typename Scalar>
260 {
261  using Teuchos::RCP;
262  using Teuchos::rcp_dynamic_cast;
263 
266  model_factory = sens_solve_model_->get_W_factory();
267  if (model_factory == Teuchos::null)
268  return Teuchos::null; // model_ doesn't support W_factory
269  if (po_ != Teuchos::null) {
270  factory = Thyra::reuseLinearOpWithSolveFactory<Scalar>(model_factory, po_);
271  }
272  else
273  factory = model_factory;
274  return Thyra::multiVectorLinearOpWithSolveFactory(factory, dfdp_space_,
275  dxdp_space_);
276 }
277 
278 template <typename Scalar>
281 {
282  return prototypeInArgs_;
283 }
284 
285 template <typename Scalar>
288 {
289  typedef Thyra::ModelEvaluatorBase MEB;
291  using Teuchos::RCP;
292  using Teuchos::rcp_dynamic_cast;
293 
294  MEB::InArgs<Scalar> me_nominal = model_->getNominalValues();
295  MEB::InArgs<Scalar> nominal = this->createInArgs();
296 
297  const Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
298 
299  // Set initial x. If dxdp_init == null, set initial dx/dp = 0
300  RCP<Thyra::VectorBase<Scalar> > x = Thyra::createMember(*dxdp_space_);
301  RCP<DMVPV> dxdp = rcp_dynamic_cast<DMVPV>(x, true);
302  if (dxdp_init_ == Teuchos::null)
303  Thyra::assign(dxdp->getNonconstMultiVector().ptr(), zero);
304  else
305  Thyra::assign(dxdp->getNonconstMultiVector().ptr(), *dxdp_init_);
306  nominal.set_x(x);
307 
308  // Set initial xdot. If dx_dotdp_init == null, set initial dxdot/dp = 0
309  if (me_nominal.supports(MEB::IN_ARG_x_dot)) {
310  RCP<Thyra::VectorBase<Scalar> > xdot = Thyra::createMember(*dxdp_space_);
311  RCP<DMVPV> dxdotdp = rcp_dynamic_cast<DMVPV>(xdot, true);
312  if (dx_dotdp_init_ == Teuchos::null)
313  Thyra::assign(dxdotdp->getNonconstMultiVector().ptr(), zero);
314  else
315  Thyra::assign(dxdotdp->getNonconstMultiVector().ptr(), *dx_dotdp_init_);
316  nominal.set_x_dot(xdot);
317  }
318 
319  // Set initial xdotdot. If dx_dotdotdp_init == null, set initial dxdotdot/dp
320  // = 0
321  if (me_nominal.supports(MEB::IN_ARG_x_dot_dot)) {
322  RCP<Thyra::VectorBase<Scalar> > xdotdot = Thyra::createMember(*dxdp_space_);
323  RCP<DMVPV> dxdotdotdp = rcp_dynamic_cast<DMVPV>(xdotdot, true);
324  if (dx_dotdotdp_init_ == Teuchos::null)
325  Thyra::assign(dxdotdotdp->getNonconstMultiVector().ptr(), zero);
326  else
327  Thyra::assign(dxdotdotdp->getNonconstMultiVector().ptr(),
328  *dx_dotdotdp_init_);
329  nominal.set_x_dot_dot(xdotdot);
330  }
331 
332  const int np = model_->Np();
333  for (int i = 0; i < np; ++i) nominal.set_p(i, me_nominal.get_p(i));
334  return nominal;
335 }
336 
337 template <typename Scalar>
340 {
341  return prototypeOutArgs_;
342 }
343 
344 template <typename Scalar>
347  const Thyra::ModelEvaluatorBase::OutArgs<Scalar>& outArgs) const
348 {
350  typedef Thyra::ModelEvaluatorBase MEB;
351  using Teuchos::RCP;
352  using Teuchos::rcp_dynamic_cast;
353 
354  // Interpolate forward solution at supplied time, reusing previous
355  // interpolation if possible
356  Scalar forward_t;
357  if (sh_ != Teuchos::null) {
358  forward_t = inArgs.get_t();
359  if (t_interp_ != forward_t) {
360  if (nc_forward_state_ == Teuchos::null)
361  nc_forward_state_ = sh_->interpolateState(forward_t);
362  else
363  sh_->interpolateState(forward_t, nc_forward_state_.get());
364  forward_state_ = nc_forward_state_;
365  t_interp_ = forward_t;
366  }
367  }
368  else {
369  TEUCHOS_ASSERT(forward_state_ != Teuchos::null);
370  forward_t = forward_state_->getTime();
371  }
372 
373  const bool use_tangent = use_dfdp_as_tangent_ || use_dgdp_as_tangent_;
374 
375  // setup input arguments for model
376  RCP<const Thyra::MultiVectorBase<Scalar> > dxdp, dxdotdp, dxdotdotdp;
377  MEB::InArgs<Scalar> me_inArgs = model_->getNominalValues();
378  dxdp = rcp_dynamic_cast<const DMVPV>(inArgs.get_x(), true)->getMultiVector();
379  me_inArgs.set_x(forward_state_->getX());
380  if (use_dfdp_as_tangent_) me_inArgs.set_p(x_tangent_index_, inArgs.get_x());
381  if (me_inArgs.supports(MEB::IN_ARG_x_dot)) {
382  if (inArgs.get_x_dot() != Teuchos::null) {
383  dxdotdp = rcp_dynamic_cast<const DMVPV>(inArgs.get_x_dot(), true)
384  ->getMultiVector();
385  me_inArgs.set_x_dot(forward_state_->getXDot());
386  if (use_dfdp_as_tangent_)
387  me_inArgs.set_p(xdot_tangent_index_, inArgs.get_x_dot());
388  }
389  else {
390  // clear out xdot if it was set in nominalValues to get to ensure we
391  // get the explicit form
392  me_inArgs.set_x_dot(Teuchos::null);
393  }
394  }
395  if (me_inArgs.supports(MEB::IN_ARG_x_dot_dot)) {
396  if (inArgs.get_x_dot_dot() != Teuchos::null) {
397  dxdotdotdp = rcp_dynamic_cast<const DMVPV>(inArgs.get_x_dot_dot(), true)
398  ->getMultiVector();
399  me_inArgs.set_x_dot_dot(forward_state_->getXDotDot());
400  if (use_dfdp_as_tangent_)
401  me_inArgs.set_p(xdotdot_tangent_index_, inArgs.get_x_dot_dot());
402  }
403  else // clear out xdotdot if it was set in nominalValues
404  me_inArgs.set_x_dot_dot(Teuchos::null);
405  }
406  if (me_inArgs.supports(MEB::IN_ARG_t)) me_inArgs.set_t(forward_t);
407  if (me_inArgs.supports(MEB::IN_ARG_alpha))
408  me_inArgs.set_alpha(inArgs.get_alpha());
409  if (me_inArgs.supports(MEB::IN_ARG_beta))
410  me_inArgs.set_beta(inArgs.get_beta());
411  if (me_inArgs.supports(MEB::IN_ARG_W_x_dot_dot_coeff))
412  me_inArgs.set_W_x_dot_dot_coeff(inArgs.get_W_x_dot_dot_coeff());
413 
414  // Set parameters -- be careful to only set ones that were set in our
415  // inArgs to not null out any specified through nominalValues or
416  // dx/dp above
417  const int np = me_inArgs.Np();
418  for (int i = 0; i < np; ++i) {
419  if (inArgs.get_p(i) != Teuchos::null)
420  if (!use_tangent ||
421  (use_tangent && i != x_tangent_index_ && i != xdot_tangent_index_ &&
422  i != xdotdot_tangent_index_))
423  me_inArgs.set_p(i, inArgs.get_p(i));
424  }
425 
426  // setup output arguments for model
427  RCP<Thyra::MultiVectorBase<Scalar> > dfdp;
428  MEB::OutArgs<Scalar> me_outArgs = model_->createOutArgs();
429  if (outArgs.get_f() != Teuchos::null) {
430  dfdp = rcp_dynamic_cast<DMVPV>(outArgs.get_f(), true)
431  ->getNonconstMultiVector();
432  if (!residual_sensitivity_is_computed_) {
433  if (!use_dfdp_as_tangent_ && is_pseudotransient_) {
434  if (my_dfdp_ == Teuchos::null)
435  my_dfdp_ = Thyra::createMembers(model_->get_f_space(),
436  model_->get_p_space(p_index_)->dim());
437  me_outArgs.set_DfDp(p_index_, my_dfdp_);
438  }
439  else
440  me_outArgs.set_DfDp(p_index_, dfdp);
441  }
442  }
443  if (outArgs.supports(MEB::OUT_ARG_W_op) && lo_ == Teuchos::null &&
444  outArgs.get_W_op() != Teuchos::null &&
445  model_.ptr() == sens_solve_model_.ptr()) {
446  RCP<Thyra::LinearOpBase<Scalar> > op = outArgs.get_W_op();
447  RCP<Thyra::MultiVectorLinearOp<Scalar> > mv_op =
448  rcp_dynamic_cast<Thyra::MultiVectorLinearOp<Scalar> >(op, true);
449  me_outArgs.set_W_op(mv_op->getNonconstLinearOp());
450  }
451  for (int j = g_offset_; j < outArgs.Ng(); ++j) {
452  me_outArgs.set_g(j - g_offset_, outArgs.get_g(j));
453  if (!me_outArgs.supports(MEB::OUT_ARG_DgDx_dot, j - g_offset_).none())
454  me_outArgs.set_DgDx_dot(j - g_offset_, outArgs.get_DgDx_dot(j));
455  if (!me_outArgs.supports(MEB::OUT_ARG_DgDx, j - g_offset_).none())
456  me_outArgs.set_DgDx(j - g_offset_, outArgs.get_DgDx(j));
457  for (int l = 0; l < outArgs.Np(); ++l)
458  if (!me_outArgs.supports(MEB::OUT_ARG_DgDp, j - g_offset_, l).none())
459  me_outArgs.set_DgDp(j - g_offset_, l, outArgs.get_DgDp(j, l));
460  }
461  if (g_index_ >= 0 && outArgs.get_g(1) != Teuchos::null)
462  me_outArgs.set_g(g_index_, outArgs.get_g(1));
463 
464  // build residual and jacobian
465  model_->evalModel(me_inArgs, me_outArgs);
466 
467  // Compute W_op separately if we have a separate sensitivity solve ME
468  if (outArgs.supports(MEB::OUT_ARG_W_op) && lo_ == Teuchos::null &&
469  outArgs.get_W_op() != Teuchos::null &&
470  model_.ptr() != sens_solve_model_.ptr()) {
471  MEB::OutArgs<Scalar> sens_me_outArgs = sens_solve_model_->createOutArgs();
472  RCP<Thyra::LinearOpBase<Scalar> > op = outArgs.get_W_op();
473  RCP<Thyra::MultiVectorLinearOp<Scalar> > mv_op =
474  rcp_dynamic_cast<Thyra::MultiVectorLinearOp<Scalar> >(op, true);
475  sens_me_outArgs.set_W_op(mv_op->getNonconstLinearOp());
476  sens_solve_model_->evalModel(me_inArgs, sens_me_outArgs);
477  }
478 
479  // Compute (df/dx) * (dx/dp) + (df/dxdot) * (dxdot/dp) + (df/dxdotdot) *
480  // (dxdotdot/dp) + (df/dp) if the underlying ME doesn't already do this. This
481  // requires computing df/dx, df/dxdot, df/dxdotdot as separate operators. For
482  // pseudo-transient, we would like to reuse these operators, but this is
483  // complicated when steppers use both implicit and explicit forms.
484  if (!use_dfdp_as_tangent_) {
485  if (dfdp != Teuchos::null && is_pseudotransient_) {
486  residual_sensitivity_is_computed_ = true;
487  Thyra::assign(dfdp.ptr(), *my_dfdp_);
488  }
489  if (dxdp != Teuchos::null && dfdp != Teuchos::null) {
490  if (my_dfdx_ == Teuchos::null)
491  my_dfdx_ = sens_residual_model_->create_W_op();
492  if (!jacobian_matrix_is_computed_) {
493  MEB::OutArgs<Scalar> meo = sens_residual_model_->createOutArgs();
494  meo.set_W_op(my_dfdx_);
495  if (me_inArgs.supports(MEB::IN_ARG_alpha)) me_inArgs.set_alpha(0.0);
496  if (me_inArgs.supports(MEB::IN_ARG_beta)) me_inArgs.set_beta(1.0);
497  if (me_inArgs.supports(MEB::IN_ARG_W_x_dot_dot_coeff))
498  me_inArgs.set_W_x_dot_dot_coeff(0.0);
499  sens_residual_model_->evalModel(me_inArgs, meo);
500  if (is_pseudotransient_) jacobian_matrix_is_computed_ = true;
501  }
502  my_dfdx_->apply(Thyra::NOTRANS, *dxdp, dfdp.ptr(), Scalar(1.0),
503  Scalar(1.0));
504  }
505  if (dxdotdp != Teuchos::null && dfdp != Teuchos::null) {
506  if (my_dfdxdot_ == Teuchos::null)
507  my_dfdxdot_ = sens_residual_model_->create_W_op();
508  if (!mass_matrix_is_computed_) {
509  MEB::OutArgs<Scalar> meo = sens_residual_model_->createOutArgs();
510  meo.set_W_op(my_dfdxdot_);
511  if (me_inArgs.supports(MEB::IN_ARG_alpha)) me_inArgs.set_alpha(1.0);
512  if (me_inArgs.supports(MEB::IN_ARG_beta)) me_inArgs.set_beta(0.0);
513  if (me_inArgs.supports(MEB::IN_ARG_W_x_dot_dot_coeff))
514  me_inArgs.set_W_x_dot_dot_coeff(0.0);
515  sens_residual_model_->evalModel(me_inArgs, meo);
516  if (is_pseudotransient_) mass_matrix_is_computed_ = true;
517  }
518  my_dfdxdot_->apply(Thyra::NOTRANS, *dxdotdp, dfdp.ptr(), Scalar(1.0),
519  Scalar(1.0));
520  }
521  if (dxdotdotdp != Teuchos::null && dfdp != Teuchos::null) {
522  if (my_dfdxdotdot_ == Teuchos::null)
523  my_dfdxdotdot_ = sens_residual_model_->create_W_op();
524  if (!acceleration_matrix_is_computed_) {
525  MEB::OutArgs<Scalar> meo = sens_residual_model_->createOutArgs();
526  meo.set_W_op(my_dfdxdotdot_);
527  if (me_inArgs.supports(MEB::IN_ARG_alpha)) me_inArgs.set_alpha(0.0);
528  if (me_inArgs.supports(MEB::IN_ARG_beta)) me_inArgs.set_beta(0.0);
529  if (me_inArgs.supports(MEB::IN_ARG_W_x_dot_dot_coeff))
530  me_inArgs.set_W_x_dot_dot_coeff(1.0);
531  sens_residual_model_->evalModel(me_inArgs, meo);
532  if (is_pseudotransient_) acceleration_matrix_is_computed_ = true;
533  }
534  my_dfdxdotdot_->apply(Thyra::NOTRANS, *dxdotdotdp, dfdp.ptr(),
535  Scalar(1.0), Scalar(1.0));
536  }
537  }
538 
539  if (g_index_ >= 0 && outArgs.get_g(0) != Teuchos::null) {
540  MEB::OutArgs<Scalar> meo = model_->createOutArgs();
541  RCP<Thyra::MultiVectorBase<Scalar> > dgdp =
542  rcp_dynamic_cast<DMVPV>(outArgs.get_g(0), true)
543  ->getNonconstMultiVector();
544  RCP<Thyra::MultiVectorBase<Scalar> > dgdp_trans;
545  MEB::DerivativeSupport dgdp_support =
546  meo.supports(MEB::OUT_ARG_DgDp, g_index_, p_index_);
547  if (dgdp_support.supports(MEB::DERIV_MV_JACOBIAN_FORM))
548  meo.set_DgDp(g_index_, p_index_,
549  MEB::Derivative<Scalar>(dgdp, MEB::DERIV_MV_JACOBIAN_FORM));
550  else if (dgdp_support.supports(MEB::DERIV_MV_GRADIENT_FORM)) {
551  dgdp_trans = createMembers(model_->get_p_space(p_index_), num_response_);
552  meo.set_DgDp(
553  g_index_, p_index_,
554  MEB::Derivative<Scalar>(dgdp_trans, MEB::DERIV_MV_GRADIENT_FORM));
555  }
556  else
558  true, std::logic_error,
559  "Operator form of dg/dp not supported for reduced sensitivity");
560 
561  if (!use_dgdp_as_tangent_ || dxdp == Teuchos::null) {
562  MEB::DerivativeSupport dgdx_support =
563  meo.supports(MEB::OUT_ARG_DgDx, g_index_);
564  if (dgdx_support.supports(MEB::DERIV_LINEAR_OP)) {
565  if (my_dgdx_ == Teuchos::null)
566  my_dgdx_ = model_->create_DgDx_op(g_index_);
567  meo.set_DgDx(g_index_, my_dgdx_);
568  }
569  else if (dgdx_support.supports(MEB::DERIV_MV_GRADIENT_FORM)) {
570  if (my_dgdx_mv_ == Teuchos::null)
571  my_dgdx_mv_ =
572  createMembers(model_->get_g_space(g_index_), num_param_);
573  meo.set_DgDx(g_index_, MEB::Derivative<Scalar>(
574  my_dgdx_mv_, MEB::DERIV_MV_GRADIENT_FORM));
575  }
576  else
578  true, std::logic_error,
579  "Jacobian form of dg/dx not supported for reduced sensitivity");
580  }
581 
582  // Clear dx/dp, dxdot/dp, dxdotdot/dp from inArgs if set and we are not
583  // using dg/dp as the tangent. Note, even though dxdp, dxdotdp, and
584  // dxdotdotdp are not used here, they are non-null if x/xdot/xdotdot are
585  // supported and supplied.
586  if (!use_dgdp_as_tangent_ && use_dfdp_as_tangent_ && dxdp != Teuchos::null)
587  me_inArgs.set_p(x_tangent_index_, Teuchos::null);
588  if (!use_dgdp_as_tangent_ && use_dfdp_as_tangent_ &&
589  dxdotdp != Teuchos::null)
590  me_inArgs.set_p(xdot_tangent_index_, Teuchos::null);
591  if (!use_dgdp_as_tangent_ && use_dfdp_as_tangent_ &&
592  dxdotdotdp != Teuchos::null)
593  me_inArgs.set_p(xdotdot_tangent_index_, Teuchos::null);
594 
595  // Set dx/dp, dxdot/dp, dxdodot/dp if necessary
596  if (use_dgdp_as_tangent_ && !use_dfdp_as_tangent_ && dxdp != Teuchos::null)
597  me_inArgs.set_p(x_tangent_index_, inArgs.get_x());
598  if (use_dgdp_as_tangent_ && !use_dfdp_as_tangent_ &&
599  dxdotdp != Teuchos::null)
600  me_inArgs.set_p(xdot_tangent_index_, inArgs.get_x_dot());
601  if (use_dgdp_as_tangent_ && !use_dfdp_as_tangent_ &&
602  dxdotdotdp != Teuchos::null)
603  me_inArgs.set_p(xdotdot_tangent_index_, inArgs.get_x_dot_dot());
604 
605  model_->evalModel(me_inArgs, meo);
606 
607  // transpose reduced dg/dp if necessary
608  if (dgdp_trans != Teuchos::null) {
610  Thyra::DetachedMultiVectorView<Scalar> dgdp_trans_view(*dgdp_trans);
611  for (int i = 0; i < num_param_; ++i)
612  for (int j = 0; j < num_response_; ++j)
613  dgdp_view(j, i) = dgdp_trans_view(i, j);
614  }
615 
616  // Compute (dg/dx) * (dx/dp) + (dg/dp) if the underlying ME doesn't already
617  // do this.
618  if (!use_dgdp_as_tangent_ && dxdp != Teuchos::null) {
619  MEB::DerivativeSupport dgdx_support =
620  me_outArgs.supports(MEB::OUT_ARG_DgDx, g_index_);
621  if (dgdx_support.supports(MEB::DERIV_LINEAR_OP)) {
622  my_dgdx_->apply(Thyra::NOTRANS, *dxdp, dgdp.ptr(), Scalar(1.0),
623  Scalar(1.0));
624  }
625  else if (dgdx_support.supports(MEB::DERIV_MV_GRADIENT_FORM)) {
626  my_dgdx_mv_->apply(Thyra::TRANS, *dxdp, dgdp.ptr(), Scalar(1.0),
627  Scalar(1.0));
628  }
629  else
631  true, std::logic_error,
632  "Jacobian form of dg/dx not supported for reduced sensitivity");
633  }
634  }
635 }
636 
637 template <class Scalar>
640 {
641  Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
642  pl->set<bool>("Use DfDp as Tangent", false);
643  pl->set<bool>("Use DgDp as Tangent", false);
644  pl->set<int>("Sensitivity Parameter Index", 0);
645  pl->set<int>("Response Function Index", -1);
646  pl->set<int>("Sensitivity X Tangent Index", 1);
647  pl->set<int>("Sensitivity X-Dot Tangent Index", 2);
648  pl->set<int>("Sensitivity X-Dot-Dot Tangent Index", 3);
649  return pl;
650 }
651 
652 } // namespace Tempus
653 
654 #endif
RCP< const VectorBase< Scalar > > get_x_dot() const
T & get(const std::string &name, T def_value)
ParameterList & set(std::string const &name, T const &value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
void setForwardSolutionHistory(const Teuchos::RCP< const Tempus::SolutionHistory< Scalar > > &sh)
Set solution history from forward state evaluation (for interpolation)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
virtual void setForwardSolutionState(const Teuchos::RCP< const Tempus::SolutionState< Scalar > > &s)
Set solution state from forward state evaluation (for frozen state)
Evaluation< VectorBase< Scalar > > get_g(int j) const
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int p) const
Evaluation< VectorBase< Scalar > > get_f() const
bool supports(EOutArgsMembers arg) const
static magnitudeType rmax()
Derivative< Scalar > get_DgDp(int j, int l) const
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_p_space(int p) const
Teuchos::RCP< Thyra::LinearOpBase< Scalar > > create_DgDx_dot_op(int j) const
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
virtual std::string description() const
void validateParametersAndSetDefaults(ParameterList const &validParamList, int const depth=1000)
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Teuchos::RCP< Thyra::LinearOpBase< Scalar > > create_DgDp_op(int j, int l) const
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_x_space() const
Implicit concrete LinearOpBase subclass that takes a flattended out multi-vector and performs a multi...
size_type size() const
Teuchos::RCP< const Thyra::LinearOpWithSolveFactoryBase< Scalar > > get_W_factory() const
Derivative< Scalar > get_DgDx_dot(int j) const
RCP< const VectorBase< Scalar > > get_x_dot_dot() const
Teuchos::RCP< Thyra::LinearOpBase< Scalar > > create_DgDx_op(int j) const
#define TEUCHOS_ASSERT(assertion_test)
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_f_space() const
RCP< LinearOpBase< Scalar > > get_W_op() const
RCP< const VectorBase< Scalar > > get_x() const
StaggeredForwardSensitivityModelEvaluator(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, const Teuchos::RCP< MultiVector > &dxdp_init=Teuchos::null, const Teuchos::RCP< MultiVector > &dx_dotdp_init=Teuchos::null, const Teuchos::RCP< MultiVector > &dx_dotdot_dp_init=Teuchos::null)
Constructor.
Solution state for integrators and steppers.
RCP< const VectorBase< Scalar > > get_p(int l) const
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_g_space(int j) const
void evalModelImpl(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Derivative< Scalar > get_DgDx(int j) const