Thyra Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Simple2DTpetraModelEvaluator_UnitTests.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Thyra: Interfaces and Support for Abstract Numerical Algorithms
4 //
5 // Copyright 2004 NTESS and the Thyra contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
11 
13 
14 
15 namespace {
16 
17 
18 using Teuchos::null;
19 using Teuchos::RCP;
20 typedef Thyra::ModelEvaluatorBase MEB;
21 
22 
23 //
24 // Unit tests
25 //
26 
27 
29 {
30  RCP<Simple2DTpetraModelEvaluator<Scalar> > model = simple2DTpetraModelEvaluator<Scalar>();
31  TEST_ASSERT(model != null);
32  TEST_EQUALITY(model->Np(), 0);
33  TEST_EQUALITY(model->Ng(), 0);
34  TEST_ASSERT(model->get_x_space() != null);
35  TEST_EQUALITY(model->get_x_space()->dim(), 2);
36  TEST_ASSERT(model->get_f_space() != null);
37  TEST_EQUALITY(model->get_f_space()->dim(), 2);
38  TEST_ASSERT(nonnull(model->getNominalValues().get_x()));
39  TEST_ASSERT(model->create_W_op() != null);
40  //TEST_ASSERT(model->get_W_factory() != null);
41  MEB::InArgs<Scalar> inArgs = model->createInArgs();
42  TEST_ASSERT(inArgs.supports(MEB::IN_ARG_x));
43  TEST_EQUALITY(inArgs.Np(), 0);
44 }
45 
46 #if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_DOUBLE)
48  Simple2DTpetraModelEvaluator, construct )
49 #endif
50 
51 #if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
53  Simple2DTpetraModelEvaluator, construct )
54 #endif
55 
56 
58 {
59  using Thyra::LinearOpBase;
60  using Thyra::VectorBase;
61  using Teuchos::ArrayRCP;
62  using Teuchos::ArrayView;
63  using Teuchos::as;
64  using Teuchos::rcp_dynamic_cast;
65  using Teuchos::tuple;
67  typedef Tpetra::Map<>::local_ordinal_type LO;
69  typedef typename ST::magnitudeType ScalarMag;
71 
72  RCP<Simple2DTpetraModelEvaluator<Scalar> > model = simple2DTpetraModelEvaluator<Scalar>();
73 
74  const Scalar d = 7.0;
75  model->set_d(d);
76 
77  const Scalar p_0 = 2.0, p_1 = 6.0;
78  model->set_p(tuple<Scalar>(p_0, p_1));
79 
80  const Scalar x_0 = 3.0, x_1 = 4.0;
81  model->set_x0(tuple<Scalar>(x_0, x_1));
82 
83  MEB::InArgs<Scalar> inArgs = model->getNominalValues();
84  MEB::OutArgs<Scalar> outArgs = model->createOutArgs();
85  const RCP<VectorBase<Scalar> > f = createMember(model->get_f_space());
86  const RCP<LinearOpBase<Scalar> > W_op = model->create_W_op();
87  outArgs.set_f(f);
88  outArgs.set_W_op(W_op);
89  model->evalModel(inArgs, outArgs);
90 
91  const ScalarMag tol = 100.0 * SMT::eps();
92 
93  const RCP<const Tpetra::Vector<Scalar> > f_tpetra =
94  ConverterT::getConstTpetraVector(f);
95 
96  const ArrayRCP<const Scalar> f_tpetra_vals = f_tpetra->get1dView();
97  TEST_FLOATING_EQUALITY(f_tpetra_vals[0], as<Scalar>(x_0 + x_1*x_1 - p_0), tol);
98  TEST_FLOATING_EQUALITY(f_tpetra_vals[1], as<Scalar>(d * (x_0*x_0 - x_1 - p_1)), tol);
99 
100  const RCP<const Tpetra::CrsMatrix<Scalar> > W_tpetra =
101  rcp_dynamic_cast<Tpetra::CrsMatrix<Scalar> >(
102  ConverterT::getTpetraOperator(W_op));
103  using crs_t = Tpetra::CrsMatrix<Scalar>;
104  typename crs_t::local_inds_host_view_type row_indices;
105  typename crs_t::values_host_view_type row_values;
106 
107  W_tpetra->getLocalRowView(0, row_indices, row_values);
108  // FIXME (mfh 22 Oct 2015) This test assumes that local indices
109  // occur in a particular order. Tpetra does not necessarily
110  // guarantee this.
111  TEST_EQUALITY( row_indices[0], static_cast<LO> (0) );
112  TEST_EQUALITY( row_indices[1], static_cast<LO> (1) );
113  TEST_FLOATING_EQUALITY( row_values[0], as<Scalar>(1.0), tol );
114  TEST_FLOATING_EQUALITY( row_values[1], as<Scalar>(2.0*x_1), tol );
115 
116  W_tpetra->getLocalRowView(1, row_indices, row_values);
117  TEST_EQUALITY( row_indices[0], static_cast<LO> (0) );
118  TEST_EQUALITY( row_indices[1], static_cast<LO> (1) );
119  TEST_FLOATING_EQUALITY( row_values[0], as<Scalar>(d*2.0*x_0), tol );
120  TEST_FLOATING_EQUALITY( row_values[1], as<Scalar>(-d), tol );
121 
122 }
123 
124 #if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_DOUBLE)
127 #endif
128 
129 #if !defined(HAVE_TPETRA_EXPLICIT_INSTANTIATION) || defined(HAVE_TPETRA_INST_FLOAT)
132 #endif
133 
134 
135 
136 } // namespace
#define TEST_ASSERT(v1)
void f()
Traits class that enables the extraction of Tpetra operator/vector objects wrapped in Thyra operator/...
#define TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT_FLOAT(TEST_GROUP, TEST_NAME)
TypeTo as(const TypeFrom &t)
#define TEST_FLOATING_EQUALITY(v1, v2, tol)
bool nonnull(const boost::shared_ptr< T > &p)
Simple 2d simulation only ModelEvaluator for f(x) = 0 using Tpetra objects.
TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL(TpetraThyraWrappers, convertTpetraToThyraComm, Scalar)
#define TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT_DOUBLE(TEST_GROUP, TEST_NAME)
#define TEST_EQUALITY(v1, v2)