58 #include <Teuchos_DefaultComm.hpp>
59 #include <Teuchos_RCP.hpp>
60 #include <Teuchos_Comm.hpp>
61 #include <Teuchos_CommHelpers.hpp>
65 using Teuchos::rcp_const_cast;
68 typedef Tpetra::Vector<zscalar_t, zlno_t, zgno_t, znode_t>
tvector_t;
69 typedef Xpetra::Vector<zscalar_t, zlno_t, zgno_t, znode_t>
xvector_t;
74 int rank = comm->getRank();
75 int nprocs = comm->getSize();
77 for (
int p=0; p < nprocs; p++){
79 std::cout << rank <<
":" << std::endl;
80 for (
zlno_t i=0; i < vlen; i++){
81 std::cout <<
" " << vtxIds[i] <<
": " <<
vals[i] << std::endl;
90 template <
typename User>
95 RCP<const Comm<int> > comm = vector.getMap()->getComm();
96 int fail = 0, gfail=0;
110 const zgno_t *vtxIds=NULL;
115 if (nvals != vector.getLocalLength())
121 if (!
fail && stride != 1)
135 for (
int w=0; !
fail && w < wdim; w++){
138 if (!
fail && stride != strides[w])
141 for (
size_t v=0; !
fail && v < vector.getLocalLength(); v++){
142 if (wgt[v*stride] !=
weights[w][v*stride])
153 int main(
int narg,
char *arg[])
155 Tpetra::ScopeGuard tscope(&narg, &arg);
156 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
158 int rank = comm->getRank();
159 int fail = 0, gfail=0;
165 RCP<UserInputForTests> uinput;
171 catch(std::exception &e){
173 std::cout << e.what() << std::endl;
180 tV = rcp(
new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap()));
182 size_t vlen = tV->getLocalLength();
194 memset(p, 0,
sizeof(
part_t) * vlen);
195 ArrayRCP<part_t> solnParts(p, 0, vlen,
true);
197 std::vector<const zscalar_t *> emptyWeights;
198 std::vector<int> emptyStrides;
207 std::cout <<
"Constructed with Tpetra::Vector" << std::endl;
209 RCP<const tvector_t> ctV = rcp_const_cast<const tvector_t>(tV);
210 RCP<adapter_t> tVInput;
215 emptyWeights, emptyStrides));
217 catch (std::exception &e){
219 std::cout << e.what() << std::endl;
223 fail = verifyInputAdapter<tvector_t>(*tVInput, *tV, 0, NULL, NULL);
230 tVInput->applyPartitioningSolution(*tV, vMigrate, solution);
231 newV = rcp(vMigrate);
233 catch (std::exception &e){
240 RCP<const tvector_t> cnewV = rcp_const_cast<const tvector_t>(newV);
241 RCP<Zoltan2::XpetraMultiVectorAdapter<tvector_t> > newInput;
244 emptyWeights, emptyStrides));
246 catch (std::exception &e){
248 std::cout << e.what() << std::endl;
253 std::cout <<
"Constructed with ";
254 std::cout <<
"Tpetra::Vector migrated to proc 0" << std::endl;
256 fail = verifyInputAdapter<tvector_t>(*newInput, *newV, 0, NULL, NULL);
270 std::cout <<
"Constructed with Xpetra::Vector" << std::endl;
273 rcp(
new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap()));
276 RCP<const xvector_t> cxV = rcp_const_cast<const xvector_t>(xV);
277 RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > xVInput;
282 emptyWeights, emptyStrides));
284 catch (std::exception &e){
286 std::cout << e.what() << std::endl;
290 fail = verifyInputAdapter<xvector_t>(*xVInput, *tV, 0, NULL, NULL);
297 xVInput->applyPartitioningSolution(*xV, vMigrate, solution);
299 catch (std::exception &e){
306 RCP<const xvector_t> cnewV(vMigrate);
307 RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > newInput;
311 emptyWeights, emptyStrides));
313 catch (std::exception &e){
315 std::cout << e.what() << std::endl;
320 std::cout <<
"Constructed with ";
321 std::cout <<
"Xpetra::Vector migrated to proc 0" << std::endl;
323 fail = verifyInputAdapter<xvector_t>(*newInput, *newV, 0, NULL, NULL);
333 #ifdef HAVE_EPETRA_DATA_TYPES
336 typedef Epetra_Vector evector_t;
339 std::cout <<
"Constructed with Epetra_Vector" << std::endl;
342 rcp(
new Epetra_Vector(uinput->getUIEpetraCrsGraph()->RowMap()));
344 RCP<const evector_t> ceV = rcp_const_cast<const evector_t>(eV);
345 RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > eVInput;
347 bool goodAdapter =
true;
351 emptyWeights, emptyStrides));
353 catch (std::exception &e){
354 if (std::is_same<znode_t, Xpetra::EpetraNode>::value) {
357 std::cout << e.what() << std::endl;
362 std::cout <<
"Node type is not supported by Xpetra's Epetra interface;"
363 <<
" Skipping this test." << std::endl;
364 std::cout <<
"FYI: Here's the exception message: " << std::endl
365 << e.what() << std::endl;
372 fail = verifyInputAdapter<evector_t>(*eVInput, *tV, 0, NULL, NULL);
377 evector_t *vMigrate =NULL;
379 eVInput->applyPartitioningSolution(*eV, vMigrate, solution);
381 catch (std::exception &e){
388 RCP<const evector_t> cnewV(vMigrate,
true);
389 RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > newInput;
393 emptyWeights, emptyStrides));
395 catch (std::exception &e){
397 std::cout << e.what() << std::endl;
402 std::cout <<
"Constructed with ";
403 std::cout <<
"Epetra_Vector migrated to proc 0" << std::endl;
405 fail = verifyInputAdapter<evector_t>(*newInput, *newV, 0, NULL, NULL);
421 std::cout <<
"PASS" << std::endl;