13 #include "Epetra_MpiComm.h"
14 #include "Epetra_Map.h"
15 #include "Thyra_VectorSpaceFactoryBase.hpp"
17 #include "Thyra_DefaultSpmdVectorSpace_decl.hpp"
18 #include "Thyra_DefaultSpmdVector_decl.hpp"
19 #include "Thyra_MultiVectorBase_decl.hpp"
20 #include "Thyra_ScalarProdVectorSpaceBase_decl.hpp"
21 #include "Thyra_DefaultSpmdMultiVector_decl.hpp"
34 using Teuchos::outArg;
37 using Teuchos::rcp_dynamic_cast;
39 using Teuchos::reduceAll;
47 (void) MPI_Comm_rank (MPI_COMM_WORLD, &myRank);
48 (void) MPI_Comm_size (MPI_COMM_WORLD, &numProcs);
53 std::ostringstream os;
54 os <<
"(Process " << myRank <<
") ";
59 std::ostringstream os;
60 os << prefix <<
"Creating Epetra_Comm" << endl;
64 const Epetra_MpiComm epetra_comm (MPI_COMM_WORLD);
66 const Epetra_SerialComm epetra_comm ();
69 std::ostringstream os;
70 os << prefix <<
"Creating Teuchos::Comm" << endl;
73 RCP<const Teuchos::Comm<Teuchos_Ordinal> > comm =
77 lclSuccess = success ? 1 : 0;
78 reduceAll (*comm, REDUCE_MIN, lclSuccess, outArg (gblSuccess));
80 if (gblSuccess != 1) {
81 out <<
"FAILED; some process(es) have a null Teuchos::Comm" << endl;
86 const int localDim = myRank % 2;
87 const int globalDim = numProcs / 2;
88 RCP<const Epetra_Map> epetra_map;
90 std::ostringstream os;
91 os << prefix <<
"Creating Epetra_Map: localDim=" << localDim <<
", globalDim=" << globalDim << endl;
94 epetra_map =
rcp (
new Epetra_Map (globalDim, localDim, 0, epetra_comm));
97 lclSuccess = success ? 1 : 0;
98 reduceAll (*comm, REDUCE_MIN, lclSuccess, outArg (gblSuccess));
100 if (gblSuccess != 1) {
101 out <<
"FAILED; some process(es) have a null Epetra_Map" << endl;
106 std::ostringstream os;
107 os << prefix <<
"Creating Thyra::DefaultSpmdVectorSpace" << endl;
110 RCP<Thyra::DefaultSpmdVectorSpace<double> > SPMD =
111 Thyra::DefaultSpmdVectorSpace<double>::create();
115 lclSuccess = success ? 1 : 0;
116 reduceAll (*comm, REDUCE_MIN, lclSuccess, outArg (gblSuccess));
118 if (gblSuccess != 1) {
119 out <<
"FAILED; some process(es) have a null SPMD" << endl;
123 SPMD->initialize(comm, localDim, globalDim);
126 std::ostringstream os;
127 os << prefix <<
"Creating Thyra::MultiVectorBase" << endl;
130 RCP<const Thyra::MultiVectorBase<double> > spmd =
131 rcp (
new Thyra::DefaultSpmdMultiVector<double> (
133 rcp_dynamic_cast<
const Thyra::ScalarProdVectorSpaceBase<double> > (
134 SPMD->smallVecSpcFcty()->createVecSpc(1),
true)
139 lclSuccess = success ? 1 : 0;
140 reduceAll (*comm, REDUCE_MIN, lclSuccess, outArg (gblSuccess));
142 if (gblSuccess != 1) {
143 out <<
"FAILED; some process(es) have a null Thyra::MultiVectorBase"
149 std::ostringstream os;
150 os << prefix <<
"Calling Thyra::get_Epetra_MultiVector "
151 "(const overload; see #1941)" << endl;
155 RCP<const Epetra_MultiVector> mv_c =
157 const_cast<const Thyra::MultiVectorBase<double>&
> (*spmd));
160 lclSuccess = success ? 1 : 0;
161 reduceAll (*comm, REDUCE_MIN, lclSuccess, outArg (gblSuccess));
163 if (gblSuccess != 1) {
164 out <<
"FAILED; some process(es) have a null const Epetra_MultiVector"
170 std::ostringstream os;
171 os << prefix <<
"Calling Thyra::get_Epetra_MultiVector "
172 "(nonconst overload; see #2061)" << endl;
176 RCP<Epetra_MultiVector> mv_nc =
178 const_cast<Thyra::MultiVectorBase<double>&
> (*spmd));
181 lclSuccess = success ? 1 : 0;
182 reduceAll (*comm, REDUCE_MIN, lclSuccess, outArg (gblSuccess));
184 if (gblSuccess != 1) {
185 out <<
"FAILED; some process(es) have a null nonconst Epetra_MultiVector"
191 std::ostringstream os;
192 os << prefix <<
"Done with test on this process" << endl;
RCP< Epetra_MultiVector > get_Epetra_MultiVector(const Epetra_Map &map, const RCP< MultiVectorBase< double > > &mv)
Get a non-const Epetra_MultiVector view from a non-const MultiVectorBase object if possible...
static Teuchos::RCP< const Comm< OrdinalType > > getComm()
TEUCHOS_UNIT_TEST(EpetraOperatorWrapper, basic)
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
#define TEST_EQUALITY(v1, v2)