Zoltan2
XpetraMultiVectorInput.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 // @HEADER
46 // ***********************************************************************
47 // Zoltan2: Sandia Partitioning Ordering & Coloring Library
48 // ***********************************************************************
49 
55 #include <string>
56 
58 #include <Zoltan2_InputTraits.hpp>
59 #include <Zoltan2_TestHelpers.hpp>
60 
61 #include <Teuchos_DefaultComm.hpp>
62 #include <Teuchos_RCP.hpp>
63 #include <Teuchos_Comm.hpp>
64 #include <Teuchos_CommHelpers.hpp>
65 
66 using Teuchos::RCP;
67 using Teuchos::rcp;
68 using Teuchos::rcp_const_cast;
69 using Teuchos::Comm;
70 
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;
73 
74 template <typename User>
77  int wdim, zscalar_t **weights, int *strides)
78 {
79  RCP<const Comm<int> > comm = vector.getMap()->getComm();
80  int fail = 0, gfail=0;
81 
82  if (!fail && ia.getNumEntriesPerID() !=nvec)
83  fail = 42;
84 
85  if (!fail && ia.getNumWeightsPerID() !=wdim)
86  fail = 41;
87 
88  size_t length = vector.getLocalLength();
89 
90  if (!fail && ia.getLocalNumIDs() != length)
91  fail = 4;
92 
93  gfail = globalFail(comm, fail);
94 
95  if (!gfail){
96  const zgno_t *vtxIds=NULL;
97  const zscalar_t *vals=NULL;
98  int stride;
99 
100  size_t nvals = ia.getLocalNumIDs();
101  if (nvals != vector.getLocalLength())
102  fail = 8;
103 
104  ia.getIDsView(vtxIds);
105 
106  for (int v=0; v < nvec; v++){
107  ia.getEntriesView(vals, stride, v);
108 
109  if (!fail && stride != 1)
110  fail = 10;
111 
112  // TODO check the values returned
113  }
114 
115  gfail = globalFail(comm, fail);
116  }
117 
118  if (!gfail && wdim){
119  const zscalar_t *wgt =NULL;
120  int stride;
121 
122  for (int w=0; !fail && w < wdim; w++){
123  ia.getWeightsView(wgt, stride, w);
124 
125  if (!fail && stride != strides[w])
126  fail = 101;
127 
128  for (size_t v=0; !fail && v < vector.getLocalLength(); v++){
129  if (wgt[v*stride] != weights[w][v*stride])
130  fail=102;
131  }
132  }
133 
134  gfail = globalFail(comm, fail);
135  }
136 
137  return gfail;
138 }
139 
140 int main(int narg, char *arg[])
141 {
142  Tpetra::ScopeGuard tscope(&narg, &arg);
143  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
144 
145  int rank = comm->getRank();
146  int fail = 0, gfail=0;
147  bool aok = true;
148 
149  // Create object that can give us test Tpetra, Xpetra
150  // and Epetra vectors for testing.
151 
152  RCP<UserInputForTests> uinput;
153 
154  try{
155  uinput =
156  rcp(new UserInputForTests(testDataFilePath,std::string("simple"), comm, true));
157  }
158  catch(std::exception &e){
159  aok = false;
160  std::cout << e.what() << std::endl;
161  }
162  TEST_FAIL_AND_EXIT(*comm, aok, "input ", 1);
163 
164  RCP<tvector_t> tV; // original vector (for checking)
165  RCP<tvector_t> newV; // migrated vector
166 
167  int nVec = 2;
168 
169  tV = rcp(new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), nVec));
170  tV->randomize();
171 
172  size_t vlen = tV->getLocalLength();
173 
174  // To test migration in the input adapter we need a Solution object.
175 
176  RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(comm));
177 
178  int nWeights = 1;
179 
182  typedef ia_t::part_t part_t;
183 
184  part_t *p = new part_t [vlen];
185  memset(p, 0, sizeof(part_t) * vlen);
186  ArrayRCP<part_t> solnParts(p, 0, vlen, true);
187 
188  soln_t solution(env, comm, nWeights);
189  solution.setParts(solnParts);
190 
191  std::vector<const zscalar_t *> emptyWeights;
192  std::vector<int> emptyStrides;
193 
195  // User object is Tpetra::MultiVector, no weights
196  if (!gfail){
197  if (rank==0)
198  std::cout << "Constructed with Tpetra::MultiVector" << std::endl;
199 
200  RCP<const tvector_t> ctV = rcp_const_cast<const tvector_t>(tV);
201  RCP<Zoltan2::XpetraMultiVectorAdapter<tvector_t> > tVInput;
202 
203  try {
204  tVInput =
206  emptyWeights, emptyStrides));
207  }
208  catch (std::exception &e){
209  aok = false;
210  std::cout << e.what() << std::endl;
211  }
212  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter ", 1);
213 
214  fail = verifyInputAdapter<tvector_t>(*tVInput, *tV, nVec, 0, NULL, NULL);
215 
216  gfail = globalFail(comm, fail);
217 
218  if (!gfail){
219  tvector_t *vMigrate = NULL;
220  try{
221  tVInput->applyPartitioningSolution(*tV, vMigrate, solution);
222  newV = rcp(vMigrate);
223  }
224  catch (std::exception &e){
225  fail = 11;
226  }
227 
228  gfail = globalFail(comm, fail);
229 
230  if (!gfail){
231  RCP<const tvector_t> cnewV = rcp_const_cast<const tvector_t>(newV);
232  RCP<Zoltan2::XpetraMultiVectorAdapter<tvector_t> > newInput;
233  try{
235  cnewV, emptyWeights, emptyStrides));
236  }
237  catch (std::exception &e){
238  aok = false;
239  std::cout << e.what() << std::endl;
240  }
241  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 2 ", 1);
242 
243  if (rank==0){
244  std::cout << "Constructed with ";
245  std::cout << "Tpetra::MultiVector migrated to proc 0" << std::endl;
246  }
247  fail = verifyInputAdapter<tvector_t>(*newInput, *newV, nVec, 0, NULL, NULL);
248  if (fail) fail += 100;
249  gfail = globalFail(comm, fail);
250  }
251  }
252  if (gfail){
253  printFailureCode(comm, fail);
254  }
255  }
256 
258  // User object is Xpetra::MultiVector
259  if (!gfail){
260  if (rank==0)
261  std::cout << "Constructed with Xpetra::MultiVector" << std::endl;
262 
263  RCP<tvector_t> tMV =
264  rcp(new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), nVec));
265  tMV->randomize();
266  RCP<xvector_t> xV = Zoltan2::XpetraTraits<tvector_t>::convertToXpetra(tMV);
267  RCP<const xvector_t> cxV = rcp_const_cast<const xvector_t>(xV);
268  RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > xVInput;
269 
270  try {
271  xVInput =
273  emptyWeights, emptyStrides));
274  }
275  catch (std::exception &e){
276  aok = false;
277  std::cout << e.what() << std::endl;
278  }
279  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 3 ", 1);
280 
281  fail = verifyInputAdapter<xvector_t>(*xVInput, *tV, nVec, 0, NULL, NULL);
282 
283  gfail = globalFail(comm, fail);
284 
285  if (!gfail){
286  xvector_t *vMigrate =NULL;
287  try{
288  xVInput->applyPartitioningSolution(*xV, vMigrate, solution);
289  }
290  catch (std::exception &e){
291  fail = 11;
292  }
293 
294  gfail = globalFail(comm, fail);
295 
296  if (!gfail){
297  RCP<const xvector_t> cnewV(vMigrate);
298  RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > newInput;
299  try{
300  newInput =
302  cnewV, emptyWeights, emptyStrides));
303  }
304  catch (std::exception &e){
305  aok = false;
306  std::cout << e.what() << std::endl;
307  }
308  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 4 ", 1);
309 
310  if (rank==0){
311  std::cout << "Constructed with ";
312  std::cout << "Xpetra::MultiVector migrated to proc 0" << std::endl;
313  }
314  fail = verifyInputAdapter<xvector_t>(*newInput, *newV, nVec, 0, NULL, NULL);
315  if (fail) fail += 100;
316  gfail = globalFail(comm, fail);
317  }
318  }
319  if (gfail){
320  printFailureCode(comm, fail);
321  }
322  }
323 
324 #ifdef HAVE_EPETRA_DATA_TYPES
325  // User object is Epetra_MultiVector
327  typedef Epetra_MultiVector evector_t;
328  if (!gfail){
329  if (rank==0)
330  std::cout << "Constructed with Epetra_MultiVector" << std::endl;
331 
332  RCP<evector_t> eV =
333  rcp(new Epetra_MultiVector(uinput->getUIEpetraCrsGraph()->RowMap(),
334  nVec));
335  eV->Random();
336  RCP<const evector_t> ceV = rcp_const_cast<const evector_t>(eV);
337  RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > eVInput;
338 
339  bool goodAdapter = true;
340  try {
341  eVInput =
343  emptyWeights, emptyStrides));
344  }
345  catch (std::exception &e){
346  if (std::is_same<znode_t, Xpetra::EpetraNode>::value) {
347  aok = false;
348  goodAdapter = false;
349  std::cout << e.what() << std::endl;
350  }
351  else {
352  // We expect an error from Xpetra when znode_t != Xpetra::EpetraNode
353  // Ignore it, but skip tests using this matrix adapter.
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;
358  goodAdapter = false;
359  }
360  }
361  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 5 ", 1);
362 
363  if (goodAdapter) {
364  fail = verifyInputAdapter<evector_t>(*eVInput, *tV, nVec, 0, NULL, NULL);
365 
366  gfail = globalFail(comm, fail);
367 
368  if (!gfail){
369  evector_t *vMigrate =NULL;
370  try{
371  eVInput->applyPartitioningSolution(*eV, vMigrate, solution);
372  }
373  catch (std::exception &e){
374  fail = 11;
375  }
376 
377  gfail = globalFail(comm, fail);
378 
379  if (!gfail){
380  RCP<const evector_t> cnewV(vMigrate, true);
381  RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > newInput;
382  try{
383  newInput =
385  emptyWeights, emptyStrides));
386  }
387  catch (std::exception &e){
388  aok = false;
389  std::cout << e.what() << std::endl;
390  }
391  TEST_FAIL_AND_EXIT(*comm, aok, "XpetraMultiVectorAdapter 6 ", 1);
392 
393  if (rank==0){
394  std::cout << "Constructed with ";
395  std::cout << "Epetra_MultiVector migrated to proc 0" << std::endl;
396  }
397  fail = verifyInputAdapter<evector_t>(*newInput, *newV, nVec, 0, NULL, NULL);
398  if (fail) fail += 100;
399  gfail = globalFail(comm, fail);
400  }
401  }
402  if (gfail){
403  printFailureCode(comm, fail);
404  }
405  }
406  }
407 #endif
408 
410  // DONE
411 
412  if (rank==0)
413  std::cout << "PASS" << std::endl;
414 }
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(".")
xml2dox.vals
dictionary vals
Definition: xml2dox.py:186
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
main
int main(int narg, char *arg[])
Definition: XpetraMultiVectorInput.cpp:140
Zoltan2::XpetraMultiVectorAdapter::getEntriesView
void getEntriesView(const scalar_t *&elements, int &stride, int idx=0) const
Definition: Zoltan2_XpetraMultiVectorAdapter.hpp:238
verifyInputAdapter
int verifyInputAdapter(Zoltan2::XpetraMultiVectorAdapter< User > &ia, tvector_t &vector, int nvec, int wdim, zscalar_t **weights, int *strides)
Definition: XpetraMultiVectorInput.cpp:75
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
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.