Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
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  Teuchos::ParameterList params;
128  params.set("input file", "simple");
129  params.set("file type", "Chaco");
130 
131  try{
132  uinput = rcp(new UserInputForTests(params, comm));
133  }
134  catch(std::exception &e){
135  fail=1;
136  }
137 
138  TEST_FAIL_AND_EXIT(*comm, !fail, "input constructor", 1);
139 
140  RCP<mv_t> coords;
141 
142  try{
143  coords = uinput->getUICoordinates();
144  }
145  catch(std::exception &e){
146  fail=1;
147  }
148 
149  TEST_FAIL_AND_EXIT(*comm, !fail, "getting coordinates", 1);
150 
151  int numLocalIds = coords->getLocalLength();
152  int numGlobalIds = coords->getGlobalLength();
153  int coordDim = coords->getNumVectors();
154  ArrayView<const zgno_t> idList = coords->getMap()->getNodeElementList();
155 
156  // Create global Ids, x-, y- and z-coordinates, and also arrays of weights.
157 
158  Array<zgno_t> myIds(numLocalIds);
159  zgno_t myFirstId = rank * numLocalIds;
160 
161  int wdim = 2;
162  Array<zscalar_t> weights(numLocalIds*wdim);
163  for (int i = 0; i < numLocalIds*wdim; i++) weights[i] = zscalar_t(i);
164 
165  zscalar_t *x_values= coords->getDataNonConst(0).getRawPtr();
166  zscalar_t *y_values= x_values; // fake 3 dimensions if needed
167  zscalar_t *z_values= x_values;
168 
169  if (coordDim > 1){
170  y_values= coords->getDataNonConst(1).getRawPtr();
171  if (coordDim > 2)
172  z_values= coords->getDataNonConst(2).getRawPtr();
173  }
174 
175  Array<zscalar_t> xyz_values(3*numLocalIds);
176 
177  for (zlno_t i=0; i < numLocalIds; i++) // global Ids
178  myIds[i] = myFirstId+i;
179 
180  zscalar_t *x = xyz_values.getRawPtr(); // a stride-3 coordinate array
181  zscalar_t *y = x+1;
182  zscalar_t *z = y+1;
183 
184  for (int i=0, ii=0; i < numLocalIds; i++, ii += 3){
185  x[ii] = x_values[i];
186  y[ii] = y_values[i];
187  z[ii] = z_values[i];
188  }
189 
190  RCP<Zoltan2::BasicVectorAdapter<userTypes_t> > ia;
191 
192  {
194  // 3-dimensional coordinates with stride one and no weights,
195  // using simpler constructor
196 
197  int ncoords = 3;
198  int nweights = 0;
199 
200  try{
202  numLocalIds, myIds.getRawPtr(), x_values, y_values, z_values));
203  }
204  catch (std::exception &e){
205  fail = 1;
206  }
207 
208  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 0", fail);
209 
210  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
211  myIds.getRawPtr(), xyz_values.getRawPtr(),
212  weights.getRawPtr(), ncoords, nweights);
213 
214  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 0", fail);
215  }
216 
217  {
219  // 3-dimensional coordinates with stride one and one weight
220  // using simpler constructor
221 
222  int ncoords = 3;
223  int nweights = 1;
224 
225  try{
227  numLocalIds, myIds.getRawPtr(),
228  x_values, y_values, z_values, 1, 1, 1,
229  true, weights.getRawPtr(), 1));
230  }
231  catch (std::exception &e){
232  fail = 1;
233  }
234 
235  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 0a", fail);
236 
237  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
238  myIds.getRawPtr(), xyz_values.getRawPtr(),
239  weights.getRawPtr(), ncoords, nweights);
240 
241  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 0a", fail);
242  }
243 
244  {
246  // 3-dimensional coordinates with stride one and no weights
247 
248  int ncoords = 3;
249  int nweights = 0;
250 
251  std::vector<const zscalar_t *> values, weightValues;
252  std::vector<int> valueStrides, weightStrides;
253 
254  values.push_back(x_values);
255  values.push_back(y_values);
256  values.push_back(z_values);
257  valueStrides.push_back(1);
258  valueStrides.push_back(1);
259  valueStrides.push_back(1);
260 
261  try{
263  numLocalIds, myIds.getRawPtr(), values, valueStrides,
264  weightValues, weightStrides));
265  }
266  catch (std::exception &e){
267  fail = 1;
268  }
269 
270  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 1", fail);
271 
272  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
273  myIds.getRawPtr(), xyz_values.getRawPtr(),
274  weights.getRawPtr(), ncoords, nweights);
275 
276  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 1", fail);
277 
278  // Try using the default: no strides supplied means strides are one.
279 
280  std::vector<int> emptyStrides;
281 
282  try{
284  numLocalIds, myIds.getRawPtr(), values, emptyStrides,
285  weightValues, emptyStrides));
286  }
287  catch (std::exception &e){
288  fail = 1;
289  }
290 
291  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 2", fail);
292 
293  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
294  myIds.getRawPtr(), xyz_values.getRawPtr(),
295  weights.getRawPtr(), ncoords, nweights);
296 
297  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 2", fail);
298  }
299 
300  {
302  // 2-dimensional coordinates with stride three and two weights
303 
304  int ncoords = 2;
305  int nweights = 2;
306 
307  std::vector<const zscalar_t *> values, weightValues;
308  std::vector<int> valueStrides, weightStrides;
309 
310  values.push_back(xyz_values.getRawPtr());
311  values.push_back(xyz_values.getRawPtr() + 1);
312  valueStrides.push_back(3);
313  valueStrides.push_back(3);
314 
315  weightValues.push_back(weights.getRawPtr());
316  weightValues.push_back(weights.getRawPtr() + numLocalIds);
317  weightStrides.push_back(1);
318  weightStrides.push_back(1);
319 
320  try{
322  numLocalIds, myIds.getRawPtr(), values, valueStrides,
323  weightValues, weightStrides));
324  }
325  catch (std::exception &e){
326  fail = 1;
327  }
328 
329  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 3", fail);
330 
331  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
332  myIds.getRawPtr(), xyz_values.getRawPtr(),
333  weights.getRawPtr(), ncoords, nweights);
334 
335  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 3", fail);
336 
337  // Try using default weight strides
338 
339  std::vector<int> emptyStrides;
340 
341  try{
343  numLocalIds, myIds.getRawPtr(), values, valueStrides,
344  weightValues, emptyStrides));
345  }
346  catch (std::exception &e){
347  fail = 1;
348  }
349 
350  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 4", fail);
351 
352  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
353  myIds.getRawPtr(), xyz_values.getRawPtr(),
354  weights.getRawPtr(), ncoords, nweights);
355 
356  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 4", fail);
357  }
358 
359  {
361  // 1-dimensional coordinates with stride one and two weights
362 
363  int ncoords = 1;
364  int nweights = 2;
365 
366  std::vector<const zscalar_t *> values, weightValues;
367  std::vector<int> valueStrides, weightStrides;
368 
369  values.push_back(x_values);
370  valueStrides.push_back(1);
371 
372  weightValues.push_back(weights.getRawPtr());
373  weightValues.push_back(weights.getRawPtr() + numLocalIds);
374  weightStrides.push_back(1);
375  weightStrides.push_back(1);
376 
377  try{
379  numLocalIds, myIds.getRawPtr(), values, valueStrides,
380  weightValues, weightStrides));
381  }
382  catch (std::exception &e){
383  fail = 1;
384  }
385 
386  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 4", fail);
387 
388  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
389  myIds.getRawPtr(), xyz_values.getRawPtr(),
390  weights.getRawPtr(), ncoords, nweights);
391 
392  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 4", fail);
393  }
394 
395  if (rank == 0)
396  std::cout << "PASS" << std::endl;
397 
398  return fail;
399 }
400 
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)
int main(int narg, char *arg[])
double zscalar_t
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
static ArrayRCP< ArrayRCP< zscalar_t > > weights
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.
list idList
Match up parameters to validators.
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...
static const std::string fail
#define TEST_FAIL_AND_RETURN_VALUE(comm, ok, s, rc)
int zgno_t
void getIDsView(const gno_t *&ids) const
Provide a pointer to this process&#39; identifiers.
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.