10 #ifndef Tempus_WrapperModelEvaluatorBasic_impl_hpp
11 #define Tempus_WrapperModelEvaluatorBasic_impl_hpp
15 template <
typename Scalar>
21 MEB::InArgsSetup<Scalar> inArgs(appModel_->getNominalValues());
24 if (y_ != Teuchos::null) inArgs.set_p(index_, y_);
25 if (inArgs.supports(MEB::IN_ARG_x_dot)) inArgs.set_x_dot(xDot_);
26 if (inArgs.supports(MEB::IN_ARG_t)) inArgs.set_t(time_);
27 if (inArgs.supports(MEB::IN_ARG_step_size))
28 inArgs.set_step_size(p_->timeStepSize_);
29 if (inArgs.supports(MEB::IN_ARG_alpha)) inArgs.set_alpha(p_->alpha_);
30 if (inArgs.supports(MEB::IN_ARG_beta)) inArgs.set_beta(p_->beta_);
31 if (inArgs.supports(MEB::IN_ARG_stage_number))
32 inArgs.set_stage_number(p_->stageNumber_);
34 inArgs.setModelEvalDescription(this->description());
35 return std::move(inArgs);
38 template <
typename Scalar>
43 MEB::OutArgsSetup<Scalar> outArgs(appModel_->createOutArgs());
44 outArgs.setModelEvalDescription(this->description());
45 return std::move(outArgs);
48 template <
typename Scalar>
56 MEB::InArgs<Scalar> appInArgs(inArgs);
57 MEB::OutArgs<Scalar> appOutArgs(outArgs);
60 switch (evaluationType_) {
63 appInArgs.set_x(inArgs.
get_x());
67 appOutArgs.set_f(outArgs.
get_f());
77 RCP<const Thyra::VectorBase<Scalar> > x = inArgs.
get_x();
79 RCP<Thyra::VectorBase<Scalar> > xDot =
81 appInArgs.get_x_dot());
82 timeDer_->compute(x, xDot);
83 appInArgs.set_x_dot(xDot);
91 appOutArgs.set_f(outArgs.
get_f());
92 appOutArgs.set_W_op(outArgs.
get_W_op());
93 if (outArgs.
supports(MEB::OUT_ARG_W_prec)) {
110 appInArgs.set_x_dot(inArgs.
get_x());
113 appOutArgs.set_f(outArgs.
get_f());
114 appOutArgs.set_W_op(outArgs.
get_W_op());
115 if (outArgs.
supports(MEB::OUT_ARG_W_prec)) {
127 appModel_->evalModel(appInArgs, appOutArgs);
132 #endif // Tempus_WrapperModelEvaluatorBasic_impl_hpp
Evaluate residual for the implicit ODE.
RCP< const VectorBase< Scalar > > get_x_dot() const
Evaluation< VectorBase< Scalar > > get_f() const
bool supports(EOutArgsMembers arg) const
Solve for xDot keeping x constant (for ICs).
RCP< PreconditionerBase< Scalar > > get_W_prec() const
void evalModelImpl(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Thyra::ModelEvaluatorBase::InArgs< Scalar > createInArgs() const
RCP< LinearOpBase< Scalar > > get_W_op() const
Solve for x and determine xDot from x.
RCP< const VectorBase< Scalar > > get_x() const
#define TEUCHOS_TEST_FOR_EXCEPT(throw_exception_test)
Thyra::ModelEvaluatorBase::OutArgs< Scalar > createOutArgsImpl() const