Zoltan2
XpetraVectorInput.cpp
Go to the documentation of this file.
1 // @HEADER
2 //
3 // ***********************************************************************
4 //
5 // Zoltan2: A package of combinatorial algorithms for scientific computing
6 // Copyright 2012 Sandia Corporation
7 //
8 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
9 // the U.S. Government retains certain rights in this software.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Karen Devine (kddevin@sandia.gov)
39 // Erik Boman (egboman@sandia.gov)
40 // Siva Rajamanickam (srajama@sandia.gov)
41 //
42 // ***********************************************************************
43 //
44 // @HEADER
45 //
46 
52 #include <string>
53 
55 #include <Zoltan2_InputTraits.hpp>
56 #include <Zoltan2_TestHelpers.hpp>
57 
58 #include <Teuchos_DefaultComm.hpp>
59 #include <Teuchos_RCP.hpp>
60 #include <Teuchos_Comm.hpp>
61 #include <Teuchos_CommHelpers.hpp>
62 
63 using Teuchos::RCP;
64 using Teuchos::rcp;
65 using Teuchos::rcp_const_cast;
66 using Teuchos::Comm;
67 
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;
70 
71 void printVector(RCP<const Comm<int> > &comm, zlno_t vlen,
72  const zgno_t *vtxIds, const zscalar_t *vals)
73 {
74  int rank = comm->getRank();
75  int nprocs = comm->getSize();
76  comm->barrier();
77  for (int p=0; p < nprocs; p++){
78  if (p == rank){
79  std::cout << rank << ":" << std::endl;
80  for (zlno_t i=0; i < vlen; i++){
81  std::cout << " " << vtxIds[i] << ": " << vals[i] << std::endl;
82  }
83  std::cout.flush();
84  }
85  comm->barrier();
86  }
87  comm->barrier();
88 }
89 
90 template <typename User>
92  Zoltan2::XpetraMultiVectorAdapter<User> &ia, tvector_t &vector, int wdim,
93  zscalar_t **weights, int *strides)
94 {
95  RCP<const Comm<int> > comm = vector.getMap()->getComm();
96  int fail = 0, gfail=0;
97 
98  if (!fail && ia.getNumEntriesPerID() !=1)
99  fail = 42;
100 
101  if (!fail && ia.getNumWeightsPerID() !=wdim)
102  fail = 41;
103 
104  if (!fail && ia.getLocalNumIDs() != vector.getLocalLength())
105  fail = 4;
106 
107  gfail = globalFail(comm, fail);
108 
109  if (!gfail){
110  const zgno_t *vtxIds=NULL;
111  const zscalar_t *vals=NULL;
112  int stride;
113 
114  size_t nvals = ia.getLocalNumIDs();
115  if (nvals != vector.getLocalLength())
116  fail = 8;
117 
118  ia.getIDsView(vtxIds);
119  ia.getEntriesView(vals, stride);
120 
121  if (!fail && stride != 1)
122  fail = 10;
123 
124  gfail = globalFail(comm, fail);
125 
126  if (gfail == 0){
127  printVector(comm, nvals, vtxIds, vals);
128  }
129  }
130 
131  if (!gfail && wdim){
132  const zscalar_t *wgt =NULL;
133  int stride;
134 
135  for (int w=0; !fail && w < wdim; w++){
136  ia.getWeightsView(wgt, stride, w);
137 
138  if (!fail && stride != strides[w])
139  fail = 101;
140 
141  for (size_t v=0; !fail && v < vector.getLocalLength(); v++){
142  if (wgt[v*stride] != weights[w][v*stride])
143  fail=102;
144  }
145  }
146 
147  gfail = globalFail(comm, fail);
148  }
149 
150  return gfail;
151 }
152 
153 int main(int narg, char *arg[])
154 {
155  Tpetra::ScopeGuard tscope(&narg, &arg);
156  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
157 
158  int rank = comm->getRank();
159  int fail = 0, gfail=0;
160  bool aok = true;
161 
162  // Create object that can give us test Tpetra, Xpetra
163  // and Epetra vectors for testing.
164 
165  RCP<UserInputForTests> uinput;
166 
167  try{
168  uinput =
169  rcp(new UserInputForTests(testDataFilePath,std::string("simple"), comm, true));
170  }
171  catch(std::exception &e){
172  aok = false;
173  std::cout << e.what() << std::endl;
174  }
175  TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
176 
177  RCP<tvector_t> tV; // original vector (for checking)
178  RCP<tvector_t> newV; // migrated vector
179 
180  tV = rcp(new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap()));
181  tV->randomize();
182  size_t vlen = tV->getLocalLength();
183 
184  // To test migration in the input adapter we need a Solution object.
185 
186  RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(comm));
187 
188  int nWeights = 1;
189 
191  typedef adapter_t::part_t part_t;
192 
193  part_t *p = new part_t [vlen];
194  memset(p, 0, sizeof(part_t) * vlen);
195  ArrayRCP<part_t> solnParts(p, 0, vlen, true);
196 
197  std::vector<const zscalar_t *> emptyWeights;
198  std::vector<int> emptyStrides;
199 
200  Zoltan2::PartitioningSolution<adapter_t> solution(env, comm, nWeights);
201  solution.setParts(solnParts);
202 
204  // User object is Tpetra::Vector, no weights
205  if (!gfail){
206  if (rank==0)
207  std::cout << "Constructed with Tpetra::Vector" << std::endl;
208 
209  RCP<const tvector_t> ctV = rcp_const_cast<const tvector_t>(tV);
210  RCP<adapter_t> tVInput;
211 
212  try {
213  tVInput =
215  emptyWeights, emptyStrides));
216  }
217  catch (std::exception &e){
218  aok = false;
219  std::cout << e.what() << std::endl;
220  }
221  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter ", 1);
222 
223  fail = verifyInputAdapter<tvector_t>(*tVInput, *tV, 0, NULL, NULL);
224 
225  gfail = globalFail(comm, fail);
226 
227  if (!gfail){
228  tvector_t *vMigrate = NULL;
229  try{
230  tVInput->applyPartitioningSolution(*tV, vMigrate, solution);
231  newV = rcp(vMigrate);
232  }
233  catch (std::exception &e){
234  fail = 11;
235  }
236 
237  gfail = globalFail(comm, fail);
238 
239  if (!gfail){
240  RCP<const tvector_t> cnewV = rcp_const_cast<const tvector_t>(newV);
241  RCP<Zoltan2::XpetraMultiVectorAdapter<tvector_t> > newInput;
242  try{
243  newInput = rcp(new Zoltan2::XpetraMultiVectorAdapter<tvector_t>(cnewV,
244  emptyWeights, emptyStrides));
245  }
246  catch (std::exception &e){
247  aok = false;
248  std::cout << e.what() << std::endl;
249  }
250  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 2 ", 1);
251 
252  if (rank==0){
253  std::cout << "Constructed with ";
254  std::cout << "Tpetra::Vector migrated to proc 0" << std::endl;
255  }
256  fail = verifyInputAdapter<tvector_t>(*newInput, *newV, 0, NULL, NULL);
257  if (fail) fail += 100;
258  gfail = globalFail(comm, fail);
259  }
260  }
261  if (gfail){
262  printFailureCode(comm, fail);
263  }
264  }
265 
267  // User object is Xpetra::Vector
268  if (!gfail){
269  if (rank==0)
270  std::cout << "Constructed with Xpetra::Vector" << std::endl;
271 
272  RCP<tvector_t> xtV =
273  rcp(new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap()));
274  xtV->randomize();
275  RCP<xvector_t> xV = Zoltan2::XpetraTraits<tvector_t>::convertToXpetra(xtV);
276  RCP<const xvector_t> cxV = rcp_const_cast<const xvector_t>(xV);
277  RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > xVInput;
278 
279  try {
280  xVInput =
282  emptyWeights, emptyStrides));
283  }
284  catch (std::exception &e){
285  aok = false;
286  std::cout << e.what() << std::endl;
287  }
288  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 3 ", 1);
289 
290  fail = verifyInputAdapter<xvector_t>(*xVInput, *tV, 0, NULL, NULL);
291 
292  gfail = globalFail(comm, fail);
293 
294  if (!gfail){
295  xvector_t *vMigrate =NULL;
296  try{
297  xVInput->applyPartitioningSolution(*xV, vMigrate, solution);
298  }
299  catch (std::exception &e){
300  fail = 11;
301  }
302 
303  gfail = globalFail(comm, fail);
304 
305  if (!gfail){
306  RCP<const xvector_t> cnewV(vMigrate);
307  RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > newInput;
308  try{
309  newInput =
311  emptyWeights, emptyStrides));
312  }
313  catch (std::exception &e){
314  aok = false;
315  std::cout << e.what() << std::endl;
316  }
317  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 4 ", 1);
318 
319  if (rank==0){
320  std::cout << "Constructed with ";
321  std::cout << "Xpetra::Vector migrated to proc 0" << std::endl;
322  }
323  fail = verifyInputAdapter<xvector_t>(*newInput, *newV, 0, NULL, NULL);
324  if (fail) fail += 100;
325  gfail = globalFail(comm, fail);
326  }
327  }
328  if (gfail){
329  printFailureCode(comm, fail);
330  }
331  }
332 
333 #ifdef HAVE_EPETRA_DATA_TYPES
334  // User object is Epetra_Vector
336  typedef Epetra_Vector evector_t;
337  if (!gfail){
338  if (rank==0)
339  std::cout << "Constructed with Epetra_Vector" << std::endl;
340 
341  RCP<evector_t> eV =
342  rcp(new Epetra_Vector(uinput->getUIEpetraCrsGraph()->RowMap()));
343  eV->Random();
344  RCP<const evector_t> ceV = rcp_const_cast<const evector_t>(eV);
345  RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > eVInput;
346 
347  bool goodAdapter = true;
348  try {
349  eVInput =
351  emptyWeights, emptyStrides));
352  }
353  catch (std::exception &e){
354  if (std::is_same<znode_t, Xpetra::EpetraNode>::value) {
355  aok = false;
356  goodAdapter = false;
357  std::cout << e.what() << std::endl;
358  }
359  else {
360  // We expect an error from Xpetra when znode_t != Xpetra::EpetraNode
361  // Ignore it, but skip tests using this matrix adapter.
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;
366  goodAdapter = false;
367  }
368  }
369  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 5 ", 1);
370 
371  if (goodAdapter) {
372  fail = verifyInputAdapter<evector_t>(*eVInput, *tV, 0, NULL, NULL);
373 
374  gfail = globalFail(comm, fail);
375 
376  if (!gfail){
377  evector_t *vMigrate =NULL;
378  try{
379  eVInput->applyPartitioningSolution(*eV, vMigrate, solution);
380  }
381  catch (std::exception &e){
382  fail = 11;
383  }
384 
385  gfail = globalFail(comm, fail);
386 
387  if (!gfail){
388  RCP<const evector_t> cnewV(vMigrate, true);
389  RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > newInput;
390  try{
391  newInput =
393  emptyWeights, emptyStrides));
394  }
395  catch (std::exception &e){
396  aok = false;
397  std::cout << e.what() << std::endl;
398  }
399  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 6 ", 1);
400 
401  if (rank==0){
402  std::cout << "Constructed with ";
403  std::cout << "Epetra_Vector migrated to proc 0" << std::endl;
404  }
405  fail = verifyInputAdapter<evector_t>(*newInput, *newV, 0, NULL, NULL);
406  if (fail) fail += 100;
407  gfail = globalFail(comm, fail);
408  }
409  }
410  if (gfail){
411  printFailureCode(comm, fail);
412  }
413  }
414  }
415 #endif
416 
418  // DONE
419 
420  if (rank==0)
421  std::cout << "PASS" << std::endl;
422 }
Zoltan2_InputTraits.hpp
Traits for application input objects.
globalFail
int globalFail(const RCP< const Comm< int > > &comm, int fail)
Definition: ErrorHandlingForTests.hpp:109
xvector_t
Xpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > xvector_t
Definition: XpetraMultiVectorInput.cpp:72
printFailureCode
void printFailureCode(const RCP< const Comm< int > > &comm, int fail)
Definition: ErrorHandlingForTests.hpp:116
Zoltan2::XpetraTraits::convertToXpetra
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.
Definition: Zoltan2_XpetraTraits.hpp:98
Zoltan2::PartitioningSolution
A PartitioningSolution is a solution to a partitioning problem.
Definition: Zoltan2_PartitioningSolution.hpp:55
zscalar_t
double zscalar_t
Definition: Zoltan2_TestHelpers.hpp:141
Zoltan2::XpetraMultiVectorAdapter::getWeightsView
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
Definition: Zoltan2_XpetraMultiVectorAdapter.hpp:143
testDataFilePath
std::string testDataFilePath(".")
printVector
void printVector(RCP< const Comm< int > > &comm, zlno_t vlen, const zgno_t *vtxIds, const zscalar_t *vals)
Definition: XpetraVectorInput.cpp:71
main
int main(int narg, char *arg[])
Definition: XpetraVectorInput.cpp:153
xml2dox.vals
dictionary vals
Definition: xml2dox.py:186
xvector_t
Xpetra::Vector< zscalar_t, zlno_t, zgno_t, znode_t > xvector_t
Definition: XpetraVectorInput.cpp:69
UserInputForTests
Definition: UserInputForTests.hpp:126
Zoltan2::XpetraMultiVectorAdapter
An adapter for Xpetra::MultiVector.
Definition: Zoltan2_XpetraMultiVectorAdapter.hpp:83
part_t
SparseMatrixAdapter_t::part_t part_t
Definition: partitioningTree.cpp:74
Zoltan2::XpetraMultiVectorAdapter::getLocalNumIDs
size_t getLocalNumIDs() const
Returns the number of objects on this process.
Definition: Zoltan2_XpetraMultiVectorAdapter.hpp:134
Zoltan2::Environment
The user parameters, debug, timing and memory profiling output objects, and error checking methods.
Definition: Zoltan2_Environment.hpp:83
TEST_FAIL_AND_EXIT
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
Definition: ErrorHandlingForTests.hpp:70
zgno_t
int zgno_t
Definition: Zoltan2_TestHelpers.hpp:143
Zoltan2::XpetraMultiVectorAdapter::getEntriesView
void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const
Definition: Zoltan2_XpetraMultiVectorAdapter.hpp:238
zlno_t
int zlno_t
Definition: Zoltan2_TestHelpers.hpp:142
tvector_t
Tpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > tvector_t
Definition: XpetraMultiVectorInput.cpp:71
Zoltan2::XpetraMultiVectorAdapter::getNumEntriesPerID
int getNumEntriesPerID() const
Return the number of vectors (typically one).
Definition: Zoltan2_XpetraMultiVectorAdapter.hpp:161
verifyInputAdapter
int verifyInputAdapter(Zoltan2::XpetraMultiVectorAdapter< User > &ia, tvector_t &vector, int wdim, zscalar_t **weights, int *strides)
Definition: XpetraVectorInput.cpp:91
Zoltan2::PartitioningSolution::setParts
void setParts(ArrayRCP< part_t > &partList)
The algorithm uses setParts to set the solution.
Definition: Zoltan2_PartitioningSolution.hpp:1255
Zoltan2::XpetraMultiVectorAdapter::getIDsView
void getIDsView(const gno_t *&ids) const
Provide a pointer to this process' identifiers.
Definition: Zoltan2_XpetraMultiVectorAdapter.hpp:136
Zoltan2::XpetraMultiVectorAdapter::getNumWeightsPerID
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater....
Definition: Zoltan2_XpetraMultiVectorAdapter.hpp:141
fail
static const std::string fail
Definition: findUniqueGids.cpp:80
weights
static ArrayRCP< ArrayRCP< zscalar_t > > weights
Definition: rcbPerformanceZ1.cpp:82
Zoltan2_TestHelpers.hpp
common code used by tests
Zoltan2_XpetraMultiVectorAdapter.hpp
Defines the XpetraMultiVectorAdapter.
tvector_t
Tpetra::Vector< zscalar_t, zlno_t, zgno_t, znode_t > tvector_t
Definition: XpetraVectorInput.cpp:68