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