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_GlobalMPISession.hpp>
62 #include <Teuchos_DefaultComm.hpp>
63 #include <Teuchos_RCP.hpp>
64 #include <Teuchos_Comm.hpp>
65 #include <Teuchos_CommHelpers.hpp>
66 
67 using namespace std;
68 using Teuchos::RCP;
69 using Teuchos::rcp;
70 using Teuchos::rcp_const_cast;
71 using Teuchos::Comm;
72 using Teuchos::DefaultComm;
73 
74 typedef Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> tvector_t;
75 typedef Xpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> xvector_t;
76 typedef Epetra_MultiVector evector_t;
77 
78 template <typename User>
81  int wdim, zscalar_t **weights, int *strides)
82 {
83  RCP<const Comm<int> > comm = vector.getMap()->getComm();
84  int fail = 0, gfail=0;
85 
86  if (!fail && ia.getNumEntriesPerID() !=nvec)
87  fail = 42;
88 
89  if (!fail && ia.getNumWeightsPerID() !=wdim)
90  fail = 41;
91 
92  size_t length = vector.getLocalLength();
93 
94  if (!fail && ia.getLocalNumIDs() != length)
95  fail = 4;
96 
97  gfail = globalFail(comm, fail);
98 
99  if (!gfail){
100  const zgno_t *vtxIds=NULL;
101  const zscalar_t *vals=NULL;
102  int stride;
103 
104  size_t nvals = ia.getLocalNumIDs();
105  if (nvals != vector.getLocalLength())
106  fail = 8;
107 
108  ia.getIDsView(vtxIds);
109 
110  for (int v=0; v < nvec; v++){
111  ia.getEntriesView(vals, stride, v);
112 
113  if (!fail && stride != 1)
114  fail = 10;
115 
116  // TODO check the values returned
117  }
118 
119  gfail = globalFail(comm, fail);
120  }
121 
122  if (!gfail && wdim){
123  const zscalar_t *wgt =NULL;
124  int stride;
125 
126  for (int w=0; !fail && w < wdim; w++){
127  ia.getWeightsView(wgt, stride, w);
128 
129  if (!fail && stride != strides[w])
130  fail = 101;
131 
132  for (size_t v=0; !fail && v < vector.getLocalLength(); v++){
133  if (wgt[v*stride] != weights[w][v*stride])
134  fail=102;
135  }
136  }
137 
138  gfail = globalFail(comm, fail);
139  }
140 
141  return gfail;
142 }
143 
144 int main(int argc, char *argv[])
145 {
146  Teuchos::GlobalMPISession session(&argc, &argv);
147  RCP<const Comm<int> > comm = DefaultComm<int>::getComm();
148  int rank = comm->getRank();
149  int fail = 0, gfail=0;
150 
151  // Create object that can give us test Tpetra, Xpetra
152  // and Epetra vectors for testing.
153 
154  RCP<UserInputForTests> uinput;
155 
156  try{
157  uinput =
158  rcp(new UserInputForTests(testDataFilePath,std::string("simple"), comm, true));
159  }
160  catch(std::exception &e){
161  TEST_FAIL_AND_EXIT(*comm, 0, string("input ")+e.what(), 1);
162  }
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);
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  RCP<const tvector_t> ctV = rcp_const_cast<const tvector_t>(tV);
198  RCP<Zoltan2::XpetraMultiVectorAdapter<tvector_t> > tVInput;
199 
200  try {
201  tVInput =
203  emptyWeights, emptyStrides));
204  }
205  catch (std::exception &e){
206  TEST_FAIL_AND_EXIT(*comm, 0,
207  string("XpetraMultiVectorAdapter ")+e.what(), 1);
208  }
209 
210  if (rank==0){
211  std::cout << "Constructed with ";
212  std::cout << "Tpetra::MultiVector" << std::endl;
213  }
214 
215  fail = verifyInputAdapter<tvector_t>(*tVInput, *tV, nVec, 0, NULL, NULL);
216 
217  gfail = globalFail(comm, fail);
218 
219  if (!gfail){
220  tvector_t *vMigrate = NULL;
221  try{
222  tVInput->applyPartitioningSolution(*tV, vMigrate, solution);
223  newV = rcp(vMigrate);
224  }
225  catch (std::exception &e){
226  fail = 11;
227  }
228 
229  gfail = globalFail(comm, fail);
230 
231  if (!gfail){
232  RCP<const tvector_t> cnewV = rcp_const_cast<const tvector_t>(newV);
233  RCP<Zoltan2::XpetraMultiVectorAdapter<tvector_t> > newInput;
234  try{
236  cnewV, emptyWeights, emptyStrides));
237  }
238  catch (std::exception &e){
239  TEST_FAIL_AND_EXIT(*comm, 0,
240  string("XpetraMultiVectorAdapter 2 ")+e.what(), 1);
241  }
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  RCP<tvector_t> tMV =
261  rcp(new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap(), nVec));
262  tMV->randomize();
263  RCP<xvector_t> xV = Zoltan2::XpetraTraits<tvector_t>::convertToXpetra(tMV);
264  RCP<const xvector_t> cxV = rcp_const_cast<const xvector_t>(xV);
265  RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > xVInput;
266 
267  try {
268  xVInput =
270  emptyWeights, emptyStrides));
271  }
272  catch (std::exception &e){
273  TEST_FAIL_AND_EXIT(*comm, 0,
274  string("XpetraMultiVectorAdapter 3 ")+e.what(), 1);
275  }
276 
277  if (rank==0){
278  std::cout << "Constructed with ";
279  std::cout << "Xpetra::MultiVector" << std::endl;
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  TEST_FAIL_AND_EXIT(*comm, 0,
306  string("XpetraMultiVectorAdapter 4 ")+e.what(), 1);
307  }
308 
309  if (rank==0){
310  std::cout << "Constructed with ";
311  std::cout << "Xpetra::MultiVector migrated to proc 0" << std::endl;
312  }
313  fail = verifyInputAdapter<xvector_t>(*newInput, *newV, nVec, 0, NULL, NULL);
314  if (fail) fail += 100;
315  gfail = globalFail(comm, fail);
316  }
317  }
318  if (gfail){
319  printFailureCode(comm, fail);
320  }
321  }
322 
323 #ifdef HAVE_EPETRA_DATA_TYPES
324  // User object is Epetra_MultiVector
326  if (!gfail){
327  RCP<evector_t> eV =
328  rcp(new Epetra_MultiVector(uinput->getUIEpetraCrsGraph()->RowMap(),
329  nVec));
330  eV->Random();
331  RCP<const evector_t> ceV = rcp_const_cast<const evector_t>(eV);
332  RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > eVInput;
333 
334  try {
335  eVInput =
337  emptyWeights, emptyStrides));
338  }
339  catch (std::exception &e){
340  TEST_FAIL_AND_EXIT(*comm, 0,
341  string("XpetraMultiVectorAdapter 5 ")+e.what(), 1);
342  }
343 
344  if (rank==0){
345  std::cout << "Constructed with ";
346  std::cout << "Epetra_MultiVector" << std::endl;
347  }
348  fail = verifyInputAdapter<evector_t>(*eVInput, *tV, nVec, 0, NULL, NULL);
349 
350  gfail = globalFail(comm, fail);
351 
352  if (!gfail){
353  evector_t *vMigrate =NULL;
354  try{
355  eVInput->applyPartitioningSolution(*eV, vMigrate, solution);
356  }
357  catch (std::exception &e){
358  fail = 11;
359  }
360 
361  gfail = globalFail(comm, fail);
362 
363  if (!gfail){
364  RCP<const evector_t> cnewV(vMigrate, true);
365  RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > newInput;
366  try{
367  newInput =
369  emptyWeights, emptyStrides));
370  }
371  catch (std::exception &e){
372  TEST_FAIL_AND_EXIT(*comm, 0,
373  string("XpetraMultiVectorAdapter 6 ")+e.what(), 1);
374  }
375 
376  if (rank==0){
377  std::cout << "Constructed with ";
378  std::cout << "Epetra_MultiVector migrated to proc 0" << std::endl;
379  }
380  fail = verifyInputAdapter<evector_t>(*newInput, *newV, nVec, 0, NULL, NULL);
381  if (fail) fail += 100;
382  gfail = globalFail(comm, fail);
383  }
384  }
385  if (gfail){
386  printFailureCode(comm, fail);
387  }
388  }
389 #endif
390 
392  // DONE
393 
394  if (rank==0)
395  std::cout << "PASS" << std::endl;
396 }
int globalFail(const RCP< const Comm< int > > &comm, int fail)
double zscalar_t
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
Tpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > tvector_t
common code used by tests
static RCP< User > convertToXpetra(const RCP< User > &a)
Convert the object to its Xpetra wrapped version.
int verifyInputAdapter(Zoltan2::XpetraMultiVectorAdapter< User > &ia, tvector_t &vector, int nvec, int wdim, zscalar_t **weights, int *strides)
Defines the XpetraMultiVectorAdapter.
int getNumEntriesPerID() const
Return the number of vectors (typically one).
dictionary vals
Definition: xml2dox.py:186
A PartitioningSolution is a solution to a partitioning problem.
Epetra_MultiVector evector_t
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.
int zgno_t
void printFailureCode(const RCP< const Comm< int > > &comm, int fail)
int main(int argc, char *argv[])
size_t getLocalNumIDs() const
Returns the number of objects on this process.
Xpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > xvector_t
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
Provide pointer to a weight array with stride.
std::string testDataFilePath(".")