46 #include "Thyra_VectorSpaceTester.hpp" 47 #include "Thyra_VectorStdOpsTester.hpp" 48 #include "Thyra_MultiVectorStdOpsTester.hpp" 49 #include "Thyra_VectorStdOps.hpp" 50 #include "Thyra_MultiVectorStdOps.hpp" 51 #include "Thyra_LinearOpTester.hpp" 52 #include "Thyra_DefaultProductVector.hpp" 53 #include "Thyra_TestingTools.hpp" 54 #include "Thyra_ScaledLinearOpBase.hpp" 55 #include "Thyra_RowStatLinearOpBase.hpp" 56 #include "Thyra_VectorStdOps.hpp" 57 #include "Tpetra_CrsMatrix.hpp" 58 #include "Teuchos_UnitTestHarness.hpp" 59 #include "Teuchos_DefaultComm.hpp" 60 #include "Teuchos_Tuple.hpp" 75 using Teuchos::ArrayView;
76 using Teuchos::rcp_dynamic_cast;
77 using Teuchos::inOutArg;
88 RCP<const TpetraMap_t>
91 typedef Teuchos::OrdinalTraits<Tpetra::global_size_t> OT;
92 return Teuchos::rcp(
new TpetraMap_t(OT::invalid(), localDim, 0,
93 Teuchos::DefaultComm<int>::getComm()));
100 template<
class Scalar>
101 RCP<const VectorSpaceBase<Scalar> >
108 template<
class Scalar>
109 RCP<Tpetra::Operator<Scalar> >
112 typedef Tpetra::global_size_t global_size_t;
113 typedef Tpetra::Map<>::global_ordinal_type GO;
117 const size_t numMyElements = map->getNodeNumElements();
118 const global_size_t numGlobalElements = map->getGlobalNumElements();
120 ArrayView<const GO> myGlobalElements = map->getNodeElementList();
126 Teuchos::ArrayRCP<size_t> numNz = Teuchos::arcp<size_t>(numMyElements);
131 for (
size_t i=0; i < numMyElements; ++i) {
132 if (myGlobalElements[i] == 0 || static_cast<global_size_t>(myGlobalElements[i]) == numGlobalElements-1) {
142 RCP< Tpetra::CrsMatrix<Scalar> > A =
143 Teuchos::rcp(
new Tpetra::CrsMatrix<Scalar>(map, numNz, Tpetra::StaticProfile) );
146 numNz = Teuchos::null;
150 const Scalar two =
static_cast<Scalar
>( 2.0);
151 const Scalar posOne =
static_cast<Scalar
>(+1.0);
152 const Scalar negOne =
static_cast<Scalar
>(-1.0);
153 for (
size_t i = 0; i < numMyElements; i++) {
154 if (myGlobalElements[i] == 0) {
155 A->insertGlobalValues( myGlobalElements[i],
156 tuple<GO>(myGlobalElements[i], myGlobalElements[i]+1)(),
157 tuple<Scalar> (two, posOne)()
160 else if (static_cast<global_size_t>(myGlobalElements[i]) == numGlobalElements-1) {
161 A->insertGlobalValues( myGlobalElements[i],
162 tuple<GO>(myGlobalElements[i]-1, myGlobalElements[i])(),
163 tuple<Scalar> (negOne, two)()
167 A->insertGlobalValues( myGlobalElements[i],
168 tuple<GO>(myGlobalElements[i]-1, myGlobalElements[i], myGlobalElements[i]+1)(),
169 tuple<Scalar> (negOne, two, posOne)()
189 Teuchos::UnitTestRepository::getCLP().setOption(
190 "show-all-tests",
"no-show-all-tests", &
showAllTests,
"Show all tests or not" );
191 Teuchos::UnitTestRepository::getCLP().setOption(
192 "dump-all",
"no-dump-all", &
dumpAll,
"Dump all objects being tested" );
193 Teuchos::UnitTestRepository::getCLP().setOption(
210 RCP<const Comm<int> > tpetraComm = Teuchos::DefaultComm<int>::getComm();
212 TEST_ASSERT(nonnull(thyraComm));
224 const RCP<const VectorSpaceBase<Scalar> > vs =
225 Thyra::createVectorSpace<Scalar>(tpetraMap);
226 TEST_ASSERT(nonnull(vs));
227 out <<
"vs = " << *vs;
228 const RCP<const SpmdVectorSpaceBase<Scalar> > vs_spmd =
229 rcp_dynamic_cast<
const SpmdVectorSpaceBase<Scalar> >(vs,
true);
230 TEST_EQUALITY(vs_spmd->localSubDim(),
g_localDim);
231 TEST_EQUALITY(vs->dim(), as<Ordinal>(tpetraMap->getGlobalNumElements()));
246 const RCP<const VectorSpaceBase<Scalar> > vs =
247 Thyra::createVectorSpace<Scalar>(tpetraMap);
249 const RCP<Tpetra::Vector<Scalar> > tpetraVector =
250 rcp(
new Tpetra::Vector<Scalar>(tpetraMap));
253 const RCP<VectorBase<Scalar> > thyraVector =
createVector(tpetraVector, vs);
254 TEST_EQUALITY(thyraVector->space(), vs);
255 const RCP<Tpetra::Vector<Scalar> > tpetraVector2 =
256 ConverterT::getTpetraVector(thyraVector);
257 TEST_EQUALITY(tpetraVector2, tpetraVector);
262 TEST_INEQUALITY(thyraVector->space(), vs);
263 TEST_ASSERT(thyraVector->space()->isCompatible(*vs));
264 const RCP<Tpetra::Vector<Scalar> > tpetraVector2 =
265 ConverterT::getTpetraVector(thyraVector);
266 TEST_EQUALITY(tpetraVector2, tpetraVector);
283 const RCP<const VectorSpaceBase<Scalar> > vs =
284 Thyra::createVectorSpace<Scalar>(tpetraMap);
286 const RCP<const Tpetra::Vector<Scalar> > tpetraVector =
287 rcp(
new Tpetra::Vector<Scalar>(tpetraMap));
290 const RCP<const VectorBase<Scalar> > thyraVector =
292 TEST_EQUALITY(thyraVector->space(), vs);
293 const RCP<const Tpetra::Vector<Scalar> > tpetraVector2 =
294 ConverterT::getConstTpetraVector(thyraVector);
295 TEST_EQUALITY(tpetraVector2, tpetraVector);
299 const RCP<const VectorBase<Scalar> > thyraVector =
301 TEST_INEQUALITY(thyraVector->space(), vs);
302 TEST_ASSERT(thyraVector->space()->isCompatible(*vs));
303 const RCP<const Tpetra::Vector<Scalar> > tpetraVector2 =
304 ConverterT::getConstTpetraVector(thyraVector);
305 TEST_EQUALITY(tpetraVector2, tpetraVector);
318 typedef Tpetra::Map<>::local_ordinal_type LO;
319 typedef Tpetra::Map<>::global_ordinal_type GO;
322 const int numCols = 3;
325 const RCP<const VectorSpaceBase<Scalar> > rangeVs =
326 Thyra::createVectorSpace<Scalar>(tpetraMap);
328 const RCP<const TpetraMap_t> tpetraLocRepMap =
329 Tpetra::createLocalMapWithNode<LO, GO>(
330 numCols, tpetraMap->getComm(), tpetraMap->getNode());
331 const RCP<const VectorSpaceBase<Scalar> > domainVs =
332 Thyra::createVectorSpace<Scalar>(tpetraLocRepMap);
334 const RCP<Tpetra::MultiVector<Scalar> > tpetraMultiVector =
335 rcp(
new Tpetra::MultiVector<Scalar>(tpetraMap, numCols));
338 const RCP<MultiVectorBase<Scalar> > thyraMultiVector =
340 TEST_EQUALITY(thyraMultiVector->range(), rangeVs);
341 TEST_EQUALITY(thyraMultiVector->domain(), domainVs);
342 const RCP<Tpetra::MultiVector<Scalar> > tpetraMultiVector2 =
343 ConverterT::getTpetraMultiVector(thyraMultiVector);
344 TEST_EQUALITY(tpetraMultiVector2, tpetraMultiVector);
348 const RCP<MultiVectorBase<Scalar> > thyraMultiVector =
350 TEST_INEQUALITY(thyraMultiVector->range(), rangeVs);
351 TEST_INEQUALITY(thyraMultiVector->domain(), domainVs);
352 TEST_ASSERT(thyraMultiVector->range()->isCompatible(*rangeVs));
353 TEST_ASSERT(thyraMultiVector->domain()->isCompatible(*domainVs));
354 const RCP<Tpetra::MultiVector<Scalar> > tpetraMultiVector2 =
355 ConverterT::getTpetraMultiVector(thyraMultiVector);
356 TEST_EQUALITY(tpetraMultiVector2, tpetraMultiVector);
369 typedef Tpetra::Map<>::local_ordinal_type LO;
370 typedef Tpetra::Map<>::global_ordinal_type GO;
373 const int numCols = 3;
376 const RCP<const VectorSpaceBase<Scalar> > rangeVs =
377 Thyra::createVectorSpace<Scalar>(tpetraMap);
379 const RCP<const TpetraMap_t> tpetraLocRepMap =
380 Tpetra::createLocalMapWithNode<LO,GO>(
381 numCols, tpetraMap->getComm(), tpetraMap->getNode());
382 const RCP<const VectorSpaceBase<Scalar> > domainVs =
383 Thyra::createVectorSpace<Scalar>(tpetraLocRepMap);
385 const RCP<const Tpetra::MultiVector<Scalar> > tpetraMultiVector =
386 rcp(
new Tpetra::MultiVector<Scalar>(tpetraMap, numCols));
389 const RCP<const MultiVectorBase<Scalar> > thyraMultiVector =
391 TEST_EQUALITY(thyraMultiVector->range(), rangeVs);
392 TEST_EQUALITY(thyraMultiVector->domain(), domainVs);
393 const RCP<const Tpetra::MultiVector<Scalar> > tpetraMultiVector2 =
394 ConverterT::getConstTpetraMultiVector(thyraMultiVector);
395 TEST_EQUALITY(tpetraMultiVector2, tpetraMultiVector);
399 const RCP<const MultiVectorBase<Scalar> > thyraMultiVector =
401 TEST_INEQUALITY(thyraMultiVector->range(), rangeVs);
402 TEST_INEQUALITY(thyraMultiVector->domain(), domainVs);
403 TEST_ASSERT(thyraMultiVector->range()->isCompatible(*rangeVs));
404 TEST_ASSERT(thyraMultiVector->domain()->isCompatible(*domainVs));
405 const RCP<const Tpetra::MultiVector<Scalar> > tpetraMultiVector2 =
406 ConverterT::getConstTpetraMultiVector(thyraMultiVector);
407 TEST_EQUALITY(tpetraMultiVector2, tpetraMultiVector);
421 const RCP<const VectorSpaceBase<Scalar> > vs =
423 const RCP<VectorBase<Scalar> > v = createMember(vs);
424 TEST_ASSERT(nonnull(v));
425 TEST_EQUALITY(v->space(), vs);
436 const RCP<const VectorSpaceBase<Scalar> > vs
437 = createTpetraVectorSpace<Scalar>(
g_localDim);
438 Thyra::VectorSpaceTester<Scalar> vectorSpaceTester;
440 vectorSpaceTester.dump_all(
dumpAll);
441 TEST_ASSERT(vectorSpaceTester.check(*vs, &out));
453 const RCP<const VectorSpaceBase<Scalar> > vs =
455 Thyra::VectorStdOpsTester<Scalar> vectorStdOpsTester;
456 vectorStdOpsTester.warning_tol(5.0e-13);
457 vectorStdOpsTester.error_tol(5.0e-14);
458 TEST_ASSERT(vectorStdOpsTester.checkStdOps(*vs, &out));
470 const RCP<const VectorSpaceBase<Scalar> > vs =
472 Thyra::MultiVectorStdOpsTester<Scalar> mvStdOpsTester;
473 mvStdOpsTester.warning_tol(5.0e-13);
474 mvStdOpsTester.error_tol(5.0e-14);
475 TEST_ASSERT(mvStdOpsTester.checkStdOps(*vs, &out));
488 const int numCols = 3;
489 const RCP<const VectorSpaceBase<Scalar> > vs
490 = createTpetraVectorSpace<Scalar>(
g_localDim);
493 const RCP<MultiVectorBase<Scalar> > mv = createMembers(vs, numCols);
494 const RCP<Tpetra::MultiVector<Scalar> > tmv =
495 ConverterT::getTpetraMultiVector(mv);
496 TEST_ASSERT(nonnull(tmv));
497 TEST_EQUALITY(as<Ordinal>(tmv->getMap()->getGlobalNumElements()), vs->dim());
501 const RCP<VectorBase<Scalar> > v = createMember(vs);
502 const RCP<Tpetra::MultiVector<Scalar> > tmv =
503 ConverterT::getTpetraMultiVector(v);
504 TEST_ASSERT(nonnull(tmv));
505 TEST_EQUALITY(as<Ordinal>(tmv->getMap()->getGlobalNumElements()), vs->dim());
509 const RCP<VectorBase<Scalar> > pv = Thyra::defaultProductVector<Scalar>();
510 TEST_THROW(ConverterT::getTpetraMultiVector(pv), std::logic_error);
525 const int numCols = 3;
526 const RCP<const VectorSpaceBase<Scalar> > vs
527 = createTpetraVectorSpace<Scalar>(
g_localDim);
530 const RCP<const MultiVectorBase<Scalar> > mv = createMembers(vs, numCols);
531 const RCP<const Tpetra::MultiVector<Scalar> > tmv =
532 ConverterT::getConstTpetraMultiVector(mv);
533 TEST_ASSERT(nonnull(tmv));
534 TEST_EQUALITY(as<Ordinal>(tmv->getMap()->getGlobalNumElements()), vs->dim());
538 const RCP<const VectorBase<Scalar> > v = createMember(vs);
539 const RCP<const Tpetra::MultiVector<Scalar> > tmv =
540 ConverterT::getConstTpetraMultiVector(v);
541 TEST_ASSERT(nonnull(tmv));
542 TEST_EQUALITY(as<Ordinal>(tmv->getMap()->getGlobalNumElements()), vs->dim());
546 const RCP<const VectorBase<Scalar> > pv = Thyra::defaultProductVector<Scalar>();
547 TEST_THROW(ConverterT::getConstTpetraMultiVector(pv), std::logic_error);
561 typedef Teuchos::ScalarTraits<Scalar> ST;
564 const RCP<Tpetra::Operator<Scalar> > tpetraOp =
565 createTriDiagonalTpetraOperator<Scalar>(
g_localDim);
566 out <<
"tpetraOp = " << Teuchos::describe(*tpetraOp, Teuchos::VERB_HIGH) << std::endl;
567 TEST_ASSERT(nonnull(tpetraOp));
569 const RCP<const VectorSpaceBase<Scalar> > rangeSpace =
570 Thyra::createVectorSpace<Scalar>(tpetraOp->getRangeMap());
571 const RCP<const VectorSpaceBase<Scalar> > domainSpace =
572 Thyra::createVectorSpace<Scalar>(tpetraOp->getDomainMap());
573 const RCP<const LinearOpBase<Scalar> > thyraLinearOp =
574 Thyra::tpetraLinearOp(rangeSpace, domainSpace, tpetraOp);
575 TEST_ASSERT(nonnull(thyraLinearOp));
577 out <<
"\nCheck that operator returns the right thing ...\n";
578 const RCP<VectorBase<Scalar> > x = createMember(thyraLinearOp->domain());
579 Thyra::V_S(x.ptr(), ST::one());
580 const RCP<VectorBase<Scalar> > y = createMember(thyraLinearOp->range());
581 Thyra::apply<Scalar>(*thyraLinearOp, Thyra::NOTRANS, *x, y.ptr());
582 const Scalar sum_y = sum(*y);
583 TEST_FLOATING_EQUALITY( sum_y, as<Scalar>(3+1+2*(y->space()->dim()-2)),
586 out <<
"\nCheck the general LinearOp interface ...\n";
587 Thyra::LinearOpTester<Scalar> linearOpTester;
589 linearOpTester.dump_all(
dumpAll);
591 TEST_ASSERT(linearOpTester.check(*thyraLinearOp, Teuchos::inOutArg(out)));
607 const RCP<Tpetra::Operator<Scalar> > tpetraOp =
608 createTriDiagonalTpetraOperator<Scalar>(
g_localDim);
609 out <<
"tpetraOp = " << Teuchos::describe(*tpetraOp, Teuchos::VERB_HIGH) << std::endl;
611 const RCP<const VectorSpaceBase<Scalar> > rangeSpace =
612 Thyra::createVectorSpace<Scalar>(tpetraOp->getRangeMap());
614 const RCP<const VectorSpaceBase<Scalar> > domainSpace =
615 Thyra::createVectorSpace<Scalar>(tpetraOp->getDomainMap());
618 const RCP<LinearOpBase<Scalar> > thyraOp =
620 TEST_EQUALITY(thyraOp->range(), rangeSpace);
621 TEST_EQUALITY(thyraOp->domain(), domainSpace);
622 const RCP<Tpetra::Operator<Scalar> > tpetraOp2 =
623 ConverterT::getTpetraOperator(thyraOp);
624 TEST_EQUALITY(tpetraOp2, tpetraOp);
628 const RCP<LinearOpBase<Scalar> > thyraOp =
630 TEST_INEQUALITY(thyraOp->range(), rangeSpace);
631 TEST_INEQUALITY(thyraOp->domain(), domainSpace);
632 TEST_ASSERT(thyraOp->range()->isCompatible(*rangeSpace));
633 TEST_ASSERT(thyraOp->domain()->isCompatible(*domainSpace));
634 const RCP<Tpetra::Operator<Scalar> > tpetraOp2 =
635 ConverterT::getTpetraOperator(thyraOp);
636 TEST_EQUALITY(tpetraOp2, tpetraOp);
652 const RCP<const Tpetra::Operator<Scalar> > tpetraOp =
653 createTriDiagonalTpetraOperator<Scalar>(
g_localDim);
654 out <<
"tpetraOp = " << Teuchos::describe(*tpetraOp, Teuchos::VERB_HIGH) << std::endl;
656 const RCP<const VectorSpaceBase<Scalar> > rangeSpace =
657 Thyra::createVectorSpace<Scalar>(tpetraOp->getRangeMap());
659 const RCP<const VectorSpaceBase<Scalar> > domainSpace =
660 Thyra::createVectorSpace<Scalar>(tpetraOp->getDomainMap());
663 const RCP<const LinearOpBase<Scalar> > thyraOp =
665 TEST_EQUALITY(thyraOp->range(), rangeSpace);
666 TEST_EQUALITY(thyraOp->domain(), domainSpace);
667 const RCP<const Tpetra::Operator<Scalar> > tpetraOp2 =
668 ConverterT::getConstTpetraOperator(thyraOp);
669 TEST_EQUALITY(tpetraOp2, tpetraOp);
673 const RCP<const LinearOpBase<Scalar> > thyraOp =
675 TEST_INEQUALITY(thyraOp->range(), rangeSpace);
676 TEST_INEQUALITY(thyraOp->domain(), domainSpace);
677 TEST_ASSERT(thyraOp->range()->isCompatible(*rangeSpace));
678 TEST_ASSERT(thyraOp->domain()->isCompatible(*domainSpace));
679 const RCP<const Tpetra::Operator<Scalar> > tpetraOp2 =
680 ConverterT::getConstTpetraOperator(thyraOp);
681 TEST_EQUALITY(tpetraOp2, tpetraOp);
694 Teuchos::RCP<Teuchos::Time> timer = Teuchos::TimeMonitor::lookupCounter(label);
695 TEUCHOS_TEST_FOR_EXCEPTION(timer == null,
697 "lookupAndAssertTimer(): timer \"" << label <<
"\" was not present in Teuchos::TimeMonitor." 698 " Unit test not valid.");
703 #define CHECK_TPETRA_FUNC_CALL_INCREMENT( timerStr, tpetraCode, thyraCode ) \ 705 out << "\nTesting that Thyra calls down to " << timerStr << "\n"; \ 707 const RCP<const Time> timer = lookupAndAssertTimer(timerStr); \ 708 const int countBefore = timer->numCalls(); \ 710 const int countAfter = timer->numCalls(); \ 711 TEST_EQUALITY( countAfter, countBefore+1 ); \ 719 typedef Teuchos::ScalarTraits<Scalar> ST;
720 typedef typename ST::magnitudeType Magnitude;
721 typedef VectorSpaceBase<Scalar> VectorSpace;
722 typedef MultiVectorBase<Scalar> MultiVec;
724 typedef Tpetra::MultiVector<Scalar> TpetraMultiVec;
727 const int numCols = 3;
729 const RCP<const VectorSpace> vs =
732 A = createMembers(vs, numCols),
733 B = createMembers(vs, numCols);
734 const RCP<TpetraMultiVec>
735 tA = TOVE::getTpetraMultiVector(A),
736 tB = TOVE::getTpetraMultiVector(B);
737 Array<Scalar> C(numCols*numCols,ST::zero());
739 Teuchos::Array<Magnitude> avMag(numCols);
740 Teuchos::Array<Scalar> avScal(numCols);
743 "Tpetra::MultiVector::putScalar()",
744 tA->putScalar(ST::zero()),
745 Thyra::assign(A.ptr(), ST::zero())
749 "Tpetra::MultiVector::dot()",
750 tA->dot(*tB, avScal() ),
751 Thyra::norms( *A, avMag() )
755 "Tpetra::MultiVector::dot()",
756 tA->dot(*tB, avScal() ),
757 A->range()->scalarProds(*A, *B, avScal() )
823 #ifdef TPETRA_TEUCHOS_TIME_MONITOR 824 # define TPETRA_TIMER_TESTS(SCALAR) \ 825 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, UseTpetraImplementations, SCALAR ) 827 # define TPETRA_TIMER_TESTS(SCALAR) 836 #ifdef HAVE_THYRA_TPETRA_EPETRA 844 using Teuchos::outArg;
845 using Teuchos::rcp_dynamic_cast;
846 using Teuchos::Array;
847 typedef Teuchos::ScalarTraits<Scalar> ST;
849 const RCP<Tpetra::Operator<Scalar> > tpetraOp =
850 createTriDiagonalTpetraOperator<Scalar>(
g_localDim);
852 const RCP<LinearOpBase<Scalar> > thyraOp =
855 const RCP<Thyra::TpetraLinearOp<Scalar> > thyraTpetraOp =
858 RCP<const Epetra_Operator> epetraOp;
859 Thyra::EOpTransp epetraOpTransp;
863 thyraTpetraOp->getEpetraOpView( outArg(epetraOp), outArg(epetraOpTransp),
864 outArg(epetraOpApplyAs), outArg(epetraOpAdjointSupport) );
866 if (
typeid(Scalar) ==
typeid(
double)) {
867 TEST_ASSERT(nonnull(epetraOp));
868 const RCP<const Epetra_RowMatrix> epetraRowMatrix =
869 rcp_dynamic_cast<
const Epetra_RowMatrix>(epetraOp,
true);
870 int numRowEntries = -1;
871 epetraRowMatrix->NumMyRowEntries(1, numRowEntries);
872 TEST_EQUALITY_CONST(numRowEntries, 3);
873 Array<double> row_values(numRowEntries);
874 Array<int> row_indices(numRowEntries);
875 epetraRowMatrix->ExtractMyRowCopy(1, numRowEntries, numRowEntries,
876 row_values.getRawPtr(), row_indices.getRawPtr());
877 TEST_EQUALITY_CONST(row_values[0], -1.0);
878 TEST_EQUALITY_CONST(row_values[1], 2.0);
879 TEST_EQUALITY_CONST(row_values[2], 1.0);
883 TEST_ASSERT(is_null(epetraOp));
900 const RCP<Tpetra::Operator<Scalar> > tpetraOp =
901 createTriDiagonalTpetraOperator<Scalar>(
g_localDim);
902 out <<
"tpetraOp = " << Teuchos::describe(*tpetraOp, Teuchos::VERB_HIGH) << std::endl;
903 TEST_ASSERT(nonnull(tpetraOp));
905 const RCP<Tpetra::CrsMatrix<Scalar> > tpetraCrsMatrix =
906 Teuchos::rcp_dynamic_cast<Tpetra::CrsMatrix<Scalar> >(tpetraOp,
true);
908 const RCP<const VectorSpaceBase<Scalar> > rangeSpace =
909 Thyra::createVectorSpace<Scalar>(tpetraOp->getRangeMap());
910 const RCP<const VectorSpaceBase<Scalar> > domainSpace =
911 Thyra::createVectorSpace<Scalar>(tpetraOp->getDomainMap());
912 const RCP<LinearOpBase<Scalar> > thyraLinearOp =
913 Thyra::tpetraLinearOp(rangeSpace, domainSpace, tpetraOp);
914 TEST_ASSERT(nonnull(thyraLinearOp));
916 const Teuchos::RCP<Thyra::RowStatLinearOpBase<Scalar> > rowStatOp =
917 Teuchos::rcp_dynamic_cast<Thyra::RowStatLinearOpBase<Scalar> >(thyraLinearOp,
true);
921 const RCP<VectorBase<Scalar> > inv_row_sums =
922 createMember<Scalar>(thyraLinearOp->range());
923 const RCP<VectorBase<Scalar> > row_sums =
924 createMember<Scalar>(thyraLinearOp->range());
926 rowStatOp->getRowStat(Thyra::RowStatLinearOpBaseUtils::ROW_STAT_INV_ROW_SUM,
928 rowStatOp->getRowStat(Thyra::RowStatLinearOpBaseUtils::ROW_STAT_ROW_SUM,
931 out <<
"inv_row_sums = " << *inv_row_sums;
932 out <<
"row_sums = " << *row_sums;
934 TEST_FLOATING_EQUALITY(
935 Thyra::sum<Scalar>(*row_sums),
936 Teuchos::as<Scalar>(4.0 * thyraLinearOp->domain()->dim() - 2.0),
937 Teuchos::as<Scalar>(10.0 * Teuchos::ScalarTraits<Scalar>::eps())
940 TEST_FLOATING_EQUALITY(
941 Thyra::sum<Scalar>(*inv_row_sums),
942 Teuchos::as<Scalar>( 1.0 / 4.0 * (thyraLinearOp->domain()->dim() - 2) + 2.0 / 3.0 ),
943 Teuchos::as<Scalar>(10.0 * Teuchos::ScalarTraits<Scalar>::eps())
950 const RCP<Tpetra::Operator<Scalar> > tpetraOp =
951 createTriDiagonalTpetraOperator<Scalar>(
g_localDim);
952 out <<
"tpetraOp = " << Teuchos::describe(*tpetraOp, Teuchos::VERB_HIGH) << std::endl;
953 TEST_ASSERT(nonnull(tpetraOp));
955 const RCP<Tpetra::CrsMatrix<Scalar> > tpetraCrsMatrix =
956 Teuchos::rcp_dynamic_cast<Tpetra::CrsMatrix<Scalar> >(tpetraOp,
true);
958 const RCP<const VectorSpaceBase<Scalar> > rangeSpace =
959 Thyra::createVectorSpace<Scalar>(tpetraOp->getRangeMap());
960 const RCP<const VectorSpaceBase<Scalar> > domainSpace =
961 Thyra::createVectorSpace<Scalar>(tpetraOp->getDomainMap());
962 const RCP<LinearOpBase<Scalar> > thyraLinearOp =
963 Thyra::tpetraLinearOp(rangeSpace, domainSpace, tpetraOp);
964 TEST_ASSERT(nonnull(thyraLinearOp));
966 const Teuchos::RCP<Thyra::RowStatLinearOpBase<Scalar> > rowStatOp =
967 Teuchos::rcp_dynamic_cast<Thyra::RowStatLinearOpBase<Scalar> >(thyraLinearOp,
true);
971 const RCP<VectorBase<Scalar> > inv_row_sums =
972 createMember<Scalar>(thyraLinearOp->range());
973 const RCP<VectorBase<Scalar> > row_sums =
974 createMember<Scalar>(thyraLinearOp->range());
976 rowStatOp->getRowStat(Thyra::RowStatLinearOpBaseUtils::ROW_STAT_INV_ROW_SUM,
978 rowStatOp->getRowStat(Thyra::RowStatLinearOpBaseUtils::ROW_STAT_ROW_SUM,
981 out <<
"inv_row_sums = " << *inv_row_sums;
982 out <<
"row_sums = " << *row_sums;
984 const Teuchos::RCP<Thyra::ScaledLinearOpBase<Scalar> > scaledOp =
985 Teuchos::rcp_dynamic_cast<Thyra::ScaledLinearOpBase<Scalar> >(thyraLinearOp,
true);
987 TEUCHOS_ASSERT(scaledOp->supportsScaleLeft());
989 scaledOp->scaleLeft(*inv_row_sums);
991 rowStatOp->getRowStat(Thyra::RowStatLinearOpBaseUtils::ROW_STAT_ROW_SUM,
994 out <<
"row_sums after left scaling by inv_row_sum = " << *row_sums;
997 TEST_FLOATING_EQUALITY(
998 Scalar(row_sums->space()->dim()),
999 Thyra::sum<Scalar>(*row_sums),
1000 as<Scalar>(10.0 * Teuchos::ScalarTraits<Scalar>::eps())
1006 TEUCHOS_ASSERT(scaledOp->supportsScaleRight());
1007 scaledOp->scaleRight(*inv_row_sums);
1008 rowStatOp->getRowStat(Thyra::RowStatLinearOpBaseUtils::ROW_STAT_ROW_SUM,row_sums.ptr());
1009 out <<
"row_sums after right scaling by inv_row_sum = " << *row_sums;
1012 #endif // HAVE_THYRA_TPETRA_EPETRA 1019 #define THYRA_TPETRA_THYRA_WRAPPERS_INSTANT(SCALAR) \ 1021 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1022 convertTpetraToThyraComm, SCALAR ) \ 1024 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1025 createVectorSpace, SCALAR ) \ 1027 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1028 createVector, SCALAR ) \ 1030 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1031 createConstVector, SCALAR ) \ 1033 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1034 createMultiVector, SCALAR ) \ 1036 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1037 createConstMultiVector, SCALAR ) \ 1039 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1040 TeptraVectorSpace, SCALAR ) \ 1042 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1043 vectorSpaceTester, SCALAR ) \ 1045 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1046 vectorStdOpsTester, SCALAR ) \ 1048 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1049 multiVectorStdOpsTester, SCALAR ) \ 1051 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1052 getTpetraMultiVector, SCALAR ) \ 1054 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1055 getConstTpetraMultiVector, SCALAR ) \ 1057 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1058 TpetraLinearOp, SCALAR ) \ 1060 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1061 createLinearOp, SCALAR ) \ 1063 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1064 createConstLinearOp, SCALAR ) \ 1066 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1067 TpetraLinearOp_EpetraRowMatrix, SCALAR ) \ 1069 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1070 TpetraLinearOp_RowStatLinearOpBase, SCALAR ) \ 1072 TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \ 1073 TpetraLinearOp_ScaledLinearOpBase, SCALAR ) \ RCP< MultiVectorBase< Scalar > > createMultiVector(const RCP< Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraMultiVector, const RCP< const VectorSpaceBase< Scalar > > rangeSpace=Teuchos::null, const RCP< const VectorSpaceBase< Scalar > > domainSpace=Teuchos::null)
#define THYRA_TPETRA_THYRA_WRAPPERS_INSTANT(SCALAR)
RCP< Tpetra::Operator< Scalar > > createTriDiagonalTpetraOperator(const int numLocalRows)
#define CHECK_TPETRA_FUNC_CALL_INCREMENT(timerStr, tpetraCode, thyraCode)
RCP< const VectorBase< Scalar > > createConstVector(const RCP< const Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVector, const RCP< const VectorSpaceBase< Scalar > > space=Teuchos::null)
RCP< const LinearOpBase< Scalar > > createConstLinearOp(const RCP< const Tpetra::Operator< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraOperator, const RCP< const VectorSpaceBase< Scalar > > rangeSpace=Teuchos::null, const RCP< const VectorSpaceBase< Scalar > > domainSpace=Teuchos::null)
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)
RCP< VectorBase< Scalar > > createVector(const RCP< Tpetra::Vector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraVector, const RCP< const VectorSpaceBase< Scalar > > space=Teuchos::null)
EAdjointEpetraOp
Determine if adjoints are supported on Epetra_Opeator or not.
RCP< const MultiVectorBase< Scalar > > createConstMultiVector(const RCP< const Tpetra::MultiVector< Scalar, LocalOrdinal, GlobalOrdinal, Node > > &tpetraMultiVector, const RCP< const VectorSpaceBase< Scalar > > rangeSpace=Teuchos::null, const RCP< const VectorSpaceBase< Scalar > > domainSpace=Teuchos::null)
Teuchos::RCP< Teuchos::Time > lookupAndAssertTimer(const std::string &label)
TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL(TpetraThyraWrappers, convertTpetraToThyraComm, Scalar)
EApplyEpetraOpAs
Determine how the apply an Epetra_Operator as a linear operator.
RCP< const VectorSpaceBase< Scalar > > createVectorSpace(const RCP< const Tpetra::Map< LocalOrdinal, GlobalOrdinal, Node > > &tpetraMap)
Create a Thyra::VectorSpaceBase object given a Tpetra::Map.
RCP< const Teuchos::Comm< Ordinal > > convertTpetraToThyraComm(const RCP< const Teuchos::Comm< int > > &tpetraComm)
Given an Tpetra Teuchos::Comm<int> object, return an equivalent Teuchos::Comm<Ordinal> object...
Concrete Thyra::LinearOpBase subclass for Tpetra::Operator.
RCP< const TpetraMap_t > createTpetraMap(const int localDim)
RCP< const VectorSpaceBase< Scalar > > createTpetraVectorSpace(const int localDim)