|
Thyra Package Browser (Single Doxygen Collection) Version of the Day
|
00001 00002 #include "Simple2DTpetraModelEvaluator.hpp" 00003 00004 #include "Teuchos_UnitTestHarness.hpp" 00005 00006 00007 namespace { 00008 00009 00010 using Teuchos::null; 00011 using Teuchos::RCP; 00012 typedef Thyra::ModelEvaluatorBase MEB; 00013 00014 00015 // 00016 // Unit tests 00017 // 00018 00019 00020 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( Simple2DTpetraModelEvaluator, construct, Scalar ) 00021 { 00022 RCP<Simple2DTpetraModelEvaluator<Scalar> > model = simple2DTpetraModelEvaluator<Scalar>(); 00023 TEST_ASSERT(model != null); 00024 TEST_EQUALITY(model->Np(), 0); 00025 TEST_EQUALITY(model->Ng(), 0); 00026 TEST_ASSERT(model->get_x_space() != null); 00027 TEST_EQUALITY(model->get_x_space()->dim(), 2); 00028 TEST_ASSERT(model->get_f_space() != null); 00029 TEST_EQUALITY(model->get_f_space()->dim(), 2); 00030 TEST_ASSERT(nonnull(model->getNominalValues().get_x())); 00031 TEST_ASSERT(model->create_W_op() != null); 00032 //TEST_ASSERT(model->get_W_factory() != null); 00033 MEB::InArgs<Scalar> inArgs = model->createInArgs(); 00034 TEST_ASSERT(inArgs.supports(MEB::IN_ARG_x)); 00035 TEST_EQUALITY(inArgs.Np(), 0); 00036 } 00037 00038 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT_REAL_SCALAR_TYPES( 00039 Simple2DTpetraModelEvaluator, construct ) 00040 00041 00042 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( Simple2DTpetraModelEvaluator, eval, Scalar ) 00043 { 00044 00045 using Teuchos::tuple; 00046 using Teuchos::ArrayRCP; 00047 using Teuchos::as; 00048 using Teuchos::rcp_dynamic_cast; 00049 typedef Teuchos::ScalarTraits<Scalar> ST; 00050 typedef typename ST::magnitudeType ScalarMag; 00051 typedef Teuchos::ScalarTraits<ScalarMag> SMT; 00052 using Thyra::VectorBase; 00053 using Thyra::LinearOpBase; 00054 typedef Thyra::ModelEvaluatorBase MEB; 00055 typedef Thyra::TpetraOperatorVectorExtraction<Scalar,int> ConverterT; 00056 00057 RCP<Simple2DTpetraModelEvaluator<Scalar> > model = simple2DTpetraModelEvaluator<Scalar>(); 00058 00059 const Scalar d = 7.0; 00060 model->set_d(d); 00061 00062 const Scalar p_0 = 2.0, p_1 = 6.0; 00063 model->set_p(tuple<Scalar>(p_0, p_1)); 00064 00065 const Scalar x_0 = 3.0, x_1 = 4.0; 00066 model->set_x0(tuple<Scalar>(x_0, x_1)); 00067 00068 MEB::InArgs<Scalar> inArgs = model->getNominalValues(); 00069 MEB::OutArgs<Scalar> outArgs = model->createOutArgs(); 00070 const RCP<VectorBase<Scalar> > f = createMember(model->get_f_space()); 00071 const RCP<LinearOpBase<Scalar> > W_op = model->create_W_op(); 00072 outArgs.set_f(f); 00073 outArgs.set_W_op(W_op); 00074 model->evalModel(inArgs, outArgs); 00075 00076 const ScalarMag tol = 100.0 * SMT::eps(); 00077 00078 const RCP<const Tpetra::Vector<Scalar,int> > f_tpetra = 00079 ConverterT::getConstTpetraVector(f); 00080 00081 const ArrayRCP<const Scalar> f_tpetra_vals = f_tpetra->get1dView(); 00082 TEST_FLOATING_EQUALITY(f_tpetra_vals[0], as<Scalar>(x_0 + x_1*x_1 - p_0), tol); 00083 TEST_FLOATING_EQUALITY(f_tpetra_vals[1], as<Scalar>(d * (x_0*x_0 - x_1 - p_1)), tol); 00084 00085 const RCP<const Tpetra::CrsMatrix<Scalar,int> > W_tpetra = 00086 rcp_dynamic_cast<Tpetra::CrsMatrix<Scalar,int> >( 00087 ConverterT::getTpetraOperator(W_op)); 00088 00089 ArrayRCP<const int> row_indices; 00090 ArrayRCP<const Scalar> row_values; 00091 00092 W_tpetra->getLocalRowView(0, row_indices, row_values); 00093 TEST_EQUALITY( row_indices[0], 0 ); 00094 TEST_EQUALITY( row_indices[1], 1 ); 00095 TEST_FLOATING_EQUALITY( row_values[0], as<Scalar>(1.0), tol ); 00096 TEST_FLOATING_EQUALITY( row_values[1], as<Scalar>(2.0*x_1), tol ); 00097 00098 W_tpetra->getLocalRowView(1, row_indices, row_values); 00099 TEST_EQUALITY( row_indices[0], 0 ); 00100 TEST_EQUALITY( row_indices[1], 1 ); 00101 TEST_FLOATING_EQUALITY( row_values[0], as<Scalar>(d*2.0*x_0), tol ); 00102 TEST_FLOATING_EQUALITY( row_values[1], as<Scalar>(-d), tol ); 00103 00104 } 00105 00106 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT_REAL_SCALAR_TYPES( 00107 Simple2DTpetraModelEvaluator, eval ) 00108 00109 00110 } // namespace
1.7.4