45 #ifndef SIMPLE_2D_TPETRA_MODEL_EVALUATOR_HPP
46 #define SIMPLE_2D_TPETRA_MODEL_EVALUATOR_HPP
51 #include "Tpetra_CrsMatrix.hpp"
58 template<
class Scalar>
66 using Thyra::VectorBase;
67 using Thyra::createMember;
68 typedef Thyra::ModelEvaluatorBase MEB;
77 MEB::InArgsSetup<Scalar> inArgs;
78 inArgs.setModelEvalDescription(this->description());
79 inArgs.setSupports(MEB::IN_ARG_x);
82 MEB::OutArgsSetup<Scalar> outArgs;
83 outArgs.setModelEvalDescription(this->description());
84 outArgs.setSupports(MEB::OUT_ARG_f);
85 outArgs.setSupports(MEB::OUT_ARG_W_op);
92 const RCP<const Teuchos::Comm<int> > comm =
95 const RCP<const Tpetra::Map<> > map =
rcp(
new Tpetra::Map<> (dim, 0, comm));
97 typedef Tpetra::Map<>::global_ordinal_type GO;
100 W_op_graph_->insertGlobalIndices(0, tuple<GO>(0, 1)());
101 W_op_graph_->insertGlobalIndices(1, tuple<GO>(0, 1)());
106 x0_ =
rcp(
new Tpetra::Vector<Scalar>(map));
107 x0_->putScalar(ST::zero());
113 x_space_ = Thyra::createVectorSpace<Scalar>(map);
124 set_p(Teuchos::tuple<Scalar>(2.0, 0.0)());
125 set_x0(Teuchos::tuple<Scalar>(1.0, 1.0)());
130 template<
class Scalar>
137 template<
class Scalar>
147 template<
class Scalar>
153 x0_->get1dViewNonConst()().assign(x0_in);
160 template<
class Scalar>
168 template<
class Scalar>
176 template<
class Scalar>
177 Thyra::ModelEvaluatorBase::InArgs<Scalar>
180 return nominalValues_;
184 template<
class Scalar>
190 Teuchos::rcp(
new Tpetra::CrsMatrix<Scalar,int>(W_op_graph_))
196 template<
class Scalar>
197 Thyra::ModelEvaluatorBase::InArgs<Scalar>
200 return prototypeInArgs_;
207 template<
class Scalar>
208 Thyra::ModelEvaluatorBase::OutArgs<Scalar>
211 return prototypeOutArgs_;
215 template<
class Scalar>
217 const Thyra::ModelEvaluatorBase::InArgs<Scalar> &inArgs,
218 const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs
224 using Teuchos::tuple;
225 using Teuchos::rcp_dynamic_cast;
227 typedef Tpetra::Map<>::global_ordinal_type GO;
230 const RCP<const Tpetra::Vector<Scalar, int> > x_vec =
231 ConverterT::getConstTpetraVector(inArgs.get_x());
232 const ArrayRCP<const Scalar> x = x_vec->get1dView();
234 const RCP<Tpetra::Vector<Scalar, int> > f_vec =
235 ConverterT::getTpetraVector(outArgs.get_f());
237 const RCP<Tpetra::CrsMatrix<Scalar, int> > W =
238 rcp_dynamic_cast<Tpetra::CrsMatrix<Scalar,int> >(
239 ConverterT::getTpetraOperator(outArgs.get_W_op()),
244 const ArrayRCP<Scalar>
f = f_vec->get1dViewNonConst();
245 f[0] = x[0] + x[1]*x[1] - p_[0];
246 f[1] = d_ * (x[0]*x[0] -x[1] - p_[1]);
250 W->setAllToScalar(ST::zero());
251 W->sumIntoGlobalValues(0, tuple<GO>(0, 1), tuple<Scalar>(1.0, 2.0*x[1]));
252 W->sumIntoGlobalValues(1, tuple<GO>(0, 1), tuple<Scalar>(2.0*d_*x[0], -d_));
258 #endif // SIMPLE_2D_TPETRA_MODEL_EVALUATOR_HPP