10 #ifndef Tempus_AdjointSensitivityModelEvaluator_impl_hpp
11 #define Tempus_AdjointSensitivityModelEvaluator_impl_hpp
14 #include "Tempus_config.hpp"
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"
28 template <
typename Scalar>
32 adjoint_residual_model,
35 const Scalar& t_init,
const Scalar& t_final,
const bool is_pseudotransient,
38 adjoint_residual_model_(adjoint_residual_model),
39 adjoint_solve_model_(adjoint_solve_model),
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())
53 if (pList != Teuchos::null) *pl = *pList;
57 p_index_ = pl->
get<
int>(
"Sensitivity Parameter Index", 0);
58 g_index_ = pl->
get<
int>(
"Response Function Index", 0);
64 "AdjointSensitivityModelEvaluator currently does not support "
65 <<
"non-constant mass matrix df/dx_dot!");
75 MEB::InArgs<Scalar> me_inArgs =
model_->createInArgs();
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);
91 inArgs.set_Np(me_inArgs.Np());
94 MEB::OutArgs<Scalar> me_outArgs =
model_->createOutArgs();
95 MEB::OutArgs<Scalar> adj_mer_outArgs =
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);
110 if (me_inArgs.supports(MEB::IN_ARG_x_dot)) {
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));
122 template <
typename Scalar>
124 const Scalar t_final)
129 template <
typename Scalar>
134 if (is_pseudotransient_)
135 forward_state_ = sh_->getCurrentState();
138 forward_state_ = Teuchos::null;
142 mass_matrix_is_computed_ =
false;
143 jacobian_matrix_is_computed_ =
false;
144 response_gradient_is_computed_ =
false;
147 template <
typename Scalar>
152 return model_->get_p_space(p);
155 template <
typename Scalar>
160 return model_->get_p_names(p);
163 template <
typename Scalar>
167 return adjoint_space_;
170 template <
typename Scalar>
174 return residual_space_;
177 template <
typename Scalar>
182 if (j == 0)
return response_space_;
183 return model_->get_g_space(g_index_);
186 template <
typename Scalar>
191 adjoint_solve_model_->create_W_op();
192 if (adjoint_op == Teuchos::null)
return Teuchos::null;
194 return Thyra::nonconstMultiVectorLinearOp(adjoint_op, num_adjoint_);
197 template <
typename Scalar>
202 using Teuchos::rcp_dynamic_cast;
205 RCP<const LOWSFB> alowsfb = adjoint_solve_model_->get_W_factory();
206 if (alowsfb == Teuchos::null)
207 return Teuchos::null;
209 return Thyra::multiVectorLinearOpWithSolveFactory(alowsfb, residual_space_,
213 template <
typename Scalar>
217 return prototypeInArgs_;
220 template <
typename Scalar>
226 using Teuchos::rcp_dynamic_cast;
228 MEB::InArgs<Scalar> me_nominal = model_->getNominalValues();
229 MEB::InArgs<Scalar> nominal = this->createInArgs();
234 RCP<Thyra::VectorBase<Scalar> > x = Thyra::createMember(*adjoint_space_);
235 Thyra::assign(x.ptr(), zero);
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);
245 const int np = model_->Np();
246 for (
int i = 0; i < np; ++i) nominal.set_p(i, me_nominal.get_p(i));
251 template <
typename Scalar>
255 return prototypeOutArgs_;
258 template <
typename Scalar>
265 using Teuchos::rcp_dynamic_cast;
267 TEMPUS_FUNC_TIME_MONITOR_DIFF(
268 "Tempus::AdjointSensitivityModelEvaluator::evalModel()", TEMPUS_EVAL);
278 const Scalar t = inArgs.
get_t();
280 if (is_pseudotransient_)
281 forward_t = forward_state_->getTime();
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);
288 sh_->interpolateState(forward_t, forward_state_.get());
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)) {
298 me_inArgs.set_x_dot(forward_state_->getXDot());
300 if (is_pseudotransient_) {
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));
309 me_inArgs.set_x_dot(my_x_dot_);
314 me_inArgs.set_x_dot(Teuchos::null);
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));
326 RCP<Thyra::LinearOpBase<Scalar> > op;
328 if (op != Teuchos::null) {
329 TEMPUS_FUNC_TIME_MONITOR_DIFF(
330 "Tempus::AdjointSensitivityModelEvaluator::evalModel::W",
332 if (me_inArgs.supports(MEB::IN_ARG_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());
338 me_inArgs.set_beta(-inArgs.
get_beta());
341 RCP<Thyra::MultiVectorLinearOp<Scalar> > mv_adjoint_op =
343 RCP<Thyra::LinearOpBase<Scalar> > adjoint_op =
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);
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();
358 rcp_dynamic_cast<
const DMVPV>(adjoint_x,
true)->getMultiVector();
366 if (adjoint_f != Teuchos::null) {
367 TEMPUS_FUNC_TIME_MONITOR_DIFF(
368 "Tempus::AdjointSensitivityModelEvaluator::evalModel::f",
371 RCP<Thyra::MultiVectorBase<Scalar> > adjoint_f_mv =
372 rcp_dynamic_cast<
DMVPV>(adjoint_f,
true)->getNonconstMultiVector();
374 MEB::OutArgs<Scalar> me_outArgs = model_->createOutArgs();
375 MEB::OutArgs<Scalar> adj_me_outArgs =
376 adjoint_residual_model_->createOutArgs();
381 TEMPUS_FUNC_TIME_MONITOR_DIFF(
382 "Tempus::AdjointSensitivityModelEvaluator::evalModel::dg/dx",
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_) {
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;
395 Thyra::assign(adjoint_f_mv.ptr(), *my_dgdx_mv_);
401 TEMPUS_FUNC_TIME_MONITOR_DIFF(
402 "Tempus::AdjointSensitivityModelEvaluator::evalModel::df/dx",
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;
413 my_dfdx_->apply(
Thyra::NOTRANS, *adjoint_x_mv, adjoint_f_mv.ptr(),
414 Scalar(-1.0), Scalar(1.0));
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)
430 if (mass_matrix_is_identity_) {
432 Thyra::V_StVpV(adjoint_f_mv.ptr(), Scalar(-1.0), *adjoint_f_mv,
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;
447 adjoint_f_mv.ptr(), Scalar(1.0), Scalar(-1.0));
457 if (adjoint_g != Teuchos::null) {
458 TEMPUS_FUNC_TIME_MONITOR_DIFF(
459 "Tempus::AdjointSensitivityModelEvaluator::evalModel::g",
461 RCP<Thyra::MultiVectorBase<Scalar> > adjoint_g_mv =
462 rcp_dynamic_cast<
DMVPV>(adjoint_g,
true)->getNonconstMultiVector();
464 MEB::OutArgs<Scalar> me_outArgs = model_->createOutArgs();
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)) {
472 MEB::Derivative<Scalar>(adjoint_g_mv, MEB::DERIV_MV_GRADIENT_FORM));
473 model_->evalModel(me_inArgs, me_outArgs);
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);
482 MEB::Derivative<Scalar>(dgdp_trans, MEB::DERIV_MV_JACOBIAN_FORM));
483 model_->evalModel(me_inArgs, me_outArgs);
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);
491 "Invalid dg/dp support");
492 me_outArgs.set_DgDp(g_index_, p_index_, MEB::Derivative<Scalar>());
495 MEB::DerivativeSupport dfdp_support =
496 me_outArgs.supports(MEB::OUT_ARG_DfDp, p_index_);
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_));
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());
510 MEB::Derivative<Scalar>(my_dfdp_mv_, MEB::DERIV_MV_JACOBIAN_FORM));
511 my_dfdp_op_ = my_dfdp_mv_;
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());
520 MEB::Derivative<Scalar>(my_dfdp_mv_, MEB::DERIV_MV_GRADIENT_FORM));
521 my_dfdp_op_ = my_dfdp_mv_;
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),
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);
539 template <
class Scalar>
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);
static Teuchos::RCP< const Teuchos::ParameterList > getValidParameters()
Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > adjoint_solve_model_
Teuchos::RCP< const DMVPVS > residual_space_
Teuchos::RCP< const DMVPVS > adjoint_space_
RCP< const VectorBase< Scalar > > get_x_dot() const
bool mass_matrix_is_identity_
T & get(const std::string &name, T def_value)
Thyra::ModelEvaluatorBase::InArgs< Scalar > createInArgs() const
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
Teuchos::RCP< const Thyra::ModelEvaluator< Scalar > > adjoint_residual_model_
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_x_space() const
RCP< LinearOpBase< Scalar > > getNonconstLinearOp()
Evaluation< VectorBase< Scalar > > get_g(int j) const
Evaluation< VectorBase< Scalar > > get_f() const
Teuchos::RCP< const Teuchos::Array< std::string > > get_p_names(int p) const
Teuchos::RCP< const DMVPVS > response_space_
ParameterList & set(std::string const &name, T &&value, std::string const &docString="", RCP< const ParameterEntryValidator > const &validator=null)
bool supports(EOutArgsMembers arg) const
Teuchos::RCP< const Thyra::LinearOpWithSolveFactoryBase< Scalar > > get_W_factory() const
static magnitudeType rmax()
bool mass_matrix_is_constant_
void evalModelImpl(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Thyra::ModelEvaluatorBase::OutArgs< Scalar > createOutArgsImpl() const
Thyra::ModelEvaluatorBase::InArgs< Scalar > prototypeInArgs_
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Thyra::ModelEvaluatorBase::OutArgs< Scalar > prototypeOutArgs_
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_f_space() const
virtual std::string description() const
void validateParametersAndSetDefaults(ParameterList const &validParamList, int const depth=1000)
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_g_space(int j) const
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
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.
void setForwardSolutionHistory(const Teuchos::RCP< const Tempus::SolutionHistory< Scalar > > &sh)
Set solution history from forward evaluation.
Implicit concrete LinearOpBase subclass that takes a flattended out multi-vector and performs a multi...
Thyra::ModelEvaluatorBase::InArgs< Scalar > getNominalValues() const
void setFinalTime(const Scalar t_final)
Set the final time from the forward evaluation.
#define TEUCHOS_ASSERT(assertion_test)
Teuchos::RCP< Thyra::LinearOpBase< Scalar > > create_W_op() const
RCP< LinearOpBase< Scalar > > get_W_op() const
RCP< const VectorBase< Scalar > > get_x() const
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_p_space(int p) const
RCP< const VectorBase< Scalar > > get_p(int l) const