Thyra Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
TpetraThyraWrappers_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 #include "Thyra_VectorSpaceTester.hpp"
12 #include "Thyra_VectorStdOpsTester.hpp"
13 #include "Thyra_MultiVectorStdOpsTester.hpp"
14 #include "Thyra_VectorStdOps.hpp"
15 #include "Thyra_MultiVectorStdOps.hpp"
16 #include "Thyra_LinearOpTester.hpp"
17 #include "Thyra_DefaultProductVector.hpp"
18 #include "Thyra_TestingTools.hpp"
19 #include "Thyra_ScaledLinearOpBase.hpp"
20 #include "Thyra_RowStatLinearOpBase.hpp"
21 #include "Thyra_VectorStdOps.hpp"
22 #include "Tpetra_CrsMatrix.hpp"
24 #include "Teuchos_DefaultComm.hpp"
25 #include "Teuchos_Tuple.hpp"
26 
27 
28 namespace Thyra {
29 
30 
31 //
32 // Helper code and declarations
33 //
34 
35 
36 using Teuchos::as;
37 using Teuchos::null;
38 using Teuchos::RCP;
39 using Teuchos::rcp;
40 using Teuchos::ArrayView;
41 using Teuchos::rcp_dynamic_cast;
42 using Teuchos::inOutArg;
43 using Teuchos::Comm;
44 using Teuchos::tuple;
45 
46 
47 const int g_localDim = 4; // ToDo: Make variable!
48 
49 
50 typedef Tpetra::Map<> TpetraMap_t;
51 
52 
54 createTpetraMap(const int localDim)
55 {
57  return Teuchos::rcp(new TpetraMap_t(OT::invalid(), localDim, 0,
59  // ToDo: Pass in the comm?
60 }
61 
62 // ToDo: Above: Vary the LocalOrdinal and GlobalOrdinal types?
63 
64 
65 template<class Scalar>
66 RCP<const VectorSpaceBase<Scalar> >
67 createTpetraVectorSpace(const int localDim)
68 {
69  return Thyra::createVectorSpace<Scalar>(createTpetraMap(g_localDim));
70 }
71 
72 
73 template<class Scalar>
74 RCP<Tpetra::Operator<Scalar> >
75 createTriDiagonalTpetraOperator(const int numLocalRows)
76 {
77  typedef Tpetra::global_size_t global_size_t;
78  typedef Tpetra::Map<>::global_ordinal_type GO;
79 
80  RCP<const Tpetra::Map<> > map = createTpetraMap(numLocalRows);
81 
82  const size_t numMyElements = map->getLocalNumElements();
83  const global_size_t numGlobalElements = map->getGlobalNumElements();
84 
85  ArrayView<const GO> myGlobalElements = map->getLocalElementList();
86 
87  // Create an OTeger vector numNz that is used to build the Petra Matrix.
88  // numNz[i] is the Number of OFF-DIAGONAL term for the ith global equation
89  // on this processor
90 
91  Teuchos::Array<size_t> numNz (numMyElements);
92 
93  // We are building a tridiagonal matrix where each row has (-1 2 -1)
94  // So we need 2 off-diagonal terms (except for the first and last equation)
95 
96  for (size_t i=0; i < numMyElements; ++i) {
97  if (myGlobalElements[i] == 0 || static_cast<global_size_t>(myGlobalElements[i]) == numGlobalElements-1) {
98  // boundary
99  numNz[i] = 2;
100  }
101  else {
102  numNz[i] = 3;
103  }
104  }
105 
106  // Create a Tpetra::Matrix using the Map, with a static allocation dictated by numNz
108  Teuchos::rcp( new Tpetra::CrsMatrix<Scalar>(map, numNz ()) );
109 
110  // We are done with NumNZ
111  {
113  swap (empty, numNz); // classic idiom for freeing a container
114  }
115 
116  // Add rows one-at-a-time
117  // Off diagonal values will always be -1
118  const Scalar two = static_cast<Scalar>( 2.0);
119  const Scalar posOne = static_cast<Scalar>(+1.0);
120  const Scalar negOne = static_cast<Scalar>(-1.0);
121  for (size_t i = 0; i < numMyElements; i++) {
122  if (myGlobalElements[i] == 0) {
123  A->insertGlobalValues( myGlobalElements[i],
124  tuple<GO>(myGlobalElements[i], myGlobalElements[i]+1)(),
125  tuple<Scalar> (two, posOne)()
126  );
127  }
128  else if (static_cast<global_size_t>(myGlobalElements[i]) == numGlobalElements-1) {
129  A->insertGlobalValues( myGlobalElements[i],
130  tuple<GO>(myGlobalElements[i]-1, myGlobalElements[i])(),
131  tuple<Scalar> (negOne, two)()
132  );
133  }
134  else {
135  A->insertGlobalValues( myGlobalElements[i],
136  tuple<GO>(myGlobalElements[i]-1, myGlobalElements[i], myGlobalElements[i]+1)(),
137  tuple<Scalar> (negOne, two, posOne)()
138  );
139  }
140  }
141 
142  // Finish up
143  A->fillComplete();
144 
145  return A;
146 
147 }
148 
149 
150 bool showAllTests = false;
151 bool dumpAll = false;
152 bool runLinearOpTester = true;
153 
154 
156 {
158  "show-all-tests", "no-show-all-tests", &showAllTests, "Show all tests or not" );
160  "dump-all", "no-dump-all", &dumpAll, "Dump all objects being tested" );
162  "run-linear-op-tester", "no-run-linear-op-tester", &runLinearOpTester, "..." );
163 }
164 
165 
166 //
167 // Unit Tests
168 //
169 
170 
171 //
172 // convertTpetraToThyraComm
173 //
174 
176  Scalar )
177 {
180  TEST_ASSERT(nonnull(thyraComm));
181 }
182 
183 
184 //
185 // createVectorSpace
186 //
187 
189  Scalar )
190 {
193  Thyra::createVectorSpace<Scalar>(tpetraMap);
194  TEST_ASSERT(nonnull(vs));
195  out << "vs = " << *vs;
196  const RCP<const SpmdVectorSpaceBase<Scalar> > vs_spmd =
197  rcp_dynamic_cast<const SpmdVectorSpaceBase<Scalar> >(vs, true);
198  TEST_EQUALITY(vs_spmd->localSubDim(), g_localDim);
199  TEST_EQUALITY(vs->dim(), as<Ordinal>(tpetraMap->getGlobalNumElements()));
200 
202  RCP<const TpetraMap_t> tpetraMap2 = ConverterT::getTpetraMap(vs);
203  TEST_EQUALITY(tpetraMap2, tpetraMap);
204 }
205 
206 
207 //
208 // createVector
209 //
210 
212  Scalar )
213 {
214 
216 
219  Thyra::createVectorSpace<Scalar>(tpetraMap);
220 
221  const RCP<Tpetra::Vector<Scalar> > tpetraVector =
222  rcp(new Tpetra::Vector<Scalar>(tpetraMap));
223 
224  {
225  const RCP<VectorBase<Scalar> > thyraVector = createVector(tpetraVector, vs);
226  TEST_EQUALITY(thyraVector->space(), vs);
227  const RCP<Tpetra::Vector<Scalar> > tpetraVector2 =
228  ConverterT::getTpetraVector(thyraVector);
229  TEST_EQUALITY(tpetraVector2, tpetraVector);
230  }
231 
232  {
233  const RCP<VectorBase<Scalar> > thyraVector = Thyra::createVector(tpetraVector);
234  TEST_INEQUALITY(thyraVector->space(), vs);
235  TEST_ASSERT(thyraVector->space()->isCompatible(*vs));
236  const RCP<Tpetra::Vector<Scalar> > tpetraVector2 =
237  ConverterT::getTpetraVector(thyraVector);
238  TEST_EQUALITY(tpetraVector2, tpetraVector);
239  }
240 
241 }
242 
243 
244 //
245 // createConstVector
246 //
247 
249  Scalar )
250 {
251 
253 
256  Thyra::createVectorSpace<Scalar>(tpetraMap);
257 
258  const RCP<const Tpetra::Vector<Scalar> > tpetraVector =
259  rcp(new Tpetra::Vector<Scalar>(tpetraMap));
260 
261  {
262  const RCP<const VectorBase<Scalar> > thyraVector =
263  createConstVector(tpetraVector, vs);
264  TEST_EQUALITY(thyraVector->space(), vs);
265  const RCP<const Tpetra::Vector<Scalar> > tpetraVector2 =
266  ConverterT::getConstTpetraVector(thyraVector);
267  TEST_EQUALITY(tpetraVector2, tpetraVector);
268  }
269 
270  {
271  const RCP<const VectorBase<Scalar> > thyraVector =
272  Thyra::createConstVector(tpetraVector);
273  TEST_INEQUALITY(thyraVector->space(), vs);
274  TEST_ASSERT(thyraVector->space()->isCompatible(*vs));
275  const RCP<const Tpetra::Vector<Scalar> > tpetraVector2 =
276  ConverterT::getConstTpetraVector(thyraVector);
277  TEST_EQUALITY(tpetraVector2, tpetraVector);
278  }
279 
280 }
281 
282 
283 //
284 // createMultiVector
285 //
286 
288  Scalar )
289 {
290  typedef Tpetra::Map<>::local_ordinal_type LO;
291  typedef Tpetra::Map<>::global_ordinal_type GO;
292  typedef Tpetra::Map<>::node_type NODE;
294 
295  const int numCols = 3;
296 
298  const RCP<const VectorSpaceBase<Scalar> > rangeVs =
299  Thyra::createVectorSpace<Scalar>(tpetraMap);
300 
301  const RCP<const TpetraMap_t> tpetraLocRepMap =
302  Tpetra::createLocalMapWithNode<LO,GO,NODE>(
303  numCols, tpetraMap->getComm());
304  const RCP<const VectorSpaceBase<Scalar> > domainVs =
305  Thyra::createVectorSpace<Scalar>(tpetraLocRepMap);
306 
307  const RCP<Tpetra::MultiVector<Scalar> > tpetraMultiVector =
308  rcp(new Tpetra::MultiVector<Scalar>(tpetraMap, numCols));
309 
310  {
311  const RCP<MultiVectorBase<Scalar> > thyraMultiVector =
312  createMultiVector(tpetraMultiVector, rangeVs, domainVs);
313  TEST_EQUALITY(thyraMultiVector->range(), rangeVs);
314  TEST_EQUALITY(thyraMultiVector->domain(), domainVs);
315  const RCP<Tpetra::MultiVector<Scalar> > tpetraMultiVector2 =
316  ConverterT::getTpetraMultiVector(thyraMultiVector);
317  TEST_EQUALITY(tpetraMultiVector2, tpetraMultiVector);
318  }
319 
320  {
321  const RCP<MultiVectorBase<Scalar> > thyraMultiVector =
322  Thyra::createMultiVector(tpetraMultiVector);
323  TEST_INEQUALITY(thyraMultiVector->range(), rangeVs);
324  TEST_INEQUALITY(thyraMultiVector->domain(), domainVs);
325  TEST_ASSERT(thyraMultiVector->range()->isCompatible(*rangeVs));
326  TEST_ASSERT(thyraMultiVector->domain()->isCompatible(*domainVs));
327  const RCP<Tpetra::MultiVector<Scalar> > tpetraMultiVector2 =
328  ConverterT::getTpetraMultiVector(thyraMultiVector);
329  TEST_EQUALITY(tpetraMultiVector2, tpetraMultiVector);
330  }
331 
332 }
333 
334 
335 //
336 // createConstMultiVector
337 //
338 
340  Scalar )
341 {
342  typedef Tpetra::Map<>::local_ordinal_type LO;
343  typedef Tpetra::Map<>::global_ordinal_type GO;
344  typedef Tpetra::Map<>::node_type NODE;
346 
347  const int numCols = 3;
348 
350  const RCP<const VectorSpaceBase<Scalar> > rangeVs =
351  Thyra::createVectorSpace<Scalar>(tpetraMap);
352 
353  const RCP<const TpetraMap_t> tpetraLocRepMap =
354  Tpetra::createLocalMapWithNode<LO,GO,NODE>(
355  numCols, tpetraMap->getComm());
356  const RCP<const VectorSpaceBase<Scalar> > domainVs =
357  Thyra::createVectorSpace<Scalar>(tpetraLocRepMap);
358 
359  const RCP<const Tpetra::MultiVector<Scalar> > tpetraMultiVector =
360  rcp(new Tpetra::MultiVector<Scalar>(tpetraMap, numCols));
361 
362  {
363  const RCP<const MultiVectorBase<Scalar> > thyraMultiVector =
364  createConstMultiVector(tpetraMultiVector, rangeVs, domainVs);
365  TEST_EQUALITY(thyraMultiVector->range(), rangeVs);
366  TEST_EQUALITY(thyraMultiVector->domain(), domainVs);
367  const RCP<const Tpetra::MultiVector<Scalar> > tpetraMultiVector2 =
368  ConverterT::getConstTpetraMultiVector(thyraMultiVector);
369  TEST_EQUALITY(tpetraMultiVector2, tpetraMultiVector);
370  }
371 
372  {
373  const RCP<const MultiVectorBase<Scalar> > thyraMultiVector =
374  Thyra::createConstMultiVector(tpetraMultiVector);
375  TEST_INEQUALITY(thyraMultiVector->range(), rangeVs);
376  TEST_INEQUALITY(thyraMultiVector->domain(), domainVs);
377  TEST_ASSERT(thyraMultiVector->range()->isCompatible(*rangeVs));
378  TEST_ASSERT(thyraMultiVector->domain()->isCompatible(*domainVs));
379  const RCP<const Tpetra::MultiVector<Scalar> > tpetraMultiVector2 =
380  ConverterT::getConstTpetraMultiVector(thyraMultiVector);
381  TEST_EQUALITY(tpetraMultiVector2, tpetraMultiVector);
382  }
383 
384 }
385 
386 
387 //
388 // TeptraVectorSpace
389 //
390 
391 
392 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, TeptraVectorSpace,
393  Scalar )
394 {
396  Thyra::createVectorSpace<Scalar>(createTpetraMap(g_localDim));
397  const RCP<VectorBase<Scalar> > v = createMember(vs);
398  TEST_ASSERT(nonnull(v));
399  TEST_EQUALITY(v->space(), vs);
400 }
401 
402 
403 //
404 // vectorSpaceTester
405 //
406 
407 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, vectorSpaceTester,
408  Scalar )
409 {
411  = createTpetraVectorSpace<Scalar>(g_localDim);
412  Thyra::VectorSpaceTester<Scalar> vectorSpaceTester;
413  vectorSpaceTester.show_all_tests(showAllTests);
414  vectorSpaceTester.dump_all(dumpAll);
415  TEST_ASSERT(vectorSpaceTester.check(*vs, &out));
416 }
417 
418 
419 //
420 // vectorStdOpsTester
421 //
422 
423 
424 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, vectorStdOpsTester,
425  Scalar )
426 {
428  Thyra::createVectorSpace<Scalar>(createTpetraMap(g_localDim));
429  Thyra::VectorStdOpsTester<Scalar> vectorStdOpsTester;
430  vectorStdOpsTester.warning_tol(5.0e-13);
431  vectorStdOpsTester.error_tol(5.0e-14);
432  TEST_ASSERT(vectorStdOpsTester.checkStdOps(*vs, &out));
433 }
434 
435 
436 //
437 // multiVectorStdOpsTester
438 //
439 
440 
441 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, multiVectorStdOpsTester,
442  Scalar )
443 {
445  Thyra::createVectorSpace<Scalar>(createTpetraMap(g_localDim));
446  Thyra::MultiVectorStdOpsTester<Scalar> mvStdOpsTester;
447  mvStdOpsTester.warning_tol(5.0e-13);
448  mvStdOpsTester.error_tol(5.0e-14);
449  TEST_ASSERT(mvStdOpsTester.checkStdOps(*vs, &out));
450 }
451 
452 
453 //
454 // getTpetraMultiVector
455 //
456 
457 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, getTpetraMultiVector,
458  Scalar )
459 {
461 
462  const int numCols = 3;
464  = createTpetraVectorSpace<Scalar>(g_localDim);
465 
466  {
467  const RCP<MultiVectorBase<Scalar> > mv = createMembers(vs, numCols);
468  const RCP<Tpetra::MultiVector<Scalar> > tmv =
469  ConverterT::getTpetraMultiVector(mv);
470  TEST_ASSERT(nonnull(tmv));
471  TEST_EQUALITY(as<Ordinal>(tmv->getMap()->getGlobalNumElements()), vs->dim());
472  }
473 
474  {
475  const RCP<VectorBase<Scalar> > v = createMember(vs);
476  const RCP<Tpetra::MultiVector<Scalar> > tmv =
477  ConverterT::getTpetraMultiVector(v);
478  TEST_ASSERT(nonnull(tmv));
479  TEST_EQUALITY(as<Ordinal>(tmv->getMap()->getGlobalNumElements()), vs->dim());
480  }
481 
482 #ifdef THYRA_DEBUG
483  const RCP<VectorBase<Scalar> > pv = Thyra::defaultProductVector<Scalar>();
484  TEST_THROW(ConverterT::getTpetraMultiVector(pv), std::logic_error);
485 #endif
486 
487 }
488 
489 
490 //
491 // getConstTpetraMultiVector
492 //
493 
494 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, getConstTpetraMultiVector,
495  Scalar )
496 {
498 
499  const int numCols = 3;
501  = createTpetraVectorSpace<Scalar>(g_localDim);
502 
503  {
504  const RCP<const MultiVectorBase<Scalar> > mv = createMembers(vs, numCols);
506  ConverterT::getConstTpetraMultiVector(mv);
507  TEST_ASSERT(nonnull(tmv));
508  TEST_EQUALITY(as<Ordinal>(tmv->getMap()->getGlobalNumElements()), vs->dim());
509  }
510 
511  {
512  const RCP<const VectorBase<Scalar> > v = createMember(vs);
514  ConverterT::getConstTpetraMultiVector(v);
515  TEST_ASSERT(nonnull(tmv));
516  TEST_EQUALITY(as<Ordinal>(tmv->getMap()->getGlobalNumElements()), vs->dim());
517  }
518 
519 #ifdef THYRA_DEBUG
520  const RCP<const VectorBase<Scalar> > pv = Thyra::defaultProductVector<Scalar>();
521  TEST_THROW(ConverterT::getConstTpetraMultiVector(pv), std::logic_error);
522 #endif
523 
524 }
525 
526 
527 //
528 // TpetraLinearOp
529 //
530 
532  Scalar )
533 {
534 
536  using Teuchos::as;
537 
538  const RCP<Tpetra::Operator<Scalar> > tpetraOp =
539  createTriDiagonalTpetraOperator<Scalar>(g_localDim);
540  out << "tpetraOp = " << Teuchos::describe(*tpetraOp, Teuchos::VERB_HIGH) << std::endl;
541  TEST_ASSERT(nonnull(tpetraOp));
542 
543  const RCP<const VectorSpaceBase<Scalar> > rangeSpace =
544  Thyra::createVectorSpace<Scalar>(tpetraOp->getRangeMap());
545  const RCP<const VectorSpaceBase<Scalar> > domainSpace =
546  Thyra::createVectorSpace<Scalar>(tpetraOp->getDomainMap());
547  const RCP<const LinearOpBase<Scalar> > thyraLinearOp =
548  Thyra::tpetraLinearOp(rangeSpace, domainSpace, tpetraOp);
549  TEST_ASSERT(nonnull(thyraLinearOp));
550 
551  out << "\nCheck that operator returns the right thing ...\n";
552  const RCP<VectorBase<Scalar> > x = createMember(thyraLinearOp->domain());
553  Thyra::V_S(x.ptr(), ST::one());
554  const RCP<VectorBase<Scalar> > y = createMember(thyraLinearOp->range());
555  Thyra::apply<Scalar>(*thyraLinearOp, Thyra::NOTRANS, *x, y.ptr());
556  const Scalar sum_y = sum(*y);
557  TEST_FLOATING_EQUALITY( sum_y, as<Scalar>(3+1+2*(y->space()->dim()-2)),
558  100.0 * ST::eps() );
559 
560  out << "\nCheck the general LinearOp interface ...\n";
561  Thyra::LinearOpTester<Scalar> linearOpTester;
562  linearOpTester.show_all_tests(showAllTests);
563  linearOpTester.dump_all(dumpAll);
564  if (runLinearOpTester) {
565  TEST_ASSERT(linearOpTester.check(*thyraLinearOp, Teuchos::inOutArg(out)));
566  }
567 
568 }
569 
570 
571 //
572 // createLinearOp
573 //
574 
576  Scalar )
577 {
578 
580 
581  const RCP<Tpetra::Operator<Scalar> > tpetraOp =
582  createTriDiagonalTpetraOperator<Scalar>(g_localDim);
583  out << "tpetraOp = " << Teuchos::describe(*tpetraOp, Teuchos::VERB_HIGH) << std::endl;
584 
585  const RCP<const VectorSpaceBase<Scalar> > rangeSpace =
586  Thyra::createVectorSpace<Scalar>(tpetraOp->getRangeMap());
587 
588  const RCP<const VectorSpaceBase<Scalar> > domainSpace =
589  Thyra::createVectorSpace<Scalar>(tpetraOp->getDomainMap());
590 
591  {
592  const RCP<LinearOpBase<Scalar> > thyraOp =
593  createLinearOp(tpetraOp, rangeSpace, domainSpace);
594  TEST_EQUALITY(thyraOp->range(), rangeSpace);
595  TEST_EQUALITY(thyraOp->domain(), domainSpace);
596  const RCP<Tpetra::Operator<Scalar> > tpetraOp2 =
597  ConverterT::getTpetraOperator(thyraOp);
598  TEST_EQUALITY(tpetraOp2, tpetraOp);
599  }
600 
601  {
602  const RCP<LinearOpBase<Scalar> > thyraOp =
603  Thyra::createLinearOp(tpetraOp);
604  TEST_INEQUALITY(thyraOp->range(), rangeSpace);
605  TEST_INEQUALITY(thyraOp->domain(), domainSpace);
606  TEST_ASSERT(thyraOp->range()->isCompatible(*rangeSpace));
607  TEST_ASSERT(thyraOp->domain()->isCompatible(*domainSpace));
608  const RCP<Tpetra::Operator<Scalar> > tpetraOp2 =
609  ConverterT::getTpetraOperator(thyraOp);
610  TEST_EQUALITY(tpetraOp2, tpetraOp);
611  }
612 
613 }
614 
615 
616 //
617 // createConstLinearOp
618 //
619 
621  Scalar )
622 {
623 
625 
626  const RCP<const Tpetra::Operator<Scalar> > tpetraOp =
627  createTriDiagonalTpetraOperator<Scalar>(g_localDim);
628  out << "tpetraOp = " << Teuchos::describe(*tpetraOp, Teuchos::VERB_HIGH) << std::endl;
629 
630  const RCP<const VectorSpaceBase<Scalar> > rangeSpace =
631  Thyra::createVectorSpace<Scalar>(tpetraOp->getRangeMap());
632 
633  const RCP<const VectorSpaceBase<Scalar> > domainSpace =
634  Thyra::createVectorSpace<Scalar>(tpetraOp->getDomainMap());
635 
636  {
637  const RCP<const LinearOpBase<Scalar> > thyraOp =
638  createConstLinearOp(tpetraOp, rangeSpace, domainSpace);
639  TEST_EQUALITY(thyraOp->range(), rangeSpace);
640  TEST_EQUALITY(thyraOp->domain(), domainSpace);
641  const RCP<const Tpetra::Operator<Scalar> > tpetraOp2 =
642  ConverterT::getConstTpetraOperator(thyraOp);
643  TEST_EQUALITY(tpetraOp2, tpetraOp);
644  }
645 
646  {
647  const RCP<const LinearOpBase<Scalar> > thyraOp =
648  Thyra::createConstLinearOp(tpetraOp);
649  TEST_INEQUALITY(thyraOp->range(), rangeSpace);
650  TEST_INEQUALITY(thyraOp->domain(), domainSpace);
651  TEST_ASSERT(thyraOp->range()->isCompatible(*rangeSpace));
652  TEST_ASSERT(thyraOp->domain()->isCompatible(*domainSpace));
653  const RCP<const Tpetra::Operator<Scalar> > tpetraOp2 =
654  ConverterT::getConstTpetraOperator(thyraOp);
655  TEST_EQUALITY(tpetraOp2, tpetraOp);
656  }
657 
658 }
659 
660 
661 //
662 // Tpetra-implemented methods
663 //
664 
665 
667 {
669  TEUCHOS_TEST_FOR_EXCEPTION(timer == null,
670  std::runtime_error,
671  "lookupAndAssertTimer(): timer \"" << label << "\" was not present in Teuchos::TimeMonitor."
672  " Unit test not valid.");
673  return timer;
674 }
675 
676 
677 #define CHECK_TPETRA_FUNC_CALL_INCREMENT( timerStr, tpetraCode, thyraCode ) \
678 { \
679  out << "\nTesting that Thyra calls down to " << timerStr << "\n"; \
680  ECHO(tpetraCode); \
681  const RCP<const Time> timer = lookupAndAssertTimer(timerStr); \
682  const int countBefore = timer->numCalls(); \
683  ECHO(thyraCode); \
684  const int countAfter = timer->numCalls(); \
685  TEST_EQUALITY( countAfter, countBefore+1 ); \
686 }
687 
688 
689 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, UseTpetraImplementations,
690  Scalar )
691 {
692  using Teuchos::Time;
694  typedef typename ST::magnitudeType Magnitude;
695  typedef VectorSpaceBase<Scalar> VectorSpace;
696  typedef MultiVectorBase<Scalar> MultiVec;
697  //typedef Tpetra::Map<int> TpetraMap;
698  typedef Tpetra::MultiVector<Scalar> TpetraMultiVec;
700 
701  const int numCols = 3;
702 
703  const RCP<const VectorSpace> vs =
704  createTpetraVectorSpace<Scalar>(g_localDim);
705  const RCP<MultiVec>
706  A = createMembers(vs, numCols),
707  B = createMembers(vs, numCols);
708  const RCP<TpetraMultiVec>
709  tA = TOVE::getTpetraMultiVector(A),
710  tB = TOVE::getTpetraMultiVector(B);
711  Array<Scalar> C(numCols*numCols,ST::zero());
712 
713  Teuchos::Array<Magnitude> avMag(numCols);
714  Teuchos::Array<Scalar> avScal(numCols);
715 
717  "Tpetra::MultiVector::putScalar()",
718  tA->putScalar(ST::zero()),
719  Thyra::assign(A.ptr(), ST::zero())
720  );
721 
723  "Tpetra::MultiVector::dot()",
724  tA->dot(*tB, avScal() ), // norms calls scalarProd, which calls Tpetra::dot
725  Thyra::norms( *A, avMag() )
726  );
727 
729  "Tpetra::MultiVector::dot()",
730  tA->dot(*tB, avScal() ),
731  A->range()->scalarProds(*A, *B, avScal() )
732  );
733 
734  /*
735  CHECK_TPETRA_FUNC_CALL_INCREMENT(
736  "Tpetra::MultiVector::scale(alpha)",
737  tB->scale( ST::zero() ),
738  Thyra::scale( ST::zero(), B.ptr() )
739  );
740  */
741 
742  /*
743  CHECK_TPETRA_FUNC_CALL_INCREMENT(
744  "Tpetra::MultiVector::operator=()",
745  (*tA) = *tB,
746  Thyra::assign( A.ptr(), *B )
747  );
748  */
749 
750  /*
751  {
752  RCP<MultiVec> Ctmvb = Thyra::createMembersView(
753  A->domain(),
754  RTOpPack::SubMultiVectorView<Scalar>(
755  0, numCols, 0, numCols,
756  Teuchos::arcpFromArrayView(C()), numCols
757  )
758  );
759  CHECK_TPETRA_FUNC_CALL_INCREMENT(
760  "Tpetra::multiplyOntoHost()",
761  Tpetra::multiplyOntoHost(Teuchos::CONJ_TRANS,Teuchos::NO_TRANS,ST::one(),*tA,*tB, C()),
762  A->apply(Thyra::CONJTRANS,*B,Ctmvb.ptr(),ST::one(),ST::zero())
763  );
764  }
765  */
766 
767  /*
768  {
769  RCP<const MultiVec> Ctmvb = Thyra::createMembersView(
770  A->domain(),
771  RTOpPack::ConstSubMultiVectorView<Scalar>(
772  0, numCols, 0, numCols,
773  Teuchos::arcpFromArrayView(C()), numCols
774  )
775  );
776  const RCP<const TpetraMultiVec>
777  tC = TOVE::getConstTpetraMultiVector(Ctmvb);
778  CHECK_TPETRA_FUNC_CALL_INCREMENT(
779  "Tpetra::MultiVector::multiply()",
780  tB->multiply(Teuchos::NO_TRANS,Teuchos::NO_TRANS,ST::one(),*tA,*tC,ST::zero()),
781  A->apply(Thyra::NOTRANS,*Ctmvb,B.ptr(),ST::one(),ST::zero())
782  );
783  }
784  */
785 
786 /*
787  RCP<Time>
788  timerUpdate = lookupAndAssertTimer("Tpetra::MultiVector::update(alpha,A,beta,B,gamma)"),
789  timerUpdate2 = lookupAndAssertTimer("Tpetra::MultiVector::update(alpha,A,beta)"),
790 
791  // TODO: test update(two vector)
792  // TODO: test update(three vector)
793 */
794 }
795 
796 
797 #ifdef TPETRA_TEUCHOS_TIME_MONITOR
798 # define TPETRA_TIMER_TESTS(SCALAR) \
799  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, UseTpetraImplementations, SCALAR )
800 #else
801 # define TPETRA_TIMER_TESTS(SCALAR)
802 #endif
803 
804 
805 //
806 // TpetraLinearOp_EpetraRowMatrix
807 //
808 
809 
810 #ifdef HAVE_THYRA_TPETRA_EPETRA
811 
812 
813 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, TpetraLinearOp_EpetraRowMatrix,
814  Scalar )
815 {
816 
817  using Teuchos::as;
818  using Teuchos::outArg;
819  using Teuchos::rcp_dynamic_cast;
820  using Teuchos::Array;
822 
823  const RCP<Tpetra::Operator<Scalar> > tpetraOp =
824  createTriDiagonalTpetraOperator<Scalar>(g_localDim);
825 
826  const RCP<LinearOpBase<Scalar> > thyraOp =
827  Thyra::createLinearOp(tpetraOp);
828 
829  const RCP<Thyra::TpetraLinearOp<Scalar> > thyraTpetraOp =
830  Teuchos::rcp_dynamic_cast<Thyra::TpetraLinearOp<Scalar> >(thyraOp);
831 
832  RCP<const Epetra_Operator> epetraOp;
833  Thyra::EOpTransp epetraOpTransp;
834  Thyra::EApplyEpetraOpAs epetraOpApplyAs;
835  Thyra::EAdjointEpetraOp epetraOpAdjointSupport;
836 
837  thyraTpetraOp->getEpetraOpView( outArg(epetraOp), outArg(epetraOpTransp),
838  outArg(epetraOpApplyAs), outArg(epetraOpAdjointSupport) );
839 
840  if (typeid(Scalar) == typeid(double)) {
841  TEST_ASSERT(nonnull(epetraOp));
842  const RCP<const Epetra_RowMatrix> epetraRowMatrix =
843  rcp_dynamic_cast<const Epetra_RowMatrix>(epetraOp, true);
844  int numRowEntries = -1;
845  epetraRowMatrix->NumMyRowEntries(1, numRowEntries);
846  TEST_EQUALITY_CONST(numRowEntries, 3);
847  Array<double> row_values(numRowEntries);
848  Array<int> row_indices(numRowEntries);
849  epetraRowMatrix->ExtractMyRowCopy(1, numRowEntries, numRowEntries,
850  row_values.getRawPtr(), row_indices.getRawPtr());
851  TEST_EQUALITY_CONST(row_values[0], -1.0);
852  TEST_EQUALITY_CONST(row_values[1], 2.0);
853  TEST_EQUALITY_CONST(row_values[2], 1.0);
854  // ToDo: Test column indices!
855  }
856  else {
857  TEST_ASSERT(is_null(epetraOp));
858  }
859 
860 }
861 
862 #else
863 
864 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, TpetraLinearOp_EpetraRowMatrix,
865  Scalar )
866 {
867 }
868 
869 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, TpetraLinearOp_RowStatLinearOpBase,
870  Scalar )
871 {
872  using Teuchos::as;
873 
874  const RCP<Tpetra::Operator<Scalar> > tpetraOp =
875  createTriDiagonalTpetraOperator<Scalar>(g_localDim);
876  out << "tpetraOp = " << Teuchos::describe(*tpetraOp, Teuchos::VERB_HIGH) << std::endl;
877  TEST_ASSERT(nonnull(tpetraOp));
878 
879  const RCP<Tpetra::CrsMatrix<Scalar> > tpetraCrsMatrix =
880  Teuchos::rcp_dynamic_cast<Tpetra::CrsMatrix<Scalar> >(tpetraOp,true);
881 
882  const RCP<const VectorSpaceBase<Scalar> > rangeSpace =
883  Thyra::createVectorSpace<Scalar>(tpetraOp->getRangeMap());
884  const RCP<const VectorSpaceBase<Scalar> > domainSpace =
885  Thyra::createVectorSpace<Scalar>(tpetraOp->getDomainMap());
886  const RCP<LinearOpBase<Scalar> > thyraLinearOp =
887  Thyra::tpetraLinearOp(rangeSpace, domainSpace, tpetraOp);
888  TEST_ASSERT(nonnull(thyraLinearOp));
889 
891  Teuchos::rcp_dynamic_cast<Thyra::RowStatLinearOpBase<Scalar> >(thyraLinearOp, true);
892 
893  // Get the inverse row sums
894 
895  const RCP<VectorBase<Scalar> > inv_row_sums =
896  createMember<Scalar>(thyraLinearOp->range());
897  const RCP<VectorBase<Scalar> > row_sums =
898  createMember<Scalar>(thyraLinearOp->range());
899 
900  rowStatOp->getRowStat(Thyra::RowStatLinearOpBaseUtils::ROW_STAT_INV_ROW_SUM,
901  inv_row_sums.ptr());
902  rowStatOp->getRowStat(Thyra::RowStatLinearOpBaseUtils::ROW_STAT_ROW_SUM,
903  row_sums.ptr());
904 
905  out << "inv_row_sums = " << *inv_row_sums;
906  out << "row_sums = " << *row_sums;
907 
909  Thyra::sum<Scalar>(*row_sums),
910  Teuchos::as<Scalar>(4.0 * thyraLinearOp->domain()->dim() - 2.0),
911  Teuchos::as<Scalar>(10.0 * Teuchos::ScalarTraits<Scalar>::eps())
912  );
913 
915  Thyra::sum<Scalar>(*inv_row_sums),
916  Teuchos::as<Scalar>( 1.0 / 4.0 * (thyraLinearOp->domain()->dim() - 2) + 2.0 / 3.0 ),
917  Teuchos::as<Scalar>(10.0 * Teuchos::ScalarTraits<Scalar>::eps())
918  );
919 }
920 
921 TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL( TpetraThyraWrappers, TpetraLinearOp_ScaledLinearOpBase,
922  Scalar )
923 {
924  const RCP<Tpetra::Operator<Scalar> > tpetraOp =
925  createTriDiagonalTpetraOperator<Scalar>(g_localDim);
926  out << "tpetraOp = " << Teuchos::describe(*tpetraOp, Teuchos::VERB_HIGH) << std::endl;
927  TEST_ASSERT(nonnull(tpetraOp));
928 
929  const RCP<Tpetra::CrsMatrix<Scalar> > tpetraCrsMatrix =
930  Teuchos::rcp_dynamic_cast<Tpetra::CrsMatrix<Scalar> >(tpetraOp,true);
931 
932  const RCP<const VectorSpaceBase<Scalar> > rangeSpace =
933  Thyra::createVectorSpace<Scalar>(tpetraOp->getRangeMap());
934  const RCP<const VectorSpaceBase<Scalar> > domainSpace =
935  Thyra::createVectorSpace<Scalar>(tpetraOp->getDomainMap());
936  const RCP<LinearOpBase<Scalar> > thyraLinearOp =
937  Thyra::tpetraLinearOp(rangeSpace, domainSpace, tpetraOp);
938  TEST_ASSERT(nonnull(thyraLinearOp));
939 
941  Teuchos::rcp_dynamic_cast<Thyra::RowStatLinearOpBase<Scalar> >(thyraLinearOp, true);
942 
943  // Get the inverse row sums
944 
945  const RCP<VectorBase<Scalar> > inv_row_sums =
946  createMember<Scalar>(thyraLinearOp->range());
947  const RCP<VectorBase<Scalar> > row_sums =
948  createMember<Scalar>(thyraLinearOp->range());
949 
950  rowStatOp->getRowStat(Thyra::RowStatLinearOpBaseUtils::ROW_STAT_INV_ROW_SUM,
951  inv_row_sums.ptr());
952  rowStatOp->getRowStat(Thyra::RowStatLinearOpBaseUtils::ROW_STAT_ROW_SUM,
953  row_sums.ptr());
954 
955  out << "inv_row_sums = " << *inv_row_sums;
956  out << "row_sums = " << *row_sums;
957 
959  Teuchos::rcp_dynamic_cast<Thyra::ScaledLinearOpBase<Scalar> >(thyraLinearOp, true);
960 
961  TEUCHOS_ASSERT(scaledOp->supportsScaleLeft());
962 
963  scaledOp->scaleLeft(*inv_row_sums);
964 
965  rowStatOp->getRowStat(Thyra::RowStatLinearOpBaseUtils::ROW_STAT_ROW_SUM,
966  row_sums.ptr());
967 
968  out << "row_sums after left scaling by inv_row_sum = " << *row_sums;
969 
970  // scaled row sums should be one for each entry
972  Scalar(row_sums->space()->dim()),
973  Thyra::sum<Scalar>(*row_sums),
974  as<Scalar>(10.0 * Teuchos::ScalarTraits<Scalar>::eps())
975  );
976 
977  // Don't currently check the results of right scaling. Tpetra tests
978  // already check this. Once column sums are supported in tpetra
979  // adapters, this can be checked easily.
980  TEUCHOS_ASSERT(scaledOp->supportsScaleRight());
981  scaledOp->scaleRight(*inv_row_sums);
982  rowStatOp->getRowStat(Thyra::RowStatLinearOpBaseUtils::ROW_STAT_ROW_SUM,row_sums.ptr());
983  out << "row_sums after right scaling by inv_row_sum = " << *row_sums;
984 }
985 
986 #endif // HAVE_THYRA_TPETRA_EPETRA
987 
988 
989 //
990 // Unit test instantiations
991 //
992 
993 #define THYRA_TPETRA_THYRA_WRAPPERS_INSTANT(SCALAR) \
994  \
995  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
996  convertTpetraToThyraComm, SCALAR ) \
997  \
998  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
999  createVectorSpace, SCALAR ) \
1000  \
1001  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1002  createVector, SCALAR ) \
1003  \
1004  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1005  createConstVector, SCALAR ) \
1006  \
1007  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1008  createMultiVector, SCALAR ) \
1009  \
1010  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1011  createConstMultiVector, SCALAR ) \
1012  \
1013  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1014  TeptraVectorSpace, SCALAR ) \
1015  \
1016  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1017  vectorSpaceTester, SCALAR ) \
1018  \
1019  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1020  vectorStdOpsTester, SCALAR ) \
1021  \
1022  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1023  multiVectorStdOpsTester, SCALAR ) \
1024  \
1025  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1026  getTpetraMultiVector, SCALAR ) \
1027  \
1028  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1029  getConstTpetraMultiVector, SCALAR ) \
1030  \
1031  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1032  TpetraLinearOp, SCALAR ) \
1033  \
1034  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1035  createLinearOp, SCALAR ) \
1036  \
1037  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1038  createConstLinearOp, SCALAR ) \
1039  \
1040  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1041  TpetraLinearOp_EpetraRowMatrix, SCALAR ) \
1042  \
1043  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1044  TpetraLinearOp_RowStatLinearOpBase, SCALAR ) \
1045  \
1046  TEUCHOS_UNIT_TEST_TEMPLATE_1_INSTANT( TpetraThyraWrappers, \
1047  TpetraLinearOp_ScaledLinearOpBase, SCALAR ) \
1048 
1049 
1050 // We can currently only explicitly instantiate with double support because
1051 // Tpetra only supports explicit instantaition with double. As for implicit
1052 // instantation, g++ 3.4.6 on my Linux machine was taking more than 30 minutes
1053 // to compile this file when all of the types double, float, complex<double>,
1054 // and complex<float> where enabled. Therefore, we will only test double for
1055 // now until explicit instantation with other types are supported by Tpetra.
1056 
1058 
1059 
1060 } // namespace Thyra
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)
#define TEST_INEQUALITY(v1, v2)
#define TEST_ASSERT(v1)
bool is_null(const boost::shared_ptr< T > &p)
static magnitudeType eps()
static CommandLineProcessor & getCLP()
#define TEST_EQUALITY_CONST(v1, v2)
RCP< Tpetra::Operator< Scalar > > createTriDiagonalTpetraOperator(const int numLocalRows)
#define CHECK_TPETRA_FUNC_CALL_INCREMENT(timerStr, tpetraCode, thyraCode)
static RCP< Time > lookupCounter(const std::string &name)
#define TEUCHOS_TEST_FOR_EXCEPTION(throw_exception_test, Exception, msg)
static Teuchos::RCP< const Comm< OrdinalType > > getComm()
Traits class that enables the extraction of Tpetra operator/vector objects wrapped in Thyra operator/...
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)
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)
void setOption(const char option_true[], const char option_false[], bool *option_val, const char documentation[]=NULL)
EAdjointEpetraOp
Determine if adjoints are supported on Epetra_Opeator or not.
Ptr< T > ptr() const
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)
TypeTo as(const TypeFrom &t)
Teuchos::RCP< Teuchos::Time > lookupAndAssertTimer(const std::string &label)
#define TEST_FLOATING_EQUALITY(v1, v2, tol)
bool nonnull(const boost::shared_ptr< T > &p)
TEUCHOS_UNIT_TEST_TEMPLATE_1_DECL(TpetraThyraWrappers, convertTpetraToThyraComm, Scalar)
#define TEST_THROW(code, ExceptType)
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.
#define TEST_EQUALITY(v1, v2)
RCP< const Teuchos::Comm< Ordinal > > convertTpetraToThyraComm(const RCP< const Teuchos::Comm< int > > &tpetraComm)
Given an Tpetra Teuchos::Comm&lt;int&gt; object, return an equivalent Teuchos::Comm&lt;Ordinal&gt; object...
#define TEUCHOS_ASSERT(assertion_test)
Concrete Thyra::LinearOpBase subclass for Tpetra::Operator.
RCP< const TpetraMap_t > createTpetraMap(const int localDim)
RCP< const VectorSpaceBase< Scalar > > createTpetraVectorSpace(const int localDim)