Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_AdjointSensitivityModelEvaluator_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_AdjointSensitivityModelEvaluator_impl_hpp
11#define Tempus_AdjointSensitivityModelEvaluator_impl_hpp
12
13#include "Teuchos_TimeMonitor.hpp"
14#include "Tempus_config.hpp"
15
16#include "Thyra_DefaultMultiVectorProductVectorSpace.hpp"
17#include "Thyra_DefaultMultiVectorProductVector.hpp"
20#include "Thyra_DefaultScaledAdjointLinearOp.hpp"
21#include "Thyra_DefaultAdjointLinearOpWithSolve.hpp"
23#include "Thyra_VectorStdOps.hpp"
24#include "Thyra_MultiVectorStdOps.hpp"
25
26namespace Tempus {
27
28template <typename Scalar>
30 const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >& model,
31 const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >&
32 adjoint_residual_model,
33 const Teuchos::RCP<const Thyra::ModelEvaluator<Scalar> >&
34 adjoint_solve_model,
35 const Scalar& t_init, const Scalar& t_final, const bool is_pseudotransient,
36 const Teuchos::RCP<const Teuchos::ParameterList>& pList)
37 : model_(model),
38 adjoint_residual_model_(adjoint_residual_model),
39 adjoint_solve_model_(adjoint_solve_model),
40 t_init_(t_init),
41 t_final_(t_final),
42 is_pseudotransient_(is_pseudotransient),
43 mass_matrix_is_computed_(false),
44 jacobian_matrix_is_computed_(false),
45 response_gradient_is_computed_(false),
46 t_interp_(Teuchos::ScalarTraits<Scalar>::rmax())
47{
48 typedef Thyra::ModelEvaluatorBase MEB;
49
50 // Set parameters
51 Teuchos::RCP<Teuchos::ParameterList> pl =
52 Teuchos::rcp(new Teuchos::ParameterList);
53 if (pList != Teuchos::null) *pl = *pList;
54 pl->validateParametersAndSetDefaults(*this->getValidParameters());
55 mass_matrix_is_constant_ = pl->get<bool>("Mass Matrix Is Constant");
56 mass_matrix_is_identity_ = pl->get<bool>("Mass Matrix Is Identity");
57 p_index_ = pl->get<int>("Sensitivity Parameter Index", 0);
58 g_index_ = pl->get<int>("Response Function Index", 0);
59 num_adjoint_ = model_->get_g_space(g_index_)->dim();
60
61 // We currently do not support a non-constant mass matrix
62 TEUCHOS_TEST_FOR_EXCEPTION(
63 mass_matrix_is_constant_ == false, std::logic_error,
64 "AdjointSensitivityModelEvaluator currently does not support "
65 << "non-constant mass matrix df/dx_dot!");
66
68 Thyra::multiVectorProductVectorSpace(model_->get_f_space(), num_adjoint_);
70 Thyra::multiVectorProductVectorSpace(model_->get_x_space(), num_adjoint_);
71 response_space_ = Thyra::multiVectorProductVectorSpace(
72 model_->get_p_space(p_index_), num_adjoint_);
73
74 // forward and adjoint models must support same InArgs
75 MEB::InArgs<Scalar> me_inArgs = model_->createInArgs();
76 me_inArgs.assertSameSupport(adjoint_residual_model_->createInArgs());
77 me_inArgs.assertSameSupport(adjoint_solve_model_->createInArgs());
78
79 MEB::InArgsSetup<Scalar> inArgs;
80 inArgs.setModelEvalDescription(this->description());
81 inArgs.setSupports(MEB::IN_ARG_x);
82 inArgs.setSupports(MEB::IN_ARG_t);
83 if (me_inArgs.supports(MEB::IN_ARG_x_dot))
84 inArgs.setSupports(MEB::IN_ARG_x_dot);
85 if (me_inArgs.supports(MEB::IN_ARG_alpha))
86 inArgs.setSupports(MEB::IN_ARG_alpha);
87 if (me_inArgs.supports(MEB::IN_ARG_beta))
88 inArgs.setSupports(MEB::IN_ARG_beta);
89
90 // Support additional parameters for x and xdot
91 inArgs.set_Np(me_inArgs.Np());
92 prototypeInArgs_ = inArgs;
93
94 MEB::OutArgs<Scalar> me_outArgs = model_->createOutArgs();
95 MEB::OutArgs<Scalar> adj_mer_outArgs =
96 adjoint_residual_model_->createOutArgs();
97 MEB::OutArgs<Scalar> adj_mes_outArgs = adjoint_solve_model_->createOutArgs();
98 MEB::OutArgsSetup<Scalar> outArgs;
99 outArgs.setModelEvalDescription(this->description());
100 outArgs.set_Np_Ng(me_inArgs.Np(), 2);
101 outArgs.setSupports(MEB::OUT_ARG_f);
102 if (adj_mes_outArgs.supports(MEB::OUT_ARG_W_op))
103 outArgs.setSupports(MEB::OUT_ARG_W_op);
104 prototypeOutArgs_ = outArgs;
105
106 // Adjoint residual ME must support W_op to define adjoint ODE/DAE.
107 // Must support alpha, beta if it suports x_dot
108 TEUCHOS_ASSERT(me_inArgs.supports(MEB::IN_ARG_x));
109 TEUCHOS_ASSERT(adj_mer_outArgs.supports(MEB::OUT_ARG_W_op));
110 if (me_inArgs.supports(MEB::IN_ARG_x_dot)) {
111 TEUCHOS_ASSERT(me_inArgs.supports(MEB::IN_ARG_alpha));
112 TEUCHOS_ASSERT(me_inArgs.supports(MEB::IN_ARG_beta));
113 }
114 MEB::DerivativeSupport dgdx_support =
115 me_outArgs.supports(MEB::OUT_ARG_DgDx, g_index_);
116 MEB::DerivativeSupport dgdp_support =
117 me_outArgs.supports(MEB::OUT_ARG_DgDp, g_index_, p_index_);
118 TEUCHOS_ASSERT(dgdx_support.supports(MEB::DERIV_MV_GRADIENT_FORM));
119 TEUCHOS_ASSERT(dgdp_support.supports(MEB::DERIV_MV_GRADIENT_FORM));
120}
121
122template <typename Scalar>
124 const Scalar t_final)
125{
126 t_final_ = t_final;
127}
128
129template <typename Scalar>
131 const Teuchos::RCP<const Tempus::SolutionHistory<Scalar> >& sh)
132{
133 sh_ = sh;
134 if (is_pseudotransient_)
135 forward_state_ = sh_->getCurrentState();
136 else {
137 t_interp_ = Teuchos::ScalarTraits<Scalar>::rmax();
138 forward_state_ = Teuchos::null;
139 }
140
141 // Reset computation flags because we have done a new forward integration
142 mass_matrix_is_computed_ = false;
143 jacobian_matrix_is_computed_ = false;
144 response_gradient_is_computed_ = false;
145}
146
147template <typename Scalar>
148Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> >
150{
151 TEUCHOS_ASSERT(p < model_->Np());
152 return model_->get_p_space(p);
153}
154
155template <typename Scalar>
156Teuchos::RCP<const Teuchos::Array<std::string> >
158{
159 TEUCHOS_ASSERT(p < model_->Np());
160 return model_->get_p_names(p);
161}
162
163template <typename Scalar>
164Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> >
166{
167 return adjoint_space_;
168}
169
170template <typename Scalar>
171Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> >
173{
174 return residual_space_;
175}
176
177template <typename Scalar>
178Teuchos::RCP<const Thyra::VectorSpaceBase<Scalar> >
180{
181 TEUCHOS_ASSERT(j == 0 || j == 1);
182 if (j == 0) return response_space_;
183 return model_->get_g_space(g_index_);
184}
185
186template <typename Scalar>
187Teuchos::RCP<Thyra::LinearOpBase<Scalar> >
189{
190 Teuchos::RCP<Thyra::LinearOpBase<Scalar> > adjoint_op =
191 adjoint_solve_model_->create_W_op();
192 if (adjoint_op == Teuchos::null) return Teuchos::null;
193
194 return Thyra::nonconstMultiVectorLinearOp(adjoint_op, num_adjoint_);
195}
196
197template <typename Scalar>
198Teuchos::RCP<const Thyra::LinearOpWithSolveFactoryBase<Scalar> >
200{
201 using Teuchos::RCP;
202 using Teuchos::rcp_dynamic_cast;
204
205 RCP<const LOWSFB> alowsfb = adjoint_solve_model_->get_W_factory();
206 if (alowsfb == Teuchos::null)
207 return Teuchos::null; // adjoint_solve_model_ doesn't support W_factory
208
209 return Thyra::multiVectorLinearOpWithSolveFactory(alowsfb, residual_space_,
210 adjoint_space_);
211}
212
213template <typename Scalar>
214Thyra::ModelEvaluatorBase::InArgs<Scalar>
216{
217 return prototypeInArgs_;
218}
219
220template <typename Scalar>
221Thyra::ModelEvaluatorBase::InArgs<Scalar>
223{
224 typedef Thyra::ModelEvaluatorBase MEB;
225 using Teuchos::RCP;
226 using Teuchos::rcp_dynamic_cast;
227
228 MEB::InArgs<Scalar> me_nominal = model_->getNominalValues();
229 MEB::InArgs<Scalar> nominal = this->createInArgs();
230
231 const Scalar zero = Teuchos::ScalarTraits<Scalar>::zero();
232
233 // Set initial x, x_dot
234 RCP<Thyra::VectorBase<Scalar> > x = Thyra::createMember(*adjoint_space_);
235 Thyra::assign(x.ptr(), zero);
236 nominal.set_x(x);
237
238 if (me_nominal.supports(MEB::IN_ARG_x_dot)) {
239 RCP<Thyra::VectorBase<Scalar> > x_dot =
240 Thyra::createMember(*adjoint_space_);
241 Thyra::assign(x_dot.ptr(), zero);
242 nominal.set_x_dot(x_dot);
243 }
244
245 const int np = model_->Np();
246 for (int i = 0; i < np; ++i) nominal.set_p(i, me_nominal.get_p(i));
247
248 return nominal;
249}
250
251template <typename Scalar>
252Thyra::ModelEvaluatorBase::OutArgs<Scalar>
254{
255 return prototypeOutArgs_;
256}
257
258template <typename Scalar>
260 const Thyra::ModelEvaluatorBase::InArgs<Scalar>& inArgs,
261 const Thyra::ModelEvaluatorBase::OutArgs<Scalar>& outArgs) const
262{
263 typedef Thyra::ModelEvaluatorBase MEB;
264 using Teuchos::RCP;
265 using Teuchos::rcp_dynamic_cast;
266
267 TEMPUS_FUNC_TIME_MONITOR_DIFF(
268 "Tempus::AdjointSensitivityModelEvaluator::evalModel()", TEMPUS_EVAL);
269
270 // Note: adjoint models compute the transposed W (either explicitly or
271 // implicitly. Thus we need to always call their evalModel() functions
272 // whenever computing the adjoint operators, and subsequent calls to apply()
273 // do not transpose them.
274
275 // Interpolate forward solution at supplied time, reusing previous
276 // interpolation if possible
277 TEUCHOS_ASSERT(sh_ != Teuchos::null);
278 const Scalar t = inArgs.get_t();
279 Scalar forward_t;
280 if (is_pseudotransient_)
281 forward_t = forward_state_->getTime();
282 else {
283 forward_t = t_final_ + t_init_ - t;
284 if (forward_state_ == Teuchos::null || t_interp_ != t) {
285 if (forward_state_ == Teuchos::null)
286 forward_state_ = sh_->interpolateState(forward_t);
287 else
288 sh_->interpolateState(forward_t, forward_state_.get());
289 t_interp_ = t;
290 }
291 }
292
293 // setup input arguments for model
294 MEB::InArgs<Scalar> me_inArgs = model_->getNominalValues();
295 me_inArgs.set_x(forward_state_->getX());
296 if (me_inArgs.supports(MEB::IN_ARG_x_dot)) {
297 if (inArgs.get_x_dot() != Teuchos::null)
298 me_inArgs.set_x_dot(forward_state_->getXDot());
299 else {
300 if (is_pseudotransient_) {
301 // For pseudo-transient, we have to always use the same form of the
302 // residual in order to reuse df/dx, df/dx_dot,..., so we force
303 // the underlying ME to always compute the implicit form with x_dot == 0
304 // if it wasn't provided.
305 if (my_x_dot_ == Teuchos::null) {
306 my_x_dot_ = Thyra::createMember(model_->get_x_space());
307 Thyra::assign(my_x_dot_.ptr(), Scalar(0.0));
308 }
309 me_inArgs.set_x_dot(my_x_dot_);
310 }
311 else {
312 // clear out xdot if it was set in nominalValues to get to ensure we
313 // get the explicit form
314 me_inArgs.set_x_dot(Teuchos::null);
315 }
316 }
317 }
318 if (me_inArgs.supports(MEB::IN_ARG_t)) me_inArgs.set_t(forward_t);
319 const int np = me_inArgs.Np();
320 for (int i = 0; i < np; ++i) me_inArgs.set_p(i, inArgs.get_p(i));
321
322 // compute adjoint W == model W
323 // It would be nice to not reevaluate W in the psuedo-transient case, but
324 // it isn't clear how to do this in a clean way. Probably just need to
325 // control that with the nonlinear solver.
326 RCP<Thyra::LinearOpBase<Scalar> > op;
327 if (outArgs.supports(MEB::OUT_ARG_W_op)) op = outArgs.get_W_op();
328 if (op != Teuchos::null) {
329 TEMPUS_FUNC_TIME_MONITOR_DIFF(
330 "Tempus::AdjointSensitivityModelEvaluator::evalModel::W",
331 TEMPUS_EVAL_W);
332 if (me_inArgs.supports(MEB::IN_ARG_alpha))
333 me_inArgs.set_alpha(inArgs.get_alpha());
334 if (me_inArgs.supports(MEB::IN_ARG_beta)) {
335 if (me_inArgs.supports(MEB::IN_ARG_x_dot))
336 me_inArgs.set_beta(inArgs.get_beta()); // Implicit form (see below)
337 else
338 me_inArgs.set_beta(-inArgs.get_beta()); // Explicit form (see below)
339 }
340
341 RCP<Thyra::MultiVectorLinearOp<Scalar> > mv_adjoint_op =
342 rcp_dynamic_cast<Thyra::MultiVectorLinearOp<Scalar> >(op, true);
343 RCP<Thyra::LinearOpBase<Scalar> > adjoint_op =
344 mv_adjoint_op->getNonconstLinearOp();
345 MEB::OutArgs<Scalar> adj_me_outArgs = adjoint_solve_model_->createOutArgs();
346 adj_me_outArgs.set_W_op(adjoint_op);
347 adjoint_solve_model_->evalModel(me_inArgs, adj_me_outArgs);
348 }
349
350 RCP<Thyra::VectorBase<Scalar> > adjoint_f = outArgs.get_f();
351 RCP<Thyra::VectorBase<Scalar> > adjoint_g = outArgs.get_g(0);
352 RCP<Thyra::VectorBase<Scalar> > g = outArgs.get_g(1);
353 RCP<const Thyra::MultiVectorBase<Scalar> > adjoint_x_mv;
354 RCP<const Thyra::VectorBase<Scalar> > adjoint_x;
355 if (adjoint_f != Teuchos::null || adjoint_g != Teuchos::null) {
356 adjoint_x = inArgs.get_x().assert_not_null();
357 adjoint_x_mv =
358 rcp_dynamic_cast<const DMVPV>(adjoint_x, true)->getMultiVector();
359 }
360
361 // Compute adjoint residual F(y):
362 // * For implicit form, F(y) = d/dt( df/dx_dot^T*y ) + df/dx^T*y - dg/dx^T
363 // * For explict form, F(y) = -df/dx^T*y + dg/dx^T
364 // For implicit form, we assume df/dx_dot is constant w.r.t. x, x_dot, and t,
365 // so the residual becomes F(y) = df/dx_dot^T*y_dot + df/dx^T*y - dg/dx^T
366 if (adjoint_f != Teuchos::null) {
367 TEMPUS_FUNC_TIME_MONITOR_DIFF(
368 "Tempus::AdjointSensitivityModelEvaluator::evalModel::f",
369 TEMPUS_EVAL_F);
370
371 RCP<Thyra::MultiVectorBase<Scalar> > adjoint_f_mv =
372 rcp_dynamic_cast<DMVPV>(adjoint_f, true)->getNonconstMultiVector();
373
374 MEB::OutArgs<Scalar> me_outArgs = model_->createOutArgs();
375 MEB::OutArgs<Scalar> adj_me_outArgs =
376 adjoint_residual_model_->createOutArgs();
377
378 // dg/dx^T
379 // Don't re-evaluate dg/dx for pseudotransient
380 {
381 TEMPUS_FUNC_TIME_MONITOR_DIFF(
382 "Tempus::AdjointSensitivityModelEvaluator::evalModel::dg/dx",
383 TEMPUS_EVAL_DGDX);
384 if (my_dgdx_mv_ == Teuchos::null)
385 my_dgdx_mv_ = Thyra::createMembers(
386 model_->get_x_space(), model_->get_g_space(g_index_)->dim());
387 if (!response_gradient_is_computed_) {
388 me_outArgs.set_DgDx(
389 g_index_,
390 MEB::Derivative<Scalar>(my_dgdx_mv_, MEB::DERIV_MV_GRADIENT_FORM));
391 model_->evalModel(me_inArgs, me_outArgs);
392 me_outArgs.set_DgDx(g_index_, MEB::Derivative<Scalar>());
393 if (is_pseudotransient_) response_gradient_is_computed_ = true;
394 }
395 Thyra::assign(adjoint_f_mv.ptr(), *my_dgdx_mv_);
396 }
397
398 // Explicit form of the residual F(y) = -df/dx^T*y + dg/dx^T
399 // Don't re-evaluate df/dx for pseudotransient
400 {
401 TEMPUS_FUNC_TIME_MONITOR_DIFF(
402 "Tempus::AdjointSensitivityModelEvaluator::evalModel::df/dx",
403 TEMPUS_EVAL_DFDX);
404 if (my_dfdx_ == Teuchos::null)
405 my_dfdx_ = adjoint_residual_model_->create_W_op();
406 if (!jacobian_matrix_is_computed_) {
407 adj_me_outArgs.set_W_op(my_dfdx_);
408 if (me_inArgs.supports(MEB::IN_ARG_alpha)) me_inArgs.set_alpha(0.0);
409 if (me_inArgs.supports(MEB::IN_ARG_beta)) me_inArgs.set_beta(1.0);
410 adjoint_residual_model_->evalModel(me_inArgs, adj_me_outArgs);
411 if (is_pseudotransient_) jacobian_matrix_is_computed_ = true;
412 }
413 my_dfdx_->apply(Thyra::NOTRANS, *adjoint_x_mv, adjoint_f_mv.ptr(),
414 Scalar(-1.0), Scalar(1.0));
415 }
416
417 // Implicit form residual F(y) df/dx_dot^T*y_dot + df/dx^T*y - dg/dx^T
418 // using the second scalar argument to apply() to change the explicit term
419 // above.
420 // Don't re-evaluate df/dx_dot for pseudotransient
421 if (me_inArgs.supports(MEB::IN_ARG_x_dot)) {
422 TEMPUS_FUNC_TIME_MONITOR_DIFF(
423 "Tempus::AdjointSensitivityModelEvaluator::evalModel::df/dx_dot",
424 TEMPUS_EVAL_DFDXDOT);
425 RCP<const Thyra::VectorBase<Scalar> > adjoint_x_dot = inArgs.get_x_dot();
426 if (adjoint_x_dot != Teuchos::null) {
427 RCP<const Thyra::MultiVectorBase<Scalar> > adjoint_x_dot_mv =
428 rcp_dynamic_cast<const DMVPV>(adjoint_x_dot, true)
429 ->getMultiVector();
430 if (mass_matrix_is_identity_) {
431 // F = -F + y_dot
432 Thyra::V_StVpV(adjoint_f_mv.ptr(), Scalar(-1.0), *adjoint_f_mv,
433 *adjoint_x_dot_mv);
434 }
435 else {
436 if (my_dfdxdot_ == Teuchos::null)
437 my_dfdxdot_ = adjoint_residual_model_->create_W_op();
438 if (!mass_matrix_is_computed_) {
439 adj_me_outArgs.set_W_op(my_dfdxdot_);
440 me_inArgs.set_alpha(1.0);
441 me_inArgs.set_beta(0.0);
442 adjoint_residual_model_->evalModel(me_inArgs, adj_me_outArgs);
443 if (is_pseudotransient_ || mass_matrix_is_constant_)
444 mass_matrix_is_computed_ = true;
445 }
446 my_dfdxdot_->apply(Thyra::NOTRANS, *adjoint_x_dot_mv,
447 adjoint_f_mv.ptr(), Scalar(1.0), Scalar(-1.0));
448 }
449 }
450 }
451 }
452
453 // Compute g = dg/dp^T - df/dp^T*y for computing the model parameter term in
454 // the adjoint sensitivity formula.
455 // We don't add pseudotransient logic here because this part is only
456 // evaluated once in that case anyway.
457 if (adjoint_g != Teuchos::null) {
458 TEMPUS_FUNC_TIME_MONITOR_DIFF(
459 "Tempus::AdjointSensitivityModelEvaluator::evalModel::g",
460 TEMPUS_EVAL_G);
461 RCP<Thyra::MultiVectorBase<Scalar> > adjoint_g_mv =
462 rcp_dynamic_cast<DMVPV>(adjoint_g, true)->getNonconstMultiVector();
463
464 MEB::OutArgs<Scalar> me_outArgs = model_->createOutArgs();
465
466 // dg/dp
467 MEB::DerivativeSupport dgdp_support =
468 me_outArgs.supports(MEB::OUT_ARG_DgDp, g_index_, p_index_);
469 if (dgdp_support.supports(MEB::DERIV_MV_GRADIENT_FORM)) {
470 me_outArgs.set_DgDp(
471 g_index_, p_index_,
472 MEB::Derivative<Scalar>(adjoint_g_mv, MEB::DERIV_MV_GRADIENT_FORM));
473 model_->evalModel(me_inArgs, me_outArgs);
474 }
475 else if (dgdp_support.supports(MEB::DERIV_MV_JACOBIAN_FORM)) {
476 const int num_g = model_->get_g_space(g_index_)->dim();
477 const int num_p = model_->get_p_space(p_index_)->dim();
478 RCP<Thyra::MultiVectorBase<Scalar> > dgdp_trans =
479 createMembers(model_->get_g_space(g_index_), num_p);
480 me_outArgs.set_DgDp(
481 g_index_, p_index_,
482 MEB::Derivative<Scalar>(dgdp_trans, MEB::DERIV_MV_JACOBIAN_FORM));
483 model_->evalModel(me_inArgs, me_outArgs);
484 Thyra::DetachedMultiVectorView<Scalar> dgdp_view(*adjoint_g_mv);
485 Thyra::DetachedMultiVectorView<Scalar> dgdp_trans_view(*dgdp_trans);
486 for (int i = 0; i < num_p; ++i)
487 for (int j = 0; j < num_g; ++j) dgdp_view(i, j) = dgdp_trans_view(j, i);
488 }
489 else
490 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
491 "Invalid dg/dp support");
492 me_outArgs.set_DgDp(g_index_, p_index_, MEB::Derivative<Scalar>());
493
494 // dg/dp - df/dp^T*y
495 MEB::DerivativeSupport dfdp_support =
496 me_outArgs.supports(MEB::OUT_ARG_DfDp, p_index_);
497 Thyra::EOpTransp trans = Thyra::CONJTRANS;
498 if (dfdp_support.supports(MEB::DERIV_LINEAR_OP)) {
499 if (my_dfdp_op_ == Teuchos::null)
500 my_dfdp_op_ = model_->create_DfDp_op(p_index_);
501 me_outArgs.set_DfDp(p_index_, MEB::Derivative<Scalar>(my_dfdp_op_));
502 trans = Thyra::CONJTRANS;
503 }
504 else if (dfdp_support.supports(MEB::DERIV_MV_JACOBIAN_FORM)) {
505 if (my_dfdp_mv_ == Teuchos::null)
506 my_dfdp_mv_ = createMembers(model_->get_f_space(),
507 model_->get_p_space(p_index_)->dim());
508 me_outArgs.set_DfDp(
509 p_index_,
510 MEB::Derivative<Scalar>(my_dfdp_mv_, MEB::DERIV_MV_JACOBIAN_FORM));
511 my_dfdp_op_ = my_dfdp_mv_;
512 trans = Thyra::CONJTRANS;
513 }
514 else if (dfdp_support.supports(MEB::DERIV_MV_GRADIENT_FORM)) {
515 if (my_dfdp_mv_ == Teuchos::null)
516 my_dfdp_mv_ = createMembers(model_->get_p_space(p_index_),
517 model_->get_f_space()->dim());
518 me_outArgs.set_DfDp(
519 p_index_,
520 MEB::Derivative<Scalar>(my_dfdp_mv_, MEB::DERIV_MV_GRADIENT_FORM));
521 my_dfdp_op_ = my_dfdp_mv_;
522 trans = Thyra::CONJ;
523 }
524 else
525 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
526 "Invalid df/dp support");
527 model_->evalModel(me_inArgs, me_outArgs);
528 my_dfdp_op_->apply(trans, *adjoint_x_mv, adjoint_g_mv.ptr(), Scalar(-1.0),
529 Scalar(1.0));
530 }
531
532 if (g != Teuchos::null) {
533 MEB::OutArgs<Scalar> me_outArgs = model_->createOutArgs();
534 me_outArgs.set_g(g_index_, g);
535 model_->evalModel(me_inArgs, me_outArgs);
536 }
537}
538
539template <class Scalar>
540Teuchos::RCP<const Teuchos::ParameterList>
542{
543 Teuchos::RCP<Teuchos::ParameterList> pl = Teuchos::parameterList();
544 pl->set<int>("Sensitivity Parameter Index", 0);
545 pl->set<int>("Response Function Index", 0);
546 pl->set<bool>("Mass Matrix Is Constant", true);
547 pl->set<bool>("Mass Matrix Is Identity", false);
548 return pl;
549}
550
551} // namespace Tempus
552
553#endif
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_p_space(int p) const
Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > adjoint_residual_model_
Thyra::ModelEvaluatorBase::InArgs< Scalar > getNominalValues() const
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_x_space() const
void setForwardSolutionHistory(const Teuchos::RCP< const Tempus::SolutionHistory< Scalar > > &sh)
Set solution history from forward evaluation.
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_f_space() const
Thyra::ModelEvaluatorBase::OutArgs< Scalar > createOutArgsImpl() const
static Teuchos::RCP< const Teuchos::ParameterList > getValidParameters()
Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > model_
AdjointSensitivityModelEvaluator(const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &model, const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &adjoint_residual_model, const Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > &adjoint_solve_model, const Scalar &t_init, const Scalar &t_final, const bool is_pseudotransient, const Teuchos::RCP< const Teuchos::ParameterList > &pList=Teuchos::null)
Constructor.
Teuchos::RCP< const Thyra::LinearOpWithSolveFactoryBase< Scalar > > get_W_factory() const
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int p) const
Thyra::ModelEvaluatorBase::InArgs< Scalar > createInArgs() const
void evalModelImpl(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > adjoint_solve_model_
void setFinalTime(const Scalar t_final)
Set the final time from the forward evaluation.
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_g_space(int j) const
Teuchos::RCP< Thyra::LinearOpBase< Scalar > > create_W_op() const
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...