Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
XpetraVectorInput.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Zoltan2: A package of combinatorial algorithms for scientific computing
4 //
5 // Copyright 2012 NTESS and the Zoltan2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 //
11 
17 #include <string>
18 
20 #include <Zoltan2_InputTraits.hpp>
21 #include <Zoltan2_TestHelpers.hpp>
22 
23 #include <Teuchos_DefaultComm.hpp>
24 #include <Teuchos_RCP.hpp>
25 #include <Teuchos_Comm.hpp>
26 #include <Teuchos_CommHelpers.hpp>
27 
28 using Teuchos::RCP;
29 using Teuchos::rcp;
30 using Teuchos::rcp_const_cast;
31 using Teuchos::Comm;
32 
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;
35 
36 void printVector(RCP<const Comm<int> > &comm, zlno_t vlen,
37  const zgno_t *vtxIds, const zscalar_t *vals)
38 {
39  int rank = comm->getRank();
40  int nprocs = comm->getSize();
41  comm->barrier();
42  for (int p=0; p < nprocs; p++){
43  if (p == rank){
44  std::cout << rank << ":" << std::endl;
45  for (zlno_t i=0; i < vlen; i++){
46  std::cout << " " << vtxIds[i] << ": " << vals[i] << std::endl;
47  }
48  std::cout.flush();
49  }
50  comm->barrier();
51  }
52  comm->barrier();
53 }
54 
55 template <typename User>
57  Zoltan2::XpetraMultiVectorAdapter<User> &ia, tvector_t &vector, int wdim,
58  zscalar_t **weights, int *strides)
59 {
60  RCP<const Comm<int> > comm = vector.getMap()->getComm();
61  int fail = 0, gfail=0;
62 
63  if (!fail && ia.getNumEntriesPerID() !=1)
64  fail = 42;
65 
66  if (!fail && ia.getNumWeightsPerID() !=wdim)
67  fail = 41;
68 
69  if (!fail && ia.getLocalNumIDs() != vector.getLocalLength())
70  fail = 4;
71 
72  gfail = globalFail(*comm, fail);
73 
74  if (!gfail){
75  const zgno_t *vtxIds=NULL;
76  const zscalar_t *vals=NULL;
77  int stride;
78 
79  size_t nvals = ia.getLocalNumIDs();
80  if (nvals != vector.getLocalLength())
81  fail = 8;
82 
83  ia.getIDsView(vtxIds);
84  ia.getEntriesView(vals, stride);
85 
86  if (!fail && stride != 1)
87  fail = 10;
88 
89  gfail = globalFail(*comm, fail);
90 
91  if (gfail == 0){
92  printVector(comm, nvals, vtxIds, vals);
93  }
94  }
95 
96  if (!gfail && wdim){
97  const zscalar_t *wgt =NULL;
98  int stride;
99 
100  for (int w=0; !fail && w < wdim; w++){
101  ia.getWeightsView(wgt, stride, w);
102 
103  if (!fail && stride != strides[w])
104  fail = 101;
105 
106  for (size_t v=0; !fail && v < vector.getLocalLength(); v++){
107  if (wgt[v*stride] != weights[w][v*stride])
108  fail=102;
109  }
110  }
111 
112  gfail = globalFail(*comm, fail);
113  }
114 
115  return gfail;
116 }
117 
118 
119 int main(int narg, char *arg[])
120 {
121  Tpetra::ScopeGuard tscope(&narg, &arg);
122  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
123 
124  int rank = comm->getRank();
125  int fail = 0, gfail=0;
126  bool aok = true;
127 
128  // Create object that can give us test Tpetra, Xpetra
129  // and Epetra vectors for testing.
130 
131  RCP<UserInputForTests> uinput;
132  Teuchos::ParameterList params;
133  params.set("input file", "simple");
134  params.set("file type", "Chaco");
135  try{
136  uinput = rcp(new UserInputForTests(params, comm));
137  }
138  catch(std::exception &e){
139  aok = false;
140  std::cout << e.what() << std::endl;
141  }
142  TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
143 
144  RCP<tvector_t> tV; // original vector (for checking)
145  RCP<tvector_t> newV; // migrated vector
146 
147  tV = rcp(new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap()));
148  tV->randomize();
149  size_t vlen = tV->getLocalLength();
150 
151  // To test migration in the input adapter we need a Solution object.
152 
153  RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(comm));
154 
155  int nWeights = 1;
156 
158  typedef adapter_t::part_t part_t;
159 
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);
163 
164  std::vector<const zscalar_t *> emptyWeights;
165  std::vector<int> emptyStrides;
166 
167  Zoltan2::PartitioningSolution<adapter_t> solution(env, comm, nWeights);
168  solution.setParts(solnParts);
169 
171  // User object is Tpetra::Vector, no weights
172  if (!gfail){
173  if (rank==0)
174  std::cout << "Constructed with Tpetra::Vector" << std::endl;
175 
176  RCP<const tvector_t> ctV = rcp_const_cast<const tvector_t>(tV);
177  RCP<adapter_t> tVInput;
178 
179  try {
180  tVInput =
182  emptyWeights, emptyStrides));
183  }
184  catch (std::exception &e){
185  aok = false;
186  std::cout << e.what() << std::endl;
187  }
188  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter ", 1);
189 
190  fail = verifyInputAdapter<tvector_t>(*tVInput, *tV, 0, NULL, NULL);
191 
192  gfail = globalFail(*comm, fail);
193 
194  if (!gfail){
195  tvector_t *vMigrate = NULL;
196  try{
197  tVInput->applyPartitioningSolution(*tV, vMigrate, solution);
198  newV = rcp(vMigrate);
199  }
200  catch (std::exception &e){
201  fail = 11;
202  }
203 
204  gfail = globalFail(*comm, fail);
205 
206  if (!gfail){
207  RCP<const tvector_t> cnewV = rcp_const_cast<const tvector_t>(newV);
208  RCP<Zoltan2::XpetraMultiVectorAdapter<tvector_t> > newInput;
209  try{
210  newInput = rcp(new Zoltan2::XpetraMultiVectorAdapter<tvector_t>(cnewV,
211  emptyWeights, emptyStrides));
212  }
213  catch (std::exception &e){
214  aok = false;
215  std::cout << e.what() << std::endl;
216  }
217  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 2 ", 1);
218 
219  if (rank==0){
220  std::cout << "Constructed with ";
221  std::cout << "Tpetra::Vector migrated to proc 0" << std::endl;
222  }
223  fail = verifyInputAdapter<tvector_t>(*newInput, *newV, 0, NULL, NULL);
224  if (fail) fail += 100;
225  gfail = globalFail(*comm, fail);
226  }
227  }
228  if (gfail){
229  printFailureCode(*comm, fail);
230  }
231  }
232 
234  // User object is Xpetra::Vector
235  if (!gfail){
236  if (rank==0)
237  std::cout << "Constructed with Xpetra::Vector" << std::endl;
238 
239  RCP<tvector_t> xtV =
240  rcp(new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap()));
241  xtV->randomize();
242  RCP<xvector_t> xV = Zoltan2::XpetraTraits<tvector_t>::convertToXpetra(xtV);
243  RCP<const xvector_t> cxV = rcp_const_cast<const xvector_t>(xV);
244  RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > xVInput;
245 
246  try {
247  xVInput =
249  emptyWeights, emptyStrides));
250  }
251  catch (std::exception &e){
252  aok = false;
253  std::cout << e.what() << std::endl;
254  }
255  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 3 ", 1);
256 
257  fail = verifyInputAdapter<xvector_t>(*xVInput, *tV, 0, NULL, NULL);
258 
259  gfail = globalFail(*comm, fail);
260 
261  if (!gfail){
262  xvector_t *vMigrate =NULL;
263  try{
264  xVInput->applyPartitioningSolution(*xV, vMigrate, solution);
265  }
266  catch (std::exception &e){
267  fail = 11;
268  }
269 
270  gfail = globalFail(*comm, fail);
271 
272  if (!gfail){
273  RCP<const xvector_t> cnewV(vMigrate);
274  RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > newInput;
275  try{
276  newInput =
278  emptyWeights, emptyStrides));
279  }
280  catch (std::exception &e){
281  aok = false;
282  std::cout << e.what() << std::endl;
283  }
284  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 4 ", 1);
285 
286  if (rank==0){
287  std::cout << "Constructed with ";
288  std::cout << "Xpetra::Vector migrated to proc 0" << std::endl;
289  }
290  fail = verifyInputAdapter<xvector_t>(*newInput, *newV, 0, NULL, NULL);
291  if (fail) fail += 100;
292  gfail = globalFail(*comm, fail);
293  }
294  }
295  if (gfail){
296  printFailureCode(*comm, fail);
297  }
298  }
299 
300 #ifdef HAVE_EPETRA_DATA_TYPES
301  // User object is Epetra_Vector
303  typedef Epetra_Vector evector_t;
304  if (!gfail){
305  if (rank==0)
306  std::cout << "Constructed with Epetra_Vector" << std::endl;
307 
308  RCP<evector_t> eV =
309  rcp(new Epetra_Vector(uinput->getUIEpetraCrsGraph()->RowMap()));
310  eV->Random();
311  RCP<const evector_t> ceV = rcp_const_cast<const evector_t>(eV);
312  RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > eVInput;
313 
314  bool goodAdapter = true;
315  try {
316  eVInput =
318  emptyWeights, emptyStrides));
319  }
320  catch (std::exception &e){
321  if (std::is_same<znode_t, Xpetra::EpetraNode>::value) {
322  aok = false;
323  goodAdapter = false;
324  std::cout << e.what() << std::endl;
325  }
326  else {
327  // We expect an error from Xpetra when znode_t != Xpetra::EpetraNode
328  // Ignore it, but skip tests using this matrix adapter.
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;
333  goodAdapter = false;
334  }
335  }
336  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 5 ", 1);
337 
338  if (goodAdapter) {
339  fail = verifyInputAdapter<evector_t>(*eVInput, *tV, 0, NULL, NULL);
340 
341  gfail = globalFail(*comm, fail);
342 
343  if (!gfail){
344  evector_t *vMigrate =NULL;
345  try{
346  eVInput->applyPartitioningSolution(*eV, vMigrate, solution);
347  }
348  catch (std::exception &e){
349  fail = 11;
350  }
351 
352  gfail = globalFail(*comm, fail);
353 
354  if (!gfail){
355  RCP<const evector_t> cnewV(vMigrate, true);
356  RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > newInput;
357  try{
358  newInput =
360  emptyWeights, emptyStrides));
361  }
362  catch (std::exception &e){
363  aok = false;
364  std::cout << e.what() << std::endl;
365  }
366  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 6 ", 1);
367 
368  if (rank==0){
369  std::cout << "Constructed with ";
370  std::cout << "Epetra_Vector migrated to proc 0" << std::endl;
371  }
372  fail = verifyInputAdapter<evector_t>(*newInput, *newV, 0, NULL, NULL);
373  if (fail) fail += 100;
374  gfail = globalFail(*comm, fail);
375  }
376  }
377  if (gfail){
378  printFailureCode(*comm, fail);
379  }
380  }
381  }
382 #endif
383 
385  // DONE
386 
387  if (rank==0)
388  std::cout << "PASS" << std::endl;
389 }
void printFailureCode(const Comm< int > &comm, int fail)
void verifyInputAdapter(adapter_t &ia, matrix_t &matrix)
static ArrayRCP< ArrayRCP< zscalar_t > > weights
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
Tpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > tvector_t
int main(int narg, char **arg)
Definition: coloring1.cpp:164
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.
Traits for application input objects.
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)
size_t getLocalNumIDs() const
Returns the number of objects on this process.
Xpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > xvector_t
float zscalar_t
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.