23 #include <Teuchos_DefaultComm.hpp>
24 #include <Teuchos_RCP.hpp>
25 #include <Teuchos_Comm.hpp>
26 #include <Teuchos_CommHelpers.hpp>
30 using Teuchos::rcp_const_cast;
33 typedef Tpetra::Vector<zscalar_t, zlno_t, zgno_t, znode_t>
tvector_t;
34 typedef Xpetra::Vector<zscalar_t, zlno_t, zgno_t, znode_t>
xvector_t;
39 int rank = comm->getRank();
40 int nprocs = comm->getSize();
42 for (
int p=0; p < nprocs; p++){
44 std::cout << rank <<
":" << std::endl;
45 for (
zlno_t i=0; i < vlen; i++){
46 std::cout <<
" " << vtxIds[i] <<
": " << vals[i] << std::endl;
55 template <
typename User>
60 RCP<const Comm<int> > comm = vector.getMap()->getComm();
61 int fail = 0, gfail=0;
80 if (nvals != vector.getLocalLength())
86 if (!fail && stride != 1)
100 for (
int w=0; !fail && w < wdim; w++){
103 if (!fail && stride != strides[w])
106 for (
size_t v=0; !fail && v < vector.getLocalLength(); v++){
107 if (wgt[v*stride] != weights[w][v*stride])
119 int main(
int narg,
char *arg[])
121 Tpetra::ScopeGuard tscope(&narg, &arg);
122 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
124 int rank = comm->getRank();
125 int fail = 0, gfail=0;
131 RCP<UserInputForTests> uinput;
132 Teuchos::ParameterList params;
133 params.set(
"input file",
"simple");
134 params.set(
"file type",
"Chaco");
138 catch(std::exception &e){
140 std::cout << e.what() << std::endl;
147 tV = rcp(
new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap()));
149 size_t vlen = tV->getLocalLength();
160 part_t *p =
new part_t [vlen];
161 memset(p, 0,
sizeof(part_t) * vlen);
162 ArrayRCP<part_t> solnParts(p, 0, vlen,
true);
164 std::vector<const zscalar_t *> emptyWeights;
165 std::vector<int> emptyStrides;
174 std::cout <<
"Constructed with Tpetra::Vector" << std::endl;
176 RCP<const tvector_t> ctV = rcp_const_cast<
const tvector_t>(tV);
177 RCP<adapter_t> tVInput;
182 emptyWeights, emptyStrides));
184 catch (std::exception &e){
186 std::cout << e.what() << std::endl;
190 fail = verifyInputAdapter<tvector_t>(*tVInput, *tV, 0, NULL, NULL);
197 tVInput->applyPartitioningSolution(*tV, vMigrate, solution);
198 newV = rcp(vMigrate);
200 catch (std::exception &e){
207 RCP<const tvector_t> cnewV = rcp_const_cast<
const tvector_t>(newV);
208 RCP<Zoltan2::XpetraMultiVectorAdapter<tvector_t> > newInput;
211 emptyWeights, emptyStrides));
213 catch (std::exception &e){
215 std::cout << e.what() << std::endl;
220 std::cout <<
"Constructed with ";
221 std::cout <<
"Tpetra::Vector migrated to proc 0" << std::endl;
223 fail = verifyInputAdapter<tvector_t>(*newInput, *newV, 0, NULL, NULL);
224 if (fail) fail += 100;
237 std::cout <<
"Constructed with Xpetra::Vector" << std::endl;
240 rcp(
new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap()));
243 RCP<const xvector_t> cxV = rcp_const_cast<
const xvector_t>(xV);
244 RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > xVInput;
249 emptyWeights, emptyStrides));
251 catch (std::exception &e){
253 std::cout << e.what() << std::endl;
257 fail = verifyInputAdapter<xvector_t>(*xVInput, *tV, 0, NULL, NULL);
264 xVInput->applyPartitioningSolution(*xV, vMigrate, solution);
266 catch (std::exception &e){
273 RCP<const xvector_t> cnewV(vMigrate);
274 RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > newInput;
278 emptyWeights, emptyStrides));
280 catch (std::exception &e){
282 std::cout << e.what() << std::endl;
287 std::cout <<
"Constructed with ";
288 std::cout <<
"Xpetra::Vector migrated to proc 0" << std::endl;
290 fail = verifyInputAdapter<xvector_t>(*newInput, *newV, 0, NULL, NULL);
291 if (fail) fail += 100;
300 #ifdef HAVE_EPETRA_DATA_TYPES
303 typedef Epetra_Vector evector_t;
306 std::cout <<
"Constructed with Epetra_Vector" << std::endl;
309 rcp(
new Epetra_Vector(uinput->getUIEpetraCrsGraph()->RowMap()));
311 RCP<const evector_t> ceV = rcp_const_cast<
const evector_t>(eV);
312 RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > eVInput;
314 bool goodAdapter =
true;
318 emptyWeights, emptyStrides));
320 catch (std::exception &e){
321 if (std::is_same<znode_t, Xpetra::EpetraNode>::value) {
324 std::cout << e.what() << std::endl;
329 std::cout <<
"Node type is not supported by Xpetra's Epetra interface;"
330 <<
" Skipping this test." << std::endl;
331 std::cout <<
"FYI: Here's the exception message: " << std::endl
332 << e.what() << std::endl;
339 fail = verifyInputAdapter<evector_t>(*eVInput, *tV, 0, NULL, NULL);
344 evector_t *vMigrate =NULL;
346 eVInput->applyPartitioningSolution(*eV, vMigrate, solution);
348 catch (std::exception &e){
355 RCP<const evector_t> cnewV(vMigrate,
true);
356 RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > newInput;
360 emptyWeights, emptyStrides));
362 catch (std::exception &e){
364 std::cout << e.what() << std::endl;
369 std::cout <<
"Constructed with ";
370 std::cout <<
"Epetra_Vector migrated to proc 0" << std::endl;
372 fail = verifyInputAdapter<evector_t>(*newInput, *newV, 0, NULL, NULL);
373 if (fail) fail += 100;
388 std::cout <<
"PASS" << std::endl;
void printFailureCode(const Comm< int > &comm, int fail)
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
int main(int narg, char **arg)
common code used by tests
void printVector(vector_t &vec, const std::string &msg, std::ostream &ostr=std::cout)
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.
Defines the XpetraMultiVectorAdapter.
int getNumEntriesPerID() const
Return the number of vectors.
SparseMatrixAdapter_t::part_t part_t
A PartitioningSolution is a solution to a partitioning problem.
void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const
Provide a pointer to the elements of the specified vector.
The user parameters, debug, timing and memory profiling output objects, and error checking methods...
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater...
An adapter for Xpetra::MultiVector.
Tpetra::Map::local_ordinal_type zlno_t
static const std::string fail
void setParts(ArrayRCP< part_t > &partList)
The algorithm uses setParts to set the solution.
int globalFail(const Comm< int > &comm, int fail)
void getIDsView(const gno_t *&ids) const
size_t getLocalNumIDs() const
Returns the number of objects on this process.
Tpetra::Map::global_ordinal_type zgno_t
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
Provide pointer to a weight array with stride.