|
Thyra Package Browser (Single Doxygen Collection) Version of the Day
|
00001 00002 #include "Thyra_TpetraThyraWrappers.hpp" 00003 #include "Thyra_VectorSpaceTester.hpp" 00004 #include "Thyra_VectorStdOpsTester.hpp" 00005 #include "Thyra_LinearOpTester.hpp" 00006 #include "Thyra_DefaultProductVector.hpp" 00007 #include "Thyra_TestingTools.hpp" 00008 #include "Tpetra_CrsMatrix.hpp" 00009 #include "Teuchos_UnitTestHarness.hpp" 00010 #include "Teuchos_DefaultComm.hpp" 00011 00012 00013 namespace { 00014 00015 00016 // 00017 // Helper code and declarations 00018 // 00019 00020 00021 using Teuchos::as; 00022 using Teuchos::null; 00023 using Teuchos::RCP; 00024 using Teuchos::rcp; 00025 using Teuchos::ArrayView; 00026 using Teuchos::rcp_dynamic_cast; 00027 using Teuchos::inOutArg; 00028 using Teuchos::Comm; 00029 using Teuchos::tuple; 00030 typedef Thyra::Ordinal Ordinal; 00031 using Thyra::VectorSpaceBase; 00032 using Thyra::SpmdVectorSpaceBase; 00033 using Thyra::MultiVectorBase; 00034 using Thyra::VectorBase; 00035 using Thyra::LinearOpBase; 00036 using Thyra::createMember; 00037 00038 00039 const int g_localDim = 4; // ToDo: Make variable! 00040 00041 00042 typedef Tpetra::Map<int> TpetraMap_t; 00043 00044 00045 RCP<const TpetraMap_t> 00046 createTpetraMap(const int localDim) 00047 { 00048 typedef Teuchos::OrdinalTraits<Tpetra::global_size_t> OT; 00049 return Teuchos::rcp(new TpetraMap_t(OT::invalid(), localDim, 0, 00050 Teuchos::DefaultComm<int>::getComm())); 00051 // ToDo: Pass in the comm? 00052 } 00053 00054 // ToDo: Above: Vary the LocalOrdinal and GlobalOrdinal types? 00055 00056 00057 template<class Scalar> 00058 RCP<const VectorSpaceBase<Scalar> > 00059 createTpetraVectorSpace(const int localDim) 00060 { 00061 return Thyra::createVectorSpace<Scalar>(createTpetraMap(g_localDim)); 00062 } 00063 00064 00065 template<class Scalar> 00066 RCP<Tpetra::Operator<Scalar,int> > 00067 createTriDiagonalTpetraOperator(const int numLocalRows) 00068 { 00069 typedef int Ordinal; 00070 typedef Tpetra::global_size_t global_size_t; 00071 00072 RCP<const Tpetra::Map<int> > map = createTpetraMap(numLocalRows); 00073 00074 const size_t numMyElements = map->getNodeNumElements(); 00075 const global_size_t numGlobalElements = map->getGlobalNumElements(); 00076 00077 ArrayView<const Ordinal> myGlobalElements = map->getNodeElementList(); 00078 00079 // Create an OTeger vector numNz that is used to build the Petra Matrix. 00080 // numNz[i] is the Number of OFF-DIAGONAL term for the ith global equation 00081 // on this processor 00082 00083 Teuchos::ArrayRCP<size_t> numNz = Teuchos::arcp<size_t>(numMyElements); 00084 00085 // We are building a tridiagonal matrix where each row has (-1 2 -1) 00086 // So we need 2 off-diagonal terms (except for the first and last equation) 00087 00088 for (size_t i=0; i < numMyElements; ++i) { 00089 if (myGlobalElements[i] == 0 || static_cast<global_size_t>(myGlobalElements[i]) == numGlobalElements-1) { 00090 // boundary 00091 numNz[i] = 2; 00092 } 00093 else { 00094 numNz[i] = 3; 00095 } 00096 } 00097 00098 // Create a Tpetra::Matrix using the Map, with a static allocation dictated by numNz 00099 RCP< Tpetra::CrsMatrix<Scalar,Ordinal> > A = 00100 Teuchos::rcp( new Tpetra::CrsMatrix<Scalar,Ordinal>(map, numNz, Tpetra::StaticProfile) ); 00101 00102 // We are done with NumNZ 00103 numNz = Teuchos::null; 00104 00105 // Add rows one-at-a-time 00106 // Off diagonal values will always be -1 00107 const Scalar two = static_cast<Scalar>( 2.0); 00108 const Scalar posOne = static_cast<Scalar>(+1.0); 00109 const Scalar negOne = static_cast<Scalar>(-1.0); 00110 for (size_t i = 0; i < numMyElements; i++) { 00111 if (myGlobalElements[i] == 0) { 00112 A->insertGlobalValues( myGlobalElements[i], 00113 tuple<Ordinal>(myGlobalElements[i], myGlobalElements[i]+1)(), 00114 tuple<Scalar> (two, posOne)() 00115 ); 00116 } 00117 else if (static_cast<global_size_t>(myGlobalElements[i]) == numGlobalElements-1) { 00118 A->insertGlobalValues( myGlobalElements[i], 00119 tuple<Ordinal>(myGlobalElements[i]-1, myGlobalElements[i])(), 00120 tuple<Scalar> (negOne, two)() 00121 ); 00122 } 00123 else { 00124 A->insertGlobalValues( myGlobalElements[i], 00125 tuple<Ordinal>(myGlobalElements[i]-1, myGlobalElements[i], myGlobalElements[i]+1)(), 00126 tuple<Scalar> (negOne, two, posOne)() 00127 ); 00128 } 00129 } 00130 00131 // Finish up 00132 A->fillComplete(Tpetra::DoOptimizeStorage); 00133 00134 return A; 00135 00136 } 00137 00138 00139 bool showAllTests = false; 00140 bool dumpAll = false; 00141 bool runLinearOpTester = true; 00142 00143 00144 TEUCHOS_STATIC_SETUP() 00145 { 00146 Teuchos::UnitTestRepository::getCLP().setOption( 00147 "show-all-tests", "no-show-all-tests", &showAllTests, "Show all tests or not" ); 00148 Teuchos::UnitTestRepository::getCLP().setOption( 00149 "dump-all", "no-dump-all", &dumpAll, "Dump all objects being tested" ); 00150 Teuchos::UnitTestRepository::getCLP().setOption( 00151 "run-linear-op-tester", "no-run-linear-op-tester", &runLinearOpTester, "..." ); 00152 } 00153 00154 00155 // 00156 // Unit Tests 00157 // 00158 00159 00160 // 00161 // convertTpetraToThyraComm 00162 // 00163 00164 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, convertTpetraToThyraComm, 00165 Scalar ) 00166 { 00167 RCP<const Comm<int> > tpetraComm = Teuchos::DefaultComm<int>::getComm(); 00168 RCP<const Comm<Ordinal> > thyraComm = Thyra::convertTpetraToThyraComm(tpetraComm); 00169 TEST_ASSERT(nonnull(thyraComm)); 00170 } 00171 00172 00173 // 00174 // createVectorSpace 00175 // 00176 00177 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, createVectorSpace, 00178 Scalar ) 00179 { 00180 const RCP<const TpetraMap_t> tpetraMap = createTpetraMap(g_localDim); 00181 const RCP<const VectorSpaceBase<Scalar> > vs = 00182 Thyra::createVectorSpace<Scalar>(tpetraMap); 00183 TEST_ASSERT(nonnull(vs)); 00184 out << "vs = " << *vs; 00185 const RCP<const SpmdVectorSpaceBase<Scalar> > vs_spmd = 00186 rcp_dynamic_cast<const SpmdVectorSpaceBase<Scalar> >(vs, true); 00187 TEST_EQUALITY(vs_spmd->localSubDim(), g_localDim); 00188 TEST_EQUALITY(vs->dim(), as<Ordinal>(tpetraMap->getGlobalNumElements())); 00189 } 00190 00191 00192 // 00193 // createVector 00194 // 00195 00196 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, createVector, 00197 Scalar ) 00198 { 00199 00200 typedef Thyra::TpetraOperatorVectorExtraction<Scalar,int> ConverterT; 00201 00202 const RCP<const TpetraMap_t> tpetraMap = createTpetraMap(g_localDim); 00203 const RCP<const VectorSpaceBase<Scalar> > vs = 00204 Thyra::createVectorSpace<Scalar>(tpetraMap); 00205 00206 const RCP<Tpetra::Vector<Scalar,int> > tpetraVector = 00207 rcp(new Tpetra::Vector<Scalar,int>(tpetraMap)); 00208 00209 { 00210 const RCP<VectorBase<Scalar> > thyraVector = createVector(tpetraVector, vs); 00211 TEST_EQUALITY(thyraVector->space(), vs); 00212 const RCP<Tpetra::Vector<Scalar,int> > tpetraVector2 = 00213 ConverterT::getTpetraVector(thyraVector); 00214 TEST_EQUALITY(tpetraVector2, tpetraVector); 00215 } 00216 00217 { 00218 const RCP<VectorBase<Scalar> > thyraVector = Thyra::createVector(tpetraVector); 00219 TEST_INEQUALITY(thyraVector->space(), vs); 00220 TEST_ASSERT(thyraVector->space()->isCompatible(*vs)); 00221 const RCP<Tpetra::Vector<Scalar,int> > tpetraVector2 = 00222 ConverterT::getTpetraVector(thyraVector); 00223 TEST_EQUALITY(tpetraVector2, tpetraVector); 00224 } 00225 00226 } 00227 00228 00229 // 00230 // createConstVector 00231 // 00232 00233 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, createConstVector, 00234 Scalar ) 00235 { 00236 00237 typedef Thyra::TpetraOperatorVectorExtraction<Scalar,int> ConverterT; 00238 00239 const RCP<const TpetraMap_t> tpetraMap = createTpetraMap(g_localDim); 00240 const RCP<const VectorSpaceBase<Scalar> > vs = 00241 Thyra::createVectorSpace<Scalar>(tpetraMap); 00242 00243 const RCP<const Tpetra::Vector<Scalar,int> > tpetraVector = 00244 rcp(new Tpetra::Vector<Scalar,int>(tpetraMap)); 00245 00246 { 00247 const RCP<const VectorBase<Scalar> > thyraVector = 00248 createConstVector(tpetraVector, vs); 00249 TEST_EQUALITY(thyraVector->space(), vs); 00250 const RCP<const Tpetra::Vector<Scalar,int> > tpetraVector2 = 00251 ConverterT::getConstTpetraVector(thyraVector); 00252 TEST_EQUALITY(tpetraVector2, tpetraVector); 00253 } 00254 00255 { 00256 const RCP<const VectorBase<Scalar> > thyraVector = 00257 Thyra::createConstVector(tpetraVector); 00258 TEST_INEQUALITY(thyraVector->space(), vs); 00259 TEST_ASSERT(thyraVector->space()->isCompatible(*vs)); 00260 const RCP<const Tpetra::Vector<Scalar,int> > tpetraVector2 = 00261 ConverterT::getConstTpetraVector(thyraVector); 00262 TEST_EQUALITY(tpetraVector2, tpetraVector); 00263 } 00264 00265 } 00266 00267 00268 // 00269 // createMultiVector 00270 // 00271 00272 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, createMultiVector, 00273 Scalar ) 00274 { 00275 00276 typedef Thyra::TpetraOperatorVectorExtraction<Scalar,int> ConverterT; 00277 00278 const int numCols = 3; 00279 00280 const RCP<const TpetraMap_t> tpetraMap = createTpetraMap(g_localDim); 00281 const RCP<const VectorSpaceBase<Scalar> > rangeVs = 00282 Thyra::createVectorSpace<Scalar>(tpetraMap); 00283 00284 const RCP<const TpetraMap_t> tpetraLocRepMap = 00285 Tpetra::createLocalMapWithNode<int,int>( 00286 numCols, tpetraMap->getComm(), tpetraMap->getNode()); 00287 const RCP<const VectorSpaceBase<Scalar> > domainVs = 00288 Thyra::createVectorSpace<Scalar>(tpetraLocRepMap); 00289 00290 const RCP<Tpetra::MultiVector<Scalar,int> > tpetraMultiVector = 00291 rcp(new Tpetra::MultiVector<Scalar,int>(tpetraMap, numCols)); 00292 00293 { 00294 const RCP<MultiVectorBase<Scalar> > thyraMultiVector = 00295 createMultiVector(tpetraMultiVector, rangeVs, domainVs); 00296 TEST_EQUALITY(thyraMultiVector->range(), rangeVs); 00297 TEST_EQUALITY(thyraMultiVector->domain(), domainVs); 00298 const RCP<Tpetra::MultiVector<Scalar,int> > tpetraMultiVector2 = 00299 ConverterT::getTpetraMultiVector(thyraMultiVector); 00300 TEST_EQUALITY(tpetraMultiVector2, tpetraMultiVector); 00301 } 00302 00303 { 00304 const RCP<MultiVectorBase<Scalar> > thyraMultiVector = 00305 Thyra::createMultiVector(tpetraMultiVector); 00306 TEST_INEQUALITY(thyraMultiVector->range(), rangeVs); 00307 TEST_INEQUALITY(thyraMultiVector->domain(), domainVs); 00308 TEST_ASSERT(thyraMultiVector->range()->isCompatible(*rangeVs)); 00309 TEST_ASSERT(thyraMultiVector->domain()->isCompatible(*domainVs)); 00310 const RCP<Tpetra::MultiVector<Scalar,int> > tpetraMultiVector2 = 00311 ConverterT::getTpetraMultiVector(thyraMultiVector); 00312 TEST_EQUALITY(tpetraMultiVector2, tpetraMultiVector); 00313 } 00314 00315 } 00316 00317 00318 // 00319 // createConstMultiVector 00320 // 00321 00322 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, createConstMultiVector, 00323 Scalar ) 00324 { 00325 00326 typedef Thyra::TpetraOperatorVectorExtraction<Scalar,int> ConverterT; 00327 00328 const int numCols = 3; 00329 00330 const RCP<const TpetraMap_t> tpetraMap = createTpetraMap(g_localDim); 00331 const RCP<const VectorSpaceBase<Scalar> > rangeVs = 00332 Thyra::createVectorSpace<Scalar>(tpetraMap); 00333 00334 const RCP<const TpetraMap_t> tpetraLocRepMap = 00335 Tpetra::createLocalMapWithNode<int,int>( 00336 numCols, tpetraMap->getComm(), tpetraMap->getNode()); 00337 const RCP<const VectorSpaceBase<Scalar> > domainVs = 00338 Thyra::createVectorSpace<Scalar>(tpetraLocRepMap); 00339 00340 const RCP<const Tpetra::MultiVector<Scalar,int> > tpetraMultiVector = 00341 rcp(new Tpetra::MultiVector<Scalar,int>(tpetraMap, numCols)); 00342 00343 { 00344 const RCP<const MultiVectorBase<Scalar> > thyraMultiVector = 00345 createConstMultiVector(tpetraMultiVector, rangeVs, domainVs); 00346 TEST_EQUALITY(thyraMultiVector->range(), rangeVs); 00347 TEST_EQUALITY(thyraMultiVector->domain(), domainVs); 00348 const RCP<const Tpetra::MultiVector<Scalar,int> > tpetraMultiVector2 = 00349 ConverterT::getConstTpetraMultiVector(thyraMultiVector); 00350 TEST_EQUALITY(tpetraMultiVector2, tpetraMultiVector); 00351 } 00352 00353 { 00354 const RCP<const MultiVectorBase<Scalar> > thyraMultiVector = 00355 Thyra::createConstMultiVector(tpetraMultiVector); 00356 TEST_INEQUALITY(thyraMultiVector->range(), rangeVs); 00357 TEST_INEQUALITY(thyraMultiVector->domain(), domainVs); 00358 TEST_ASSERT(thyraMultiVector->range()->isCompatible(*rangeVs)); 00359 TEST_ASSERT(thyraMultiVector->domain()->isCompatible(*domainVs)); 00360 const RCP<const Tpetra::MultiVector<Scalar,int> > tpetraMultiVector2 = 00361 ConverterT::getConstTpetraMultiVector(thyraMultiVector); 00362 TEST_EQUALITY(tpetraMultiVector2, tpetraMultiVector); 00363 } 00364 00365 } 00366 00367 00368 // 00369 // TeptraVectorSpace 00370 // 00371 00372 00373 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, TeptraVectorSpace, 00374 Scalar ) 00375 { 00376 const RCP<const VectorSpaceBase<Scalar> > vs = 00377 Thyra::createVectorSpace<Scalar>(createTpetraMap(g_localDim)); 00378 const RCP<VectorBase<Scalar> > v = createMember(vs); 00379 TEST_ASSERT(nonnull(v)); 00380 TEST_EQUALITY(v->space(), vs); 00381 } 00382 00383 00384 // 00385 // vectorSpaceTester 00386 // 00387 00388 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, vectorSpaceTester, 00389 Scalar ) 00390 { 00391 const RCP<const VectorSpaceBase<Scalar> > vs 00392 = createTpetraVectorSpace<Scalar>(g_localDim); 00393 Thyra::VectorSpaceTester<Scalar> vectorSpaceTester; 00394 vectorSpaceTester.show_all_tests(showAllTests); 00395 vectorSpaceTester.dump_all(dumpAll); 00396 TEST_ASSERT(vectorSpaceTester.check(*vs, &out)); 00397 } 00398 00399 00400 // ToDo: Fix the default tolerances for below 00401 00402 /* 00403 00404 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraVectorSpace, vectorStdOpsTester, 00405 Scalar ) 00406 { 00407 const RCP<const VectorSpaceBase<Scalar> > vs = 00408 Thyra::createVectorSpace<Scalar>(createTpetraMap(g_localDim)); 00409 Thyra::VectorStdOpsTester<Scalar> vectorStdOpsTester; 00410 //vectorStdOpsTester.show_all_tests(showAllTests); 00411 //vectorStdOpsTester.dump_all(dumpAll); 00412 TEST_ASSERT(vectorStdOpsTester.checkStdOps(*vs, &out)); 00413 } 00414 00415 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT_SCALAR_TYPES( TpetraVectorSpace, 00416 vectorStdOpsTester ) 00417 00418 */ 00419 00420 00421 // 00422 // getTpetraMultiVector 00423 // 00424 00425 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, getTpetraMultiVector, 00426 Scalar ) 00427 { 00428 typedef Thyra::TpetraOperatorVectorExtraction<Scalar,int> ConverterT; 00429 00430 const int numCols = 3; 00431 const RCP<const VectorSpaceBase<Scalar> > vs 00432 = createTpetraVectorSpace<Scalar>(g_localDim); 00433 00434 { 00435 const RCP<MultiVectorBase<Scalar> > mv = createMembers(vs, numCols); 00436 const RCP<Tpetra::MultiVector<Scalar,int> > tmv = 00437 ConverterT::getTpetraMultiVector(mv); 00438 TEST_ASSERT(nonnull(tmv)); 00439 TEST_EQUALITY(as<Ordinal>(tmv->getMap()->getGlobalNumElements()), vs->dim()); 00440 } 00441 00442 { 00443 const RCP<VectorBase<Scalar> > v = createMember(vs); 00444 const RCP<Tpetra::MultiVector<Scalar,int> > tmv = 00445 ConverterT::getTpetraMultiVector(v); 00446 TEST_ASSERT(nonnull(tmv)); 00447 TEST_EQUALITY(as<Ordinal>(tmv->getMap()->getGlobalNumElements()), vs->dim()); 00448 } 00449 00450 #ifdef THYRA_DEBUG 00451 const RCP<VectorBase<Scalar> > pv = Thyra::defaultProductVector<Scalar>(); 00452 TEST_THROW(ConverterT::getTpetraMultiVector(pv), std::logic_error); 00453 #endif 00454 00455 } 00456 00457 00458 // 00459 // getConstTpetraMultiVector 00460 // 00461 00462 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, getConstTpetraMultiVector, 00463 Scalar ) 00464 { 00465 typedef Thyra::TpetraOperatorVectorExtraction<Scalar,int> ConverterT; 00466 00467 const int numCols = 3; 00468 const RCP<const VectorSpaceBase<Scalar> > vs 00469 = createTpetraVectorSpace<Scalar>(g_localDim); 00470 00471 { 00472 const RCP<const MultiVectorBase<Scalar> > mv = createMembers(vs, numCols); 00473 const RCP<const Tpetra::MultiVector<Scalar,int> > tmv = 00474 ConverterT::getConstTpetraMultiVector(mv); 00475 TEST_ASSERT(nonnull(tmv)); 00476 TEST_EQUALITY(as<Ordinal>(tmv->getMap()->getGlobalNumElements()), vs->dim()); 00477 } 00478 00479 { 00480 const RCP<const VectorBase<Scalar> > v = createMember(vs); 00481 const RCP<const Tpetra::MultiVector<Scalar,int> > tmv = 00482 ConverterT::getConstTpetraMultiVector(v); 00483 TEST_ASSERT(nonnull(tmv)); 00484 TEST_EQUALITY(as<Ordinal>(tmv->getMap()->getGlobalNumElements()), vs->dim()); 00485 } 00486 00487 #ifdef THYRA_DEBUG 00488 const RCP<const VectorBase<Scalar> > pv = Thyra::defaultProductVector<Scalar>(); 00489 TEST_THROW(ConverterT::getConstTpetraMultiVector(pv), std::logic_error); 00490 #endif 00491 00492 } 00493 00494 00495 // 00496 // TpetraLinearOp 00497 // 00498 00499 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, TpetraLinearOp, 00500 Scalar ) 00501 { 00502 00503 typedef Teuchos::ScalarTraits<Scalar> ST; 00504 using Teuchos::as; 00505 00506 const RCP<Tpetra::Operator<Scalar,int> > tpetraOp = 00507 createTriDiagonalTpetraOperator<Scalar>(g_localDim); 00508 out << "tpetraOp = " << Teuchos::describe(*tpetraOp, Teuchos::VERB_HIGH) << std::endl; 00509 TEST_ASSERT(nonnull(tpetraOp)); 00510 00511 const RCP<const VectorSpaceBase<Scalar> > rangeSpace = 00512 Thyra::createVectorSpace<Scalar>(tpetraOp->getRangeMap()); 00513 const RCP<const VectorSpaceBase<Scalar> > domainSpace = 00514 Thyra::createVectorSpace<Scalar>(tpetraOp->getDomainMap()); 00515 const RCP<const LinearOpBase<Scalar> > thyraLinearOp = 00516 Thyra::tpetraLinearOp(rangeSpace, domainSpace, tpetraOp); 00517 TEST_ASSERT(nonnull(thyraLinearOp)); 00518 00519 out << "\nCheck that operator returns the right thing ...\n"; 00520 const RCP<VectorBase<Scalar> > x = createMember(thyraLinearOp->domain()); 00521 Thyra::V_S(x.ptr(), ST::one()); 00522 const RCP<VectorBase<Scalar> > y = createMember(thyraLinearOp->range()); 00523 Thyra::apply<Scalar>(*thyraLinearOp, Thyra::NOTRANS, *x, y.ptr()); 00524 const Scalar sum_y = sum(*y); 00525 TEST_FLOATING_EQUALITY( sum_y, as<Scalar>(3+1+2*(y->space()->dim()-2)), 00526 100.0 * ST::eps() ); 00527 00528 out << "\nCheck the general LinearOp interface ...\n"; 00529 Thyra::LinearOpTester<Scalar> linearOpTester; 00530 linearOpTester.show_all_tests(showAllTests); 00531 linearOpTester.dump_all(dumpAll); 00532 if (runLinearOpTester) { 00533 TEST_ASSERT(linearOpTester.check(*thyraLinearOp, &out)); 00534 } 00535 00536 } 00537 00538 00539 // 00540 // createLinearOp 00541 // 00542 00543 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, createLinearOp, 00544 Scalar ) 00545 { 00546 00547 typedef Thyra::TpetraOperatorVectorExtraction<Scalar,int> ConverterT; 00548 00549 const RCP<Tpetra::Operator<Scalar,int> > tpetraOp = 00550 createTriDiagonalTpetraOperator<Scalar>(g_localDim); 00551 out << "tpetraOp = " << Teuchos::describe(*tpetraOp, Teuchos::VERB_HIGH) << std::endl; 00552 00553 const RCP<const VectorSpaceBase<Scalar> > rangeSpace = 00554 Thyra::createVectorSpace<Scalar>(tpetraOp->getRangeMap()); 00555 00556 const RCP<const VectorSpaceBase<Scalar> > domainSpace = 00557 Thyra::createVectorSpace<Scalar>(tpetraOp->getDomainMap()); 00558 00559 { 00560 const RCP<LinearOpBase<Scalar> > thyraOp = 00561 createLinearOp(tpetraOp, rangeSpace, domainSpace); 00562 TEST_EQUALITY(thyraOp->range(), rangeSpace); 00563 TEST_EQUALITY(thyraOp->domain(), domainSpace); 00564 const RCP<Tpetra::Operator<Scalar,int> > tpetraOp2 = 00565 ConverterT::getTpetraOperator(thyraOp); 00566 TEST_EQUALITY(tpetraOp2, tpetraOp); 00567 } 00568 00569 { 00570 const RCP<LinearOpBase<Scalar> > thyraOp = 00571 Thyra::createLinearOp(tpetraOp); 00572 TEST_INEQUALITY(thyraOp->range(), rangeSpace); 00573 TEST_INEQUALITY(thyraOp->domain(), domainSpace); 00574 TEST_ASSERT(thyraOp->range()->isCompatible(*rangeSpace)); 00575 TEST_ASSERT(thyraOp->domain()->isCompatible(*domainSpace)); 00576 const RCP<Tpetra::Operator<Scalar,int> > tpetraOp2 = 00577 ConverterT::getTpetraOperator(thyraOp); 00578 TEST_EQUALITY(tpetraOp2, tpetraOp); 00579 } 00580 00581 } 00582 00583 00584 // 00585 // createConstLinearOp 00586 // 00587 00588 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, createConstLinearOp, 00589 Scalar ) 00590 { 00591 00592 typedef Thyra::TpetraOperatorVectorExtraction<Scalar,int> ConverterT; 00593 00594 const RCP<const Tpetra::Operator<Scalar,int> > tpetraOp = 00595 createTriDiagonalTpetraOperator<Scalar>(g_localDim); 00596 out << "tpetraOp = " << Teuchos::describe(*tpetraOp, Teuchos::VERB_HIGH) << std::endl; 00597 00598 const RCP<const VectorSpaceBase<Scalar> > rangeSpace = 00599 Thyra::createVectorSpace<Scalar>(tpetraOp->getRangeMap()); 00600 00601 const RCP<const VectorSpaceBase<Scalar> > domainSpace = 00602 Thyra::createVectorSpace<Scalar>(tpetraOp->getDomainMap()); 00603 00604 { 00605 const RCP<const LinearOpBase<Scalar> > thyraOp = 00606 createConstLinearOp(tpetraOp, rangeSpace, domainSpace); 00607 TEST_EQUALITY(thyraOp->range(), rangeSpace); 00608 TEST_EQUALITY(thyraOp->domain(), domainSpace); 00609 const RCP<const Tpetra::Operator<Scalar,int> > tpetraOp2 = 00610 ConverterT::getConstTpetraOperator(thyraOp); 00611 TEST_EQUALITY(tpetraOp2, tpetraOp); 00612 } 00613 00614 { 00615 const RCP<const LinearOpBase<Scalar> > thyraOp = 00616 Thyra::createConstLinearOp(tpetraOp); 00617 TEST_INEQUALITY(thyraOp->range(), rangeSpace); 00618 TEST_INEQUALITY(thyraOp->domain(), domainSpace); 00619 TEST_ASSERT(thyraOp->range()->isCompatible(*rangeSpace)); 00620 TEST_ASSERT(thyraOp->domain()->isCompatible(*domainSpace)); 00621 const RCP<const Tpetra::Operator<Scalar,int> > tpetraOp2 = 00622 ConverterT::getConstTpetraOperator(thyraOp); 00623 TEST_EQUALITY(tpetraOp2, tpetraOp); 00624 } 00625 00626 } 00627 00628 00629 // 00630 // TpetraLinearOp_EpetraRowMatrix 00631 // 00632 00633 00634 #ifdef HAVE_THYRA_TPETRA_EPETRA 00635 00636 00637 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, TpetraLinearOp_EpetraRowMatrix, 00638 Scalar ) 00639 { 00640 00641 using Teuchos::as; 00642 using Teuchos::outArg; 00643 using Teuchos::rcp_dynamic_cast; 00644 using Teuchos::Array; 00645 typedef Teuchos::ScalarTraits<Scalar> ST; 00646 00647 const RCP<Tpetra::Operator<Scalar,int> > tpetraOp = 00648 createTriDiagonalTpetraOperator<Scalar>(g_localDim); 00649 00650 const RCP<LinearOpBase<Scalar> > thyraOp = 00651 Thyra::createLinearOp(tpetraOp); 00652 00653 const RCP<Thyra::TpetraLinearOp<Scalar, int> > thyraTpetraOp = 00654 Teuchos::rcp_dynamic_cast<Thyra::TpetraLinearOp<Scalar, int> >(thyraOp); 00655 00656 RCP<const Epetra_Operator> epetraOp; 00657 Thyra::EOpTransp epetraOpTransp; 00658 Thyra::EApplyEpetraOpAs epetraOpApplyAs; 00659 Thyra::EAdjointEpetraOp epetraOpAdjointSupport; 00660 00661 thyraTpetraOp->getEpetraOpView( outArg(epetraOp), outArg(epetraOpTransp), 00662 outArg(epetraOpApplyAs), outArg(epetraOpAdjointSupport) ); 00663 00664 if (typeid(Scalar) == typeid(double)) { 00665 TEST_ASSERT(nonnull(epetraOp)); 00666 const RCP<const Epetra_RowMatrix> epetraRowMatrix = 00667 rcp_dynamic_cast<const Epetra_RowMatrix>(epetraOp, true); 00668 int numRowEntries = -1; 00669 epetraRowMatrix->NumMyRowEntries(1, numRowEntries); 00670 TEST_EQUALITY_CONST(numRowEntries, 3); 00671 Array<double> row_values(numRowEntries); 00672 Array<int> row_indices(numRowEntries); 00673 epetraRowMatrix->ExtractMyRowCopy(1, numRowEntries, numRowEntries, 00674 row_values.getRawPtr(), row_indices.getRawPtr()); 00675 TEST_EQUALITY_CONST(row_values[0], -1.0); 00676 TEST_EQUALITY_CONST(row_values[1], 2.0); 00677 TEST_EQUALITY_CONST(row_values[2], 1.0); 00678 // ToDo: Test column indices! 00679 } 00680 else { 00681 TEST_ASSERT(is_null(epetraOp)); 00682 } 00683 00684 } 00685 00686 00687 #endif // HAVE_THYRA_TPETRA_EPETRA 00688 00689 00690 // 00691 // Unit test instantiations 00692 // 00693 00694 #define THYRA_TPETRA_THYRA_WRAPPERS_INSTANT(SCALAR) \ 00695 \ 00696 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 00697 convertTpetraToThyraComm, SCALAR ) \ 00698 \ 00699 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 00700 createVectorSpace, SCALAR ) \ 00701 \ 00702 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 00703 createVector, SCALAR ) \ 00704 \ 00705 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 00706 createConstVector, SCALAR ) \ 00707 \ 00708 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 00709 createMultiVector, SCALAR ) \ 00710 \ 00711 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 00712 createConstMultiVector, SCALAR ) \ 00713 \ 00714 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 00715 TeptraVectorSpace, SCALAR ) \ 00716 \ 00717 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 00718 vectorSpaceTester, SCALAR ) \ 00719 \ 00720 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 00721 getTpetraMultiVector, SCALAR ) \ 00722 \ 00723 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 00724 getConstTpetraMultiVector, SCALAR ) \ 00725 \ 00726 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 00727 TpetraLinearOp, SCALAR ) \ 00728 \ 00729 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 00730 createLinearOp, SCALAR ) \ 00731 \ 00732 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 00733 createConstLinearOp, SCALAR ) \ 00734 \ 00735 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 00736 TpetraLinearOp_EpetraRowMatrix, SCALAR ) \ 00737 00738 00739 // We can currently only explicitly instantiate with double support because 00740 // Tpetra only supports explicit instantaition with double. As for implicit 00741 // instantation, g++ 3.4.6 on my Linux machine was taking more than 30 minutes 00742 // to compile this file when all of the types double, float, complex<double>, 00743 // and complex<float> where enabled. Therefore, we will only test double for 00744 // now until explicit instantation with other types are supported by Tpetra. 00745 00746 THYRA_TPETRA_THYRA_WRAPPERS_INSTANT(double) 00747 00748 00749 } // namespace
1.7.4