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