Zoltan2
BasicCoordinateInput.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 // Test for Zoltan2::BasicVectorAdapter for coordinate-based problems
47 
49 #include <Zoltan2_TestHelpers.hpp>
50 
51 #include <Teuchos_DefaultComm.hpp>
52 #include <Teuchos_RCP.hpp>
53 #include <Teuchos_CommHelpers.hpp>
54 
55 using Teuchos::RCP;
56 using Teuchos::Comm;
57 using Teuchos::Array;
58 
60 
63  int len, int glen, zgno_t *ids,
64  zscalar_t *xyz,
66  int nCoords, int nWeights)
67 {
68  int fail = 0;
69 
70  if (ia->getNumEntriesPerID() != nCoords)
71  fail = 100;
72 
73  if (!fail && ia->getNumWeightsPerID() != nWeights)
74  fail = 101;
75 
76  if (!fail && ia->getLocalNumIDs() != size_t(len))
77  fail = 102;
78 
79  for (int x=0; !fail && x < nCoords; x++){
80  const zgno_t *idList;
81  const zscalar_t *vals;
82  int stride;
83 
84  ia->getIDsView(idList);
85  ia->getEntriesView(vals, stride, x);
86 
87  zscalar_t *coordVal = xyz + x;
88  for (int i=0; !fail && i < len; i++, coordVal += 3){
89 
90  if (idList[i] != ids[i])
91  fail = 105;
92 
93  if (!fail && vals[stride*i] != *coordVal)
94  fail = 106;
95  }
96  }
97 
98  for (int w=0; !fail && w < nWeights; w++){
99  const zscalar_t *wgts;
100  int stride;
101 
102  ia->getWeightsView(wgts, stride, w);
103 
104  zscalar_t *weightVal = weights + len*w;
105  for (int i=0; !fail && i < len; i++, weightVal++){
106  if (wgts[stride*i] != *weightVal)
107  fail = 110;
108  }
109  }
110 
111  return fail;
112 }
113 
114 
115 int main(int narg, char *arg[])
116 {
117  Tpetra::ScopeGuard tscope(&narg, &arg);
118  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
119 
120  int rank = comm->getRank();
121  int fail = 0;
122 
123  // Get some coordinates
124 
125  typedef Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> mv_t;
126  RCP<UserInputForTests> uinput;
127  std::string fname("simple");
128 
129  try{
130  uinput = rcp(new UserInputForTests(testDataFilePath, fname, comm, true));
131  }
132  catch(std::exception &e){
133  fail=1;
134  }
135 
136  TEST_FAIL_AND_EXIT(*comm, !fail, "input constructor", 1);
137 
138  RCP<mv_t> coords;
139 
140  try{
141  coords = uinput->getUICoordinates();
142  }
143  catch(std::exception &e){
144  fail=1;
145  }
146 
147  TEST_FAIL_AND_EXIT(*comm, !fail, "getting coordinates", 1);
148 
149  int numLocalIds = coords->getLocalLength();
150  int numGlobalIds = coords->getGlobalLength();
151  int coordDim = coords->getNumVectors();
152  ArrayView<const zgno_t> idList = coords->getMap()->getNodeElementList();
153 
154  // Create global Ids, x-, y- and z-coordinates, and also arrays of weights.
155 
156  Array<zgno_t> myIds(numLocalIds);
157  zgno_t myFirstId = rank * numLocalIds;
158 
159  int wdim = 2;
160  Array<zscalar_t> weights(numLocalIds*wdim);
161  for (int i = 0; i < numLocalIds*wdim; i++) weights[i] = zscalar_t(i);
162 
163  zscalar_t *x_values= coords->getDataNonConst(0).getRawPtr();
164  zscalar_t *y_values= x_values; // fake 3 dimensions if needed
165  zscalar_t *z_values= x_values;
166 
167  if (coordDim > 1){
168  y_values= coords->getDataNonConst(1).getRawPtr();
169  if (coordDim > 2)
170  z_values= coords->getDataNonConst(2).getRawPtr();
171  }
172 
173  Array<zscalar_t> xyz_values(3*numLocalIds);
174 
175  for (zlno_t i=0; i < numLocalIds; i++) // global Ids
176  myIds[i] = myFirstId+i;
177 
178  zscalar_t *x = xyz_values.getRawPtr(); // a stride-3 coordinate array
179  zscalar_t *y = x+1;
180  zscalar_t *z = y+1;
181 
182  for (int i=0, ii=0; i < numLocalIds; i++, ii += 3){
183  x[ii] = x_values[i];
184  y[ii] = y_values[i];
185  z[ii] = z_values[i];
186  }
187 
188  RCP<Zoltan2::BasicVectorAdapter<userTypes_t> > ia;
189 
190  {
192  // 3-dimensional coordinates with stride one and no weights,
193  // using simpler constructor
194 
195  int ncoords = 3;
196  int nweights = 0;
197 
198  try{
200  numLocalIds, myIds.getRawPtr(), x_values, y_values, z_values));
201  }
202  catch (std::exception &e){
203  fail = 1;
204  }
205 
206  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 0", fail);
207 
208  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
209  myIds.getRawPtr(), xyz_values.getRawPtr(),
210  weights.getRawPtr(), ncoords, nweights);
211 
212  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 0", fail);
213  }
214 
215  {
217  // 3-dimensional coordinates with stride one and one weight
218  // using simpler constructor
219 
220  int ncoords = 3;
221  int nweights = 1;
222 
223  try{
225  numLocalIds, myIds.getRawPtr(),
226  x_values, y_values, z_values, 1, 1, 1,
227  true, weights.getRawPtr(), 1));
228  }
229  catch (std::exception &e){
230  fail = 1;
231  }
232 
233  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 0a", fail);
234 
235  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
236  myIds.getRawPtr(), xyz_values.getRawPtr(),
237  weights.getRawPtr(), ncoords, nweights);
238 
239  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 0a", fail);
240  }
241 
242  {
244  // 3-dimensional coordinates with stride one and no weights
245 
246  int ncoords = 3;
247  int nweights = 0;
248 
249  std::vector<const zscalar_t *> values, weightValues;
250  std::vector<int> valueStrides, weightStrides;
251 
252  values.push_back(x_values);
253  values.push_back(y_values);
254  values.push_back(z_values);
255  valueStrides.push_back(1);
256  valueStrides.push_back(1);
257  valueStrides.push_back(1);
258 
259  try{
261  numLocalIds, myIds.getRawPtr(), values, valueStrides,
262  weightValues, weightStrides));
263  }
264  catch (std::exception &e){
265  fail = 1;
266  }
267 
268  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 1", fail);
269 
270  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
271  myIds.getRawPtr(), xyz_values.getRawPtr(),
272  weights.getRawPtr(), ncoords, nweights);
273 
274  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 1", fail);
275 
276  // Try using the default: no strides supplied means strides are one.
277 
278  std::vector<int> emptyStrides;
279 
280  try{
282  numLocalIds, myIds.getRawPtr(), values, emptyStrides,
283  weightValues, emptyStrides));
284  }
285  catch (std::exception &e){
286  fail = 1;
287  }
288 
289  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 2", fail);
290 
291  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
292  myIds.getRawPtr(), xyz_values.getRawPtr(),
293  weights.getRawPtr(), ncoords, nweights);
294 
295  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 2", fail);
296  }
297 
298  {
300  // 2-dimensional coordinates with stride three and two weights
301 
302  int ncoords = 2;
303  int nweights = 2;
304 
305  std::vector<const zscalar_t *> values, weightValues;
306  std::vector<int> valueStrides, weightStrides;
307 
308  values.push_back(xyz_values.getRawPtr());
309  values.push_back(xyz_values.getRawPtr() + 1);
310  valueStrides.push_back(3);
311  valueStrides.push_back(3);
312 
313  weightValues.push_back(weights.getRawPtr());
314  weightValues.push_back(weights.getRawPtr() + numLocalIds);
315  weightStrides.push_back(1);
316  weightStrides.push_back(1);
317 
318  try{
320  numLocalIds, myIds.getRawPtr(), values, valueStrides,
321  weightValues, weightStrides));
322  }
323  catch (std::exception &e){
324  fail = 1;
325  }
326 
327  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 3", fail);
328 
329  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
330  myIds.getRawPtr(), xyz_values.getRawPtr(),
331  weights.getRawPtr(), ncoords, nweights);
332 
333  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 3", fail);
334 
335  // Try using default weight strides
336 
337  std::vector<int> emptyStrides;
338 
339  try{
341  numLocalIds, myIds.getRawPtr(), values, valueStrides,
342  weightValues, emptyStrides));
343  }
344  catch (std::exception &e){
345  fail = 1;
346  }
347 
348  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 4", fail);
349 
350  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
351  myIds.getRawPtr(), xyz_values.getRawPtr(),
352  weights.getRawPtr(), ncoords, nweights);
353 
354  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 4", fail);
355  }
356 
357  {
359  // 1-dimensional coordinates with stride one and two weights
360 
361  int ncoords = 1;
362  int nweights = 2;
363 
364  std::vector<const zscalar_t *> values, weightValues;
365  std::vector<int> valueStrides, weightStrides;
366 
367  values.push_back(x_values);
368  valueStrides.push_back(1);
369 
370  weightValues.push_back(weights.getRawPtr());
371  weightValues.push_back(weights.getRawPtr() + numLocalIds);
372  weightStrides.push_back(1);
373  weightStrides.push_back(1);
374 
375  try{
377  numLocalIds, myIds.getRawPtr(), values, valueStrides,
378  weightValues, weightStrides));
379  }
380  catch (std::exception &e){
381  fail = 1;
382  }
383 
384  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 4", fail);
385 
386  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
387  myIds.getRawPtr(), xyz_values.getRawPtr(),
388  weights.getRawPtr(), ncoords, nweights);
389 
390  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 4", fail);
391  }
392 
393  if (rank == 0)
394  std::cout << "PASS" << std::endl;
395 
396  return fail;
397 }
398 
userTypes_t
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > userTypes_t
Definition: BasicCoordinateInput.cpp:59
Zoltan2::BasicVectorAdapter::getLocalNumIDs
size_t getLocalNumIDs() const
Returns the number of objects on this process.
Definition: Zoltan2_BasicVectorAdapter.hpp:237
Zoltan2_BasicVectorAdapter.hpp
Defines the BasicVectorAdapter class.
zscalar_t
double zscalar_t
Definition: Zoltan2_TestHelpers.hpp:141
testDataFilePath
std::string testDataFilePath(".")
checkXMLParameters.idList
list idList
Match up parameters to validators.
Definition: checkXMLParameters.py:50
xml2dox.vals
dictionary vals
Definition: xml2dox.py:186
UserInputForTests
Definition: UserInputForTests.hpp:126
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::BasicVectorAdapter::getNumWeightsPerID
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater....
Definition: Zoltan2_BasicVectorAdapter.hpp:241
Zoltan2::BasicVectorAdapter::getEntriesView
void getEntriesView(const scalar_t *&entries, int &stride, int idx=0) const
Definition: Zoltan2_BasicVectorAdapter.hpp:261
Zoltan2::BasicVectorAdapter::getWeightsView
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
Definition: Zoltan2_BasicVectorAdapter.hpp:243
Zoltan2::BasicUserTypes
A simple class that can be the User template argument for an InputAdapter.
Definition: Zoltan2_InputTraits.hpp:139
main
int main(int narg, char *arg[])
Definition: BasicCoordinateInput.cpp:115
zlno_t
int zlno_t
Definition: Zoltan2_TestHelpers.hpp:142
Zoltan2::BasicVectorAdapter::getNumEntriesPerID
int getNumEntriesPerID() const
Return the number of vectors (typically one).
Definition: Zoltan2_BasicVectorAdapter.hpp:259
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
TEST_FAIL_AND_RETURN_VALUE
#define TEST_FAIL_AND_RETURN_VALUE(comm, ok, s, rc)
Definition: ErrorHandlingForTests.hpp:96
checkBasicCoordinate
int checkBasicCoordinate(Zoltan2::BasicVectorAdapter< userTypes_t > *ia, int len, int glen, zgno_t *ids, zscalar_t *xyz, zscalar_t *weights, int nCoords, int nWeights)
Definition: BasicCoordinateInput.cpp:61
Zoltan2::BasicVectorAdapter::getIDsView
void getIDsView(const gno_t *&ids) const
Provide a pointer to this process' identifiers.
Definition: Zoltan2_BasicVectorAdapter.hpp:239
Zoltan2::BasicVectorAdapter
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
Definition: Zoltan2_BasicVectorAdapter.hpp:81