11 #include "Thyra_VectorStdOps.hpp"
14 #include "Thyra_DefaultBlockedLinearOp.hpp"
15 #include "Thyra_ProductVectorBase.hpp"
16 #include "Thyra_SpmdVectorSpaceBase.hpp"
17 #include "Thyra_DetachedSpmdVectorView.hpp"
18 #include "Thyra_TestingTools.hpp"
19 #include "Thyra_LinearOpTester.hpp"
20 #include "Trilinos_Util_CrsMatrixGallery.h"
28 # include "Epetra_MpiComm.h"
30 # include "Epetra_SerialComm.h"
32 #include "Epetra_Vector.h"
33 #include "Epetra_CrsMatrix.h"
44 using Teuchos::rcp_dynamic_cast;
45 using Teuchos::inOutArg;
53 Epetra_MpiComm comm(MPI_COMM_WORLD);
55 Epetra_SerialComm comm;
58 out <<
"\nRunning on " << comm.NumProc() <<
" processors\n";
63 out <<
"Using Trilinos_Util to create test matrices\n";
66 Trilinos_Util::CrsMatrixGallery FGallery(
"recirc_2d",comm,
false);
67 FGallery.Set(
"nx",nx);
68 FGallery.Set(
"ny",ny);
71 Trilinos_Util::CrsMatrixGallery CGallery(
"laplace_2d",comm,
false);
72 CGallery.Set(
"nx",nx);
73 CGallery.Set(
"ny",ny);
76 Trilinos_Util::CrsMatrixGallery BGallery(
"diag",comm,
false);
77 BGallery.Set(
"nx",nx*ny);
78 BGallery.Set(
"a",5.0);
81 Trilinos_Util::CrsMatrixGallery BtGallery(
"diag",comm,
false);
82 BtGallery.Set(
"nx",nx*ny);
83 BtGallery.Set(
"a",3.0);
87 out <<
"Building block2x2 Thyra matrix ... wrapping in EpetraOperatorWrapper\n";
89 Thyra::block2x2<double>(
90 Thyra::epetraLinearOp(F),
91 Thyra::epetraLinearOp(Bt),
92 Thyra::epetraLinearOp(B),
93 Thyra::epetraLinearOp(C),
101 const Epetra_Map & rangeMap = epetra_A->OperatorRangeMap();
102 const Epetra_Map & domainMap = epetra_A->OperatorDomainMap();
109 TEST_EQUALITY(rangeMap.NumGlobalElements()-1, rangeMap.MaxAllGID());
110 TEST_EQUALITY(domainMap.NumGlobalElements()-1, domainMap.MaxAllGID());
115 Thyra::randomize(-100.0, 100.0, tv.
ptr());
117 Thyra::productVectorBase<double>(tv)->getVectorBlock(0);
119 Thyra::productVectorBase<double>(tv)->getVectorBlock(1);
120 const Thyra::ConstDetachedSpmdVectorView<double> vv_0(tv_0);
121 const Thyra::ConstDetachedSpmdVectorView<double> vv_1(tv_1);
123 int off_0 = vv_0.globalOffset();
124 int off_1 = vv_1.globalOffset();
127 Epetra_Vector ev(epetra_A->OperatorDomainMap());
128 epetra_A->copyThyraIntoEpetra(*tv, ev);
131 TEST_EQUALITY(tv->space()->dim(), as<Ordinal>(ev.GlobalLength()));
132 const int numMyElements = domainMap.NumMyElements();
134 for(
int i=0; i < numMyElements; i++) {
135 int gid = domainMap.GID(i);
137 tval = vv_0[gid-off_0];
139 tval = vv_1[gid-off_1-nx*ny];
147 Epetra_Vector ev(epetra_A->OperatorDomainMap());
153 Thyra::productVectorBase<double>(tv)->getVectorBlock(0);
155 Thyra::productVectorBase<double>(tv)->getVectorBlock(1);
156 const Thyra::ConstDetachedSpmdVectorView<double> vv_0(tv_0);
157 const Thyra::ConstDetachedSpmdVectorView<double> vv_1(tv_1);
159 int off_0 = rcp_dynamic_cast<
const Thyra::SpmdVectorSpaceBase<double> >(
160 tv_0->space())->localOffset();
161 int off_1 = rcp_dynamic_cast<
const Thyra::SpmdVectorSpaceBase<double> >(
162 tv_1->space())->localOffset();
164 epetra_A->copyEpetraIntoThyra(ev, tv.
ptr());
167 TEST_EQUALITY(tv->space()->dim(), as<Ordinal>(ev.GlobalLength()));
168 int numMyElements = domainMap.NumMyElements();
170 for(
int i=0;i<numMyElements;i++) {
171 int gid = domainMap.GID(i);
173 tval = vv_0[gid-off_0];
175 tval = vv_1[gid-off_1-nx*ny];
182 LinearOpTester<double> linearOpTester;
183 linearOpTester.show_all_tests(
true);
184 const bool checkResult = linearOpTester.check(*thyraEpetraOp, inOutArg(out));
TEUCHOS_UNIT_TEST(EpetraOperatorWrapper, basic)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Implements the Epetra_Operator interface with a Thyra LinearOperator.
TypeTo as(const TypeFrom &t)
#define TEST_EQUALITY(v1, v2)