61 #include <Teuchos_DefaultComm.hpp>
62 #include <Teuchos_RCP.hpp>
63 #include <Teuchos_Comm.hpp>
64 #include <Teuchos_CommHelpers.hpp>
68 using Teuchos::rcp_const_cast;
71 typedef Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t>
tvector_t;
72 typedef Xpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t>
xvector_t;
74 template <
typename User>
79 RCP<const Comm<int> > comm = vector.getMap()->getComm();
80 int fail = 0, gfail=0;
88 size_t length = vector.getLocalLength();
101 if (nvals != vector.getLocalLength())
106 for (
int v=0; v < nvec; v++){
109 if (!
fail && stride != 1)
122 for (
int w=0; !
fail && w < wdim; w++){
125 if (!
fail && stride != strides[w])
128 for (
size_t v=0; !
fail && v < vector.getLocalLength(); v++){
129 if (wgt[v*stride] !=
weights[w][v*stride])
140 int main(
int narg,
char *arg[])
142 Tpetra::ScopeGuard tscope(&narg, &arg);
143 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
145 int rank = comm->getRank();
146 int fail = 0, gfail=0;
152 RCP<UserInputForTests> uinput;
158 catch(std::exception &e){
160 std::cout << e.what() << std::endl;
169 tV = rcp(
new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), nVec));
172 size_t vlen = tV->getLocalLength();
185 memset(p, 0,
sizeof(
part_t) * vlen);
186 ArrayRCP<part_t> solnParts(p, 0, vlen,
true);
188 soln_t solution(env, comm, nWeights);
189 solution.setParts(solnParts);
191 std::vector<const zscalar_t *> emptyWeights;
192 std::vector<int> emptyStrides;
198 std::cout <<
"Constructed with Tpetra::MultiVector" << std::endl;
200 RCP<const tvector_t> ctV = rcp_const_cast<const tvector_t>(tV);
201 RCP<Zoltan2::XpetraMultiVectorAdapter<tvector_t> > tVInput;
206 emptyWeights, emptyStrides));
208 catch (std::exception &e){
210 std::cout << e.what() << std::endl;
214 fail = verifyInputAdapter<tvector_t>(*tVInput, *tV, nVec, 0, NULL, NULL);
221 tVInput->applyPartitioningSolution(*tV, vMigrate, solution);
222 newV = rcp(vMigrate);
224 catch (std::exception &e){
231 RCP<const tvector_t> cnewV = rcp_const_cast<const tvector_t>(newV);
232 RCP<Zoltan2::XpetraMultiVectorAdapter<tvector_t> > newInput;
235 cnewV, emptyWeights, emptyStrides));
237 catch (std::exception &e){
239 std::cout << e.what() << std::endl;
244 std::cout <<
"Constructed with ";
245 std::cout <<
"Tpetra::MultiVector migrated to proc 0" << std::endl;
247 fail = verifyInputAdapter<tvector_t>(*newInput, *newV, nVec, 0, NULL, NULL);
261 std::cout <<
"Constructed with Xpetra::MultiVector" << std::endl;
264 rcp(
new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), nVec));
267 RCP<const xvector_t> cxV = rcp_const_cast<const xvector_t>(xV);
268 RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > xVInput;
273 emptyWeights, emptyStrides));
275 catch (std::exception &e){
277 std::cout << e.what() << std::endl;
281 fail = verifyInputAdapter<xvector_t>(*xVInput, *tV, nVec, 0, NULL, NULL);
288 xVInput->applyPartitioningSolution(*xV, vMigrate, solution);
290 catch (std::exception &e){
297 RCP<const xvector_t> cnewV(vMigrate);
298 RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > newInput;
302 cnewV, emptyWeights, emptyStrides));
304 catch (std::exception &e){
306 std::cout << e.what() << std::endl;
311 std::cout <<
"Constructed with ";
312 std::cout <<
"Xpetra::MultiVector migrated to proc 0" << std::endl;
314 fail = verifyInputAdapter<xvector_t>(*newInput, *newV, nVec, 0, NULL, NULL);
324 #ifdef HAVE_EPETRA_DATA_TYPES
327 typedef Epetra_MultiVector evector_t;
330 std::cout <<
"Constructed with Epetra_MultiVector" << std::endl;
333 rcp(
new Epetra_MultiVector(uinput->getUIEpetraCrsGraph()->RowMap(),
336 RCP<const evector_t> ceV = rcp_const_cast<const evector_t>(eV);
337 RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > eVInput;
339 bool goodAdapter =
true;
343 emptyWeights, emptyStrides));
345 catch (std::exception &e){
346 if (std::is_same<znode_t, Xpetra::EpetraNode>::value) {
349 std::cout << e.what() << std::endl;
354 std::cout <<
"Node type is not supported by Xpetra's Epetra interface;"
355 <<
" Skipping this test." << std::endl;
356 std::cout <<
"FYI: Here's the exception message: " << std::endl
357 << e.what() << std::endl;
364 fail = verifyInputAdapter<evector_t>(*eVInput, *tV, nVec, 0, NULL, NULL);
369 evector_t *vMigrate =NULL;
371 eVInput->applyPartitioningSolution(*eV, vMigrate, solution);
373 catch (std::exception &e){
380 RCP<const evector_t> cnewV(vMigrate,
true);
381 RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > newInput;
385 emptyWeights, emptyStrides));
387 catch (std::exception &e){
389 std::cout << e.what() << std::endl;
394 std::cout <<
"Constructed with ";
395 std::cout <<
"Epetra_MultiVector migrated to proc 0" << std::endl;
397 fail = verifyInputAdapter<evector_t>(*newInput, *newV, nVec, 0, NULL, NULL);
413 std::cout <<
"PASS" << std::endl;