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_GlobalMPISession.hpp>
59 #include <Teuchos_DefaultComm.hpp>
60 #include <Teuchos_RCP.hpp>
61 #include <Teuchos_Comm.hpp>
62 #include <Teuchos_CommHelpers.hpp>
63 
64 using namespace std;
65 using Teuchos::RCP;
66 using Teuchos::rcp;
67 using Teuchos::rcp_const_cast;
68 using Teuchos::Comm;
69 using Teuchos::DefaultComm;
70 
71 typedef Tpetra::Vector<zscalar_t, zlno_t, zgno_t, znode_t> tvector_t;
72 typedef Xpetra::Vector<zscalar_t, zlno_t, zgno_t, znode_t> xvector_t;
73 typedef Epetra_Vector evector_t;
74 
75 void printVector(RCP<const Comm<int> > &comm, zlno_t vlen,
76  const zgno_t *vtxIds, const zscalar_t *vals)
77 {
78  int rank = comm->getRank();
79  int nprocs = comm->getSize();
80  comm->barrier();
81  for (int p=0; p < nprocs; p++){
82  if (p == rank){
83  std::cout << rank << ":" << std::endl;
84  for (zlno_t i=0; i < vlen; i++){
85  std::cout << " " << vtxIds[i] << ": " << vals[i] << std::endl;
86  }
87  std::cout.flush();
88  }
89  comm->barrier();
90  }
91  comm->barrier();
92 }
93 
94 template <typename User>
96  Zoltan2::XpetraMultiVectorAdapter<User> &ia, tvector_t &vector, int wdim,
97  zscalar_t **weights, int *strides)
98 {
99  RCP<const Comm<int> > comm = vector.getMap()->getComm();
100  int fail = 0, gfail=0;
101 
102  if (!fail && ia.getNumEntriesPerID() !=1)
103  fail = 42;
104 
105  if (!fail && ia.getNumWeightsPerID() !=wdim)
106  fail = 41;
107 
108  if (!fail && ia.getLocalNumIDs() != vector.getLocalLength())
109  fail = 4;
110 
111  gfail = globalFail(comm, fail);
112 
113  if (!gfail){
114  const zgno_t *vtxIds=NULL;
115  const zscalar_t *vals=NULL;
116  int stride;
117 
118  size_t nvals = ia.getLocalNumIDs();
119  if (nvals != vector.getLocalLength())
120  fail = 8;
121 
122  ia.getIDsView(vtxIds);
123  ia.getEntriesView(vals, stride);
124 
125  if (!fail && stride != 1)
126  fail = 10;
127 
128  gfail = globalFail(comm, fail);
129 
130  if (gfail == 0){
131  printVector(comm, nvals, vtxIds, vals);
132  }
133  }
134 
135  if (!gfail && wdim){
136  const zscalar_t *wgt =NULL;
137  int stride;
138 
139  for (int w=0; !fail && w < wdim; w++){
140  ia.getWeightsView(wgt, stride, w);
141 
142  if (!fail && stride != strides[w])
143  fail = 101;
144 
145  for (size_t v=0; !fail && v < vector.getLocalLength(); v++){
146  if (wgt[v*stride] != weights[w][v*stride])
147  fail=102;
148  }
149  }
150 
151  gfail = globalFail(comm, fail);
152  }
153 
154  return gfail;
155 }
156 
157 int main(int argc, char *argv[])
158 {
159  Teuchos::GlobalMPISession session(&argc, &argv);
160  RCP<const Comm<int> > comm = DefaultComm<int>::getComm();
161  int rank = comm->getRank();
162  int fail = 0, gfail=0;
163 
164  // Create object that can give us test Tpetra, Xpetra
165  // and Epetra vectors for testing.
166 
167  RCP<UserInputForTests> uinput;
168 
169  try{
170  uinput =
171  rcp(new UserInputForTests(testDataFilePath,std::string("simple"), comm, true));
172  }
173  catch(std::exception &e){
174  TEST_FAIL_AND_EXIT(*comm, 0, string("input ")+e.what(), 1);
175  }
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);
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  RCP<const tvector_t> ctV = rcp_const_cast<const tvector_t>(tV);
207  RCP<adapter_t> tVInput;
208 
209  try {
210  tVInput =
212  emptyWeights, emptyStrides));
213  }
214  catch (std::exception &e){
215  TEST_FAIL_AND_EXIT(*comm, 0,
216  string("XpetraMultiVectorAdapter ")+e.what(), 1);
217  }
218 
219  if (rank==0){
220  std::cout << "Constructed with ";
221  std::cout << "Tpetra::Vector" << std::endl;
222  }
223 
224  fail = verifyInputAdapter<tvector_t>(*tVInput, *tV, 0, NULL, NULL);
225 
226  gfail = globalFail(comm, fail);
227 
228  if (!gfail){
229  tvector_t *vMigrate = NULL;
230  try{
231  tVInput->applyPartitioningSolution(*tV, vMigrate, solution);
232  newV = rcp(vMigrate);
233  }
234  catch (std::exception &e){
235  fail = 11;
236  }
237 
238  gfail = globalFail(comm, fail);
239 
240  if (!gfail){
241  RCP<const tvector_t> cnewV = rcp_const_cast<const tvector_t>(newV);
242  RCP<Zoltan2::XpetraMultiVectorAdapter<tvector_t> > newInput;
243  try{
244  newInput = rcp(new Zoltan2::XpetraMultiVectorAdapter<tvector_t>(cnewV,
245  emptyWeights, emptyStrides));
246  }
247  catch (std::exception &e){
248  TEST_FAIL_AND_EXIT(*comm, 0,
249  string("XpetraMultiVectorAdapter 2 ")+e.what(), 1);
250  }
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  RCP<tvector_t> xtV =
270  rcp(new tvector_t(uinput->getUITpetraCrsGraph()->getRowMap()));
271  xtV->randomize();
272  RCP<xvector_t> xV = Zoltan2::XpetraTraits<tvector_t>::convertToXpetra(xtV);
273  RCP<const xvector_t> cxV = rcp_const_cast<const xvector_t>(xV);
274  RCP<Zoltan2::XpetraMultiVectorAdapter<xvector_t> > xVInput;
275 
276  try {
277  xVInput =
279  emptyWeights, emptyStrides));
280  }
281  catch (std::exception &e){
282  TEST_FAIL_AND_EXIT(*comm, 0,
283  string("XpetraMultiVectorAdapter 3 ")+e.what(), 1);
284  }
285 
286  if (rank==0){
287  std::cout << "Constructed with ";
288  std::cout << "Xpetra::Vector" << std::endl;
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  TEST_FAIL_AND_EXIT(*comm, 0,
315  string("XpetraMultiVectorAdapter 4 ")+e.what(), 1);
316  }
317 
318  if (rank==0){
319  std::cout << "Constructed with ";
320  std::cout << "Xpetra::Vector migrated to proc 0" << std::endl;
321  }
322  fail = verifyInputAdapter<xvector_t>(*newInput, *newV, 0, NULL, NULL);
323  if (fail) fail += 100;
324  gfail = globalFail(comm, fail);
325  }
326  }
327  if (gfail){
328  printFailureCode(comm, fail);
329  }
330  }
331 
332 #ifdef HAVE_EPETRA_DATA_TYPES
333  // User object is Epetra_Vector
335  if (!gfail){
336  RCP<evector_t> eV =
337  rcp(new Epetra_Vector(uinput->getUIEpetraCrsGraph()->RowMap()));
338  eV->Random();
339  RCP<const evector_t> ceV = rcp_const_cast<const evector_t>(eV);
340  RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > eVInput;
341 
342  try {
343  eVInput =
345  emptyWeights, emptyStrides));
346  }
347  catch (std::exception &e){
348  TEST_FAIL_AND_EXIT(*comm, 0,
349  string("XpetraMultiVectorAdapter 5 ")+e.what(), 1);
350  }
351 
352  if (rank==0){
353  std::cout << "Constructed with ";
354  std::cout << "Epetra_Vector" << std::endl;
355  }
356  fail = verifyInputAdapter<evector_t>(*eVInput, *tV, 0, NULL, NULL);
357 
358  gfail = globalFail(comm, fail);
359 
360  if (!gfail){
361  evector_t *vMigrate =NULL;
362  try{
363  eVInput->applyPartitioningSolution(*eV, vMigrate, solution);
364  }
365  catch (std::exception &e){
366  fail = 11;
367  }
368 
369  gfail = globalFail(comm, fail);
370 
371  if (!gfail){
372  RCP<const evector_t> cnewV(vMigrate, true);
373  RCP<Zoltan2::XpetraMultiVectorAdapter<evector_t> > newInput;
374  try{
375  newInput =
377  emptyWeights, emptyStrides));
378  }
379  catch (std::exception &e){
380  TEST_FAIL_AND_EXIT(*comm, 0,
381  string("XpetraMultiVectorAdapter 6 ")+e.what(), 1);
382  }
383 
384  if (rank==0){
385  std::cout << "Constructed with ";
386  std::cout << "Epetra_Vector migrated to proc 0" << std::endl;
387  }
388  fail = verifyInputAdapter<evector_t>(*newInput, *newV, 0, NULL, NULL);
389  if (fail) fail += 100;
390  gfail = globalFail(comm, fail);
391  }
392  }
393  if (gfail){
394  printFailureCode(comm, fail);
395  }
396  }
397 #endif
398 
400  // DONE
401 
402  if (rank==0)
403  std::cout << "PASS" << std::endl;
404 }
int globalFail(const RCP< const Comm< int > > &comm, int fail)
double zscalar_t
int main(int argc, char *argv[])
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
Tpetra::Vector< zscalar_t, zlno_t, zgno_t, znode_t > tvector_t
Epetra_Vector evector_t
Tpetra::MultiVector< zscalar_t, zlno_t, zgno_t, znode_t > tvector_t
int zlno_t
common code used by tests
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 (typically one).
int verifyInputAdapter(Zoltan2::XpetraMultiVectorAdapter< User > &ia, tvector_t &vector, int wdim, zscalar_t **weights, int *strides)
dictionary vals
Definition: xml2dox.py:186
Xpetra::Vector< zscalar_t, zlno_t, zgno_t, znode_t > xvector_t
A PartitioningSolution is a solution to a partitioning problem.
Epetra_MultiVector evector_t
void printVector(RCP< const Comm< int > > &comm, zlno_t vlen, const zgno_t *vtxIds, const zscalar_t *vals)
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.
void setParts(ArrayRCP< part_t > &partList)
The algorithm uses setParts to set the solution.
int zgno_t
void printFailureCode(const RCP< 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
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
Provide pointer to a weight array with stride.
std::string testDataFilePath(".")