10#ifndef THYRA_SIMPLE_2D_MODEL_EVALUATOR_DEF_HPP
11#define THYRA_SIMPLE_2D_MODEL_EVALUATOR_DEF_HPP
14#include "Thyra_Simple2DModelEvaluator_decl.hpp"
15#include "Thyra_SimpleDenseLinearOp.hpp"
16#include "Thyra_DefaultSpmdVectorSpace.hpp"
17#include "Thyra_DefaultSerialDenseLinearOpWithSolveFactory.hpp"
18#include "Thyra_DefaultPreconditioner.hpp"
19#include "Thyra_DetachedMultiVectorView.hpp"
20#include "Thyra_DetachedVectorView.hpp"
21#include "Thyra_MultiVectorStdOps.hpp"
22#include "Thyra_VectorStdOps.hpp"
33simple2DModelEvaluator()
73 showGetInvalidArg_ = showGetInvalidArg;
100 return nominalValues_;
104template<
class Scalar>
108 return createNonconstSimpleDenseLinearOp<Scalar>(
109 createMembers<Scalar>(f_space_, x_space_->dim())
114template<
class Scalar>
118 return nonconstUnspecifiedPrec<Scalar>(
119 createNonconstSimpleDenseLinearOp<Scalar>(
120 createMembers<Scalar>(f_space_, x_space_->dim())
126template<
class Scalar>
134template<
class Scalar>
138 return prototypeInArgs_;
145template<
class Scalar>
149 return prototypeOutArgs_;
153template<
class Scalar>
154void Simple2DModelEvaluator<Scalar>::evalModelImpl(
159 using Teuchos::rcp_dynamic_cast;
160 const Scalar one = 1.0, two = 2.0, zero = 0.0;
162 const ConstDetachedVectorView<Scalar> x(inArgs.
get_x());
164 const RCP<Thyra::VectorBase<Scalar> > f_out = outArgs.
get_f();
165 const RCP<Thyra::LinearOpBase< Scalar > > W_op_out = outArgs.
get_W_op();
166 const RCP<Thyra::PreconditionerBase< Scalar > > W_prec_out = outArgs.
get_W_prec();
168 if (nonnull(f_out)) {
169 const DetachedVectorView<Scalar> f(f_out);
170 f[0] = x[0] + x[1] * x[1] - p_[0];
171 f[1] = d_ * (x[0] * x[0] - x[1] - p_[1]);
175 const RCP<SimpleDenseLinearOp<Scalar> > W =
176 rcp_dynamic_cast<SimpleDenseLinearOp<Scalar> >(W_op_out,
true);
177 const RCP<MultiVectorBase<Scalar> > W_mv = W->getNonconstMultiVector();
180 W_dmvv(0, 1) = two * x[1];
181 W_dmvv(1, 0) = d_ * two * x[0];
186 const RCP<SimpleDenseLinearOp<Scalar> > W_prec_op =
187 rcp_dynamic_cast<SimpleDenseLinearOp<Scalar> >(
188 W_prec_out->getNonconstUnspecifiedPrecOp(),
true);
189 const RCP<MultiVectorBase<Scalar> > W_prec_mv = W_prec_op->getNonconstMultiVector();
192 W_prec_dmvv(0, 0) = one;
193 W_prec_dmvv(0, 1) = zero;
194 W_prec_dmvv(1, 0) = zero;
195 W_prec_dmvv(1, 1) = -one/d_;
204template<
class Scalar>
205Simple2DModelEvaluator<Scalar>::Simple2DModelEvaluator()
206 : x_space_(Thyra::defaultSpmdVectorSpace<Scalar>(2)),
208 W_factory_(Thyra::defaultSerialDenseLinearOpWithSolveFactory<Scalar>()),
210 p_(x_space_->dim(),
Teuchos::ScalarTraits<Scalar>::zero()),
211 showGetInvalidArg_(false)
216 using Thyra::createMember;
220 MEB::InArgsSetup<Scalar> inArgs;
221 inArgs.setModelEvalDescription(this->
description());
222 inArgs.setSupports(MEB::IN_ARG_x);
223 prototypeInArgs_ = inArgs;
225 MEB::OutArgsSetup<Scalar> outArgs;
226 outArgs.setModelEvalDescription(this->
description());
227 outArgs.setSupports(MEB::OUT_ARG_f);
228 outArgs.setSupports(MEB::OUT_ARG_W_op);
229 outArgs.setSupports(MEB::OUT_ARG_W_prec);
230 prototypeOutArgs_ = outArgs;
232 nominalValues_ = inArgs;
233 x0_ = createMember(x_space_);
234 V_S(x0_.
ptr(), ST::zero());
235 nominalValues_.set_x(x0_);
253#define SIMPLE_2D_MODEL_EVALUATOR_INSTANT(SCALAR) \
255 template class Simple2DModelEvaluator<SCALAR >; \
257 template Teuchos::RCP<Simple2DModelEvaluator<SCALAR > > \
258 simple2DModelEvaluator(); \
const ArrayRCP< Scalar > values() const
virtual std::string description() const
Create an explicit mutable (non-const) view of a MultiVectorBase object.
Create an explicit mutable (non-const) view of a VectorBase object.
const RTOpPack::SubVectorView< Scalar > & sv() const
Returns the explicit view as an RTOpPack::ConstSubVectorView<Scalar> object.
Concrete aggregate class for all input arguments computable by a ModelEvaluator subclass object.
RCP< const VectorBase< Scalar > > get_x() const
Precondition: supports(IN_ARG_x)==true.
Concrete aggregate class for all output arguments computable by a ModelEvaluator subclass object.
Evaluation< VectorBase< Scalar > > get_f() const
Precondition: supports(OUT_ARG_f)==true.
RCP< LinearOpBase< Scalar > > get_W_op() const
Precondition: supports(OUT_ARG_W_op)==true.
RCP< PreconditionerBase< Scalar > > get_W_prec() const
Precondition: supports(OUT_ARG_W_op)==true.
Base subclass for ModelEvaluator that defines some basic types.
Simple 2d simulation only ModelEvaluator for f(x) = 0.
void set_p(const Teuchos::ArrayView< const Scalar > &p)
void set_d(const Scalar &d)
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_x_space() const
Teuchos::RCP< const Thyra::LinearOpWithSolveFactoryBase< Scalar > > get_W_factory() const
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_f_space() const
Teuchos::RCP< Thyra::LinearOpBase< Scalar > > create_W_op() const
Thyra::ModelEvaluatorBase::InArgs< Scalar > createInArgs() const
Thyra::ModelEvaluatorBase::InArgs< Scalar > getNominalValues() const
void set_x0(const Teuchos::ArrayView< const Scalar > &x0)
void setShowGetInvalidArgs(bool showGetInvalidArg)
Teuchos::RCP< Thyra::PreconditionerBase< Scalar > > create_W_prec() const
Abstract interface for finite-dimensional dense vectors.
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)
bool nonnull(const std::shared_ptr< T > &p)
T_To & dyn_cast(T_From &from)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)