10 #ifndef SIMPLE_2D_TPETRA_MODEL_EVALUATOR_HPP
11 #define SIMPLE_2D_TPETRA_MODEL_EVALUATOR_HPP
16 #include "Tpetra_CrsMatrix.hpp"
23 template<
class Scalar>
31 using Thyra::VectorBase;
32 using Thyra::createMember;
33 typedef Thyra::ModelEvaluatorBase MEB;
42 MEB::InArgsSetup<Scalar> inArgs;
43 inArgs.setModelEvalDescription(this->description());
44 inArgs.setSupports(MEB::IN_ARG_x);
47 MEB::OutArgsSetup<Scalar> outArgs;
48 outArgs.setModelEvalDescription(this->description());
49 outArgs.setSupports(MEB::OUT_ARG_f);
50 outArgs.setSupports(MEB::OUT_ARG_W_op);
57 const RCP<const Teuchos::Comm<int> > comm =
60 const RCP<const Tpetra::Map<> > map =
rcp(
new Tpetra::Map<> (dim, 0, comm));
62 typedef Tpetra::Map<>::global_ordinal_type GO;
65 W_op_graph_->insertGlobalIndices(0, tuple<GO>(0, 1)());
66 W_op_graph_->insertGlobalIndices(1, tuple<GO>(0, 1)());
71 x0_ =
rcp(
new Tpetra::Vector<Scalar>(map));
72 x0_->putScalar(ST::zero());
78 x_space_ = Thyra::createVectorSpace<Scalar>(map);
89 set_p(Teuchos::tuple<Scalar>(2.0, 0.0)());
90 set_x0(Teuchos::tuple<Scalar>(1.0, 1.0)());
95 template<
class Scalar>
102 template<
class Scalar>
112 template<
class Scalar>
118 x0_->get1dViewNonConst()().assign(x0_in);
125 template<
class Scalar>
133 template<
class Scalar>
141 template<
class Scalar>
142 Thyra::ModelEvaluatorBase::InArgs<Scalar>
145 return nominalValues_;
149 template<
class Scalar>
155 Teuchos::rcp(
new Tpetra::CrsMatrix<Scalar,int>(W_op_graph_))
161 template<
class Scalar>
162 Thyra::ModelEvaluatorBase::InArgs<Scalar>
165 return prototypeInArgs_;
172 template<
class Scalar>
173 Thyra::ModelEvaluatorBase::OutArgs<Scalar>
176 return prototypeOutArgs_;
180 template<
class Scalar>
182 const Thyra::ModelEvaluatorBase::InArgs<Scalar> &inArgs,
183 const Thyra::ModelEvaluatorBase::OutArgs<Scalar> &outArgs
189 using Teuchos::tuple;
190 using Teuchos::rcp_dynamic_cast;
192 typedef Tpetra::Map<>::global_ordinal_type GO;
195 const RCP<const Tpetra::Vector<Scalar, int> > x_vec =
196 ConverterT::getConstTpetraVector(inArgs.get_x());
197 const ArrayRCP<const Scalar> x = x_vec->get1dView();
199 const RCP<Tpetra::Vector<Scalar, int> > f_vec =
200 ConverterT::getTpetraVector(outArgs.get_f());
202 const RCP<Tpetra::CrsMatrix<Scalar, int> > W =
203 rcp_dynamic_cast<Tpetra::CrsMatrix<Scalar,int> >(
204 ConverterT::getTpetraOperator(outArgs.get_W_op()),
209 const ArrayRCP<Scalar>
f = f_vec->get1dViewNonConst();
210 f[0] = x[0] + x[1]*x[1] - p_[0];
211 f[1] = d_ * (x[0]*x[0] -x[1] - p_[1]);
215 W->setAllToScalar(ST::zero());
216 W->sumIntoGlobalValues(0, tuple<GO>(0, 1), tuple<Scalar>(1.0, 2.0*x[1]));
217 W->sumIntoGlobalValues(1, tuple<GO>(0, 1), tuple<Scalar>(2.0*d_*x[0], -d_));
223 #endif // SIMPLE_2D_TPETRA_MODEL_EVALUATOR_HPP
void evalModelImpl(const Thyra::ModelEvaluatorBase::InArgs< Scalar > &inArgs, const Thyra::ModelEvaluatorBase::OutArgs< Scalar > &outArgs) const
Thyra::ModelEvaluatorBase::OutArgs< Scalar > prototypeOutArgs_
Thyra::ModelEvaluatorBase::InArgs< Scalar > prototypeInArgs_
Teuchos::RCP< Tpetra::CrsGraph<> > W_op_graph_
Thyra::ModelEvaluatorBase::InArgs< Scalar > createInArgs() const
static Teuchos::RCP< const Comm< OrdinalType > > getComm()
void set_x0(const Teuchos::ArrayView< const Scalar > &x0)
Teuchos::RCP< Thyra::LinearOpBase< Scalar > > create_W_op() const
void set_d(const Scalar &d)
RCP< LinearOpBase< Scalar > > createLinearOp(const RCP< Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraOperator, const RCP< const VectorSpaceBase< Scalar > > rangeSpace=Teuchos::null, const RCP< const VectorSpaceBase< Scalar > > domainSpace=Teuchos::null)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
RCP< VectorBase< Scalar > > createVector(const RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVector, const RCP< const VectorSpaceBase< Scalar > > space=Teuchos::null)
Teuchos::RCP< Tpetra::Vector< Scalar > > x0_
Teuchos::Array< Scalar > p_
void set_p(const Teuchos::ArrayView< const Scalar > &p)
void resize(size_type new_size, const value_type &x=value_type())
bool nonnull(const boost::shared_ptr< T > &p)
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_f_space() const
Thyra::ModelEvaluatorBase::InArgs< Scalar > getNominalValues() const
Simple2DTpetraModelEvaluator()
Thyra::ModelEvaluatorBase::OutArgs< Scalar > createOutArgsImpl() const
#define TEUCHOS_ASSERT_EQUALITY(val1, val2)
Thyra::ModelEvaluatorBase::InArgs< Scalar > nominalValues_
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > f_space_
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > get_x_space() const
Teuchos::RCP< const Thyra::VectorSpaceBase< Scalar > > x_space_