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 // Zoltan2: A package of combinatorial algorithms for scientific computing
4 //
5 // Copyright 2012 NTESS and the Zoltan2 contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 //
11 // Test for Zoltan2::BasicVectorAdapter for coordinate-based problems
12 
14 #include <Zoltan2_TestHelpers.hpp>
15 
16 #include <Teuchos_DefaultComm.hpp>
17 #include <Teuchos_RCP.hpp>
18 #include <Teuchos_CommHelpers.hpp>
19 
20 using Teuchos::RCP;
21 using Teuchos::Comm;
22 using Teuchos::Array;
23 
25 
28  int len, int glen, zgno_t *ids,
29  zscalar_t *xyz,
31  int nCoords, int nWeights)
32 {
33  int fail = 0;
34 
35  if (ia->getNumEntriesPerID() != nCoords)
36  fail = 100;
37 
38  if (!fail && ia->getNumWeightsPerID() != nWeights)
39  fail = 101;
40 
41  if (!fail && ia->getLocalNumIDs() != size_t(len))
42  fail = 102;
43 
44  for (int x=0; !fail && x < nCoords; x++){
45  const zgno_t *idList;
46  const zscalar_t *vals;
47  int stride;
48 
49  ia->getIDsView(idList);
50  ia->getEntriesView(vals, stride, x);
51 
52  zscalar_t *coordVal = xyz + x;
53  for (int i=0; !fail && i < len; i++, coordVal += 3){
54 
55  if (idList[i] != ids[i])
56  fail = 105;
57 
58  if (!fail && vals[stride*i] != *coordVal)
59  fail = 106;
60  }
61  }
62 
63  for (int w=0; !fail && w < nWeights; w++){
64  const zscalar_t *wgts;
65  int stride;
66 
67  ia->getWeightsView(wgts, stride, w);
68 
69  zscalar_t *weightVal = weights + len*w;
70  for (int i=0; !fail && i < len; i++, weightVal++){
71  if (wgts[stride*i] != *weightVal)
72  fail = 110;
73  }
74  }
75 
76  return fail;
77 }
78 
79 
80 int main(int narg, char *arg[])
81 {
82  Tpetra::ScopeGuard tscope(&narg, &arg);
83  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
84 
85  int rank = comm->getRank();
86  int fail = 0;
87 
88  // Get some coordinates
89 
90  typedef Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> mv_t;
91  RCP<UserInputForTests> uinput;
92  Teuchos::ParameterList params;
93  params.set("input file", "simple");
94  params.set("file type", "Chaco");
95 
96  try{
97  uinput = rcp(new UserInputForTests(params, comm));
98  }
99  catch(std::exception &e){
100  fail=1;
101  }
102 
103  TEST_FAIL_AND_EXIT(*comm, !fail, "input constructor", 1);
104 
105  RCP<mv_t> coords;
106 
107  try{
108  coords = uinput->getUICoordinates();
109  }
110  catch(std::exception &e){
111  fail=1;
112  }
113 
114  TEST_FAIL_AND_EXIT(*comm, !fail, "getting coordinates", 1);
115 
116  int numLocalIds = coords->getLocalLength();
117  int numGlobalIds = coords->getGlobalLength();
118  int coordDim = coords->getNumVectors();
119  ArrayView<const zgno_t> idList = coords->getMap()->getLocalElementList();
120 
121  // Create global Ids, x-, y- and z-coordinates, and also arrays of weights.
122 
123  Array<zgno_t> myIds(numLocalIds);
124  zgno_t myFirstId = rank * numLocalIds;
125 
126  int wdim = 2;
127  Array<zscalar_t> weights(numLocalIds*wdim);
128  for (int i = 0; i < numLocalIds*wdim; i++) weights[i] = zscalar_t(i);
129 
130  zscalar_t *x_values= coords->getDataNonConst(0).getRawPtr();
131  zscalar_t *y_values= x_values; // fake 3 dimensions if needed
132  zscalar_t *z_values= x_values;
133 
134  if (coordDim > 1){
135  y_values= coords->getDataNonConst(1).getRawPtr();
136  if (coordDim > 2)
137  z_values= coords->getDataNonConst(2).getRawPtr();
138  }
139 
140  Array<zscalar_t> xyz_values(3*numLocalIds);
141 
142  for (zlno_t i=0; i < numLocalIds; i++) // global Ids
143  myIds[i] = myFirstId+i;
144 
145  zscalar_t *x = xyz_values.getRawPtr(); // a stride-3 coordinate array
146  zscalar_t *y = x+1;
147  zscalar_t *z = y+1;
148 
149  for (int i=0, ii=0; i < numLocalIds; i++, ii += 3){
150  x[ii] = x_values[i];
151  y[ii] = y_values[i];
152  z[ii] = z_values[i];
153  }
154 
155  RCP<Zoltan2::BasicVectorAdapter<userTypes_t> > ia;
156 
157  {
159  // 3-dimensional coordinates with stride one and no weights,
160  // using simpler constructor
161 
162  int ncoords = 3;
163  int nweights = 0;
164 
165  try{
167  numLocalIds, myIds.getRawPtr(), x_values, y_values, z_values));
168  }
169  catch (std::exception &e){
170  fail = 1;
171  }
172 
173  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 0", fail);
174 
175  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
176  myIds.getRawPtr(), xyz_values.getRawPtr(),
177  weights.getRawPtr(), ncoords, nweights);
178 
179  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 0", fail);
180  }
181 
182  {
184  // 3-dimensional coordinates with stride one and one weight
185  // using simpler constructor
186 
187  int ncoords = 3;
188  int nweights = 1;
189 
190  try{
192  numLocalIds, myIds.getRawPtr(),
193  x_values, y_values, z_values, 1, 1, 1,
194  true, weights.getRawPtr(), 1));
195  }
196  catch (std::exception &e){
197  fail = 1;
198  }
199 
200  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 0a", fail);
201 
202  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
203  myIds.getRawPtr(), xyz_values.getRawPtr(),
204  weights.getRawPtr(), ncoords, nweights);
205 
206  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 0a", fail);
207  }
208 
209  {
211  // 3-dimensional coordinates with stride one and no weights
212 
213  int ncoords = 3;
214  int nweights = 0;
215 
216  std::vector<const zscalar_t *> values, weightValues;
217  std::vector<int> valueStrides, weightStrides;
218 
219  values.push_back(x_values);
220  values.push_back(y_values);
221  values.push_back(z_values);
222  valueStrides.push_back(1);
223  valueStrides.push_back(1);
224  valueStrides.push_back(1);
225 
226  try{
228  numLocalIds, myIds.getRawPtr(), values, valueStrides,
229  weightValues, weightStrides));
230  }
231  catch (std::exception &e){
232  fail = 1;
233  }
234 
235  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 1", 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 1", fail);
242 
243  // Try using the default: no strides supplied means strides are one.
244 
245  std::vector<int> emptyStrides;
246 
247  try{
249  numLocalIds, myIds.getRawPtr(), values, emptyStrides,
250  weightValues, emptyStrides));
251  }
252  catch (std::exception &e){
253  fail = 1;
254  }
255 
256  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 2", fail);
257 
258  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
259  myIds.getRawPtr(), xyz_values.getRawPtr(),
260  weights.getRawPtr(), ncoords, nweights);
261 
262  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 2", fail);
263  }
264 
265  {
267  // 2-dimensional coordinates with stride three and two weights
268 
269  int ncoords = 2;
270  int nweights = 2;
271 
272  std::vector<const zscalar_t *> values, weightValues;
273  std::vector<int> valueStrides, weightStrides;
274 
275  values.push_back(xyz_values.getRawPtr());
276  values.push_back(xyz_values.getRawPtr() + 1);
277  valueStrides.push_back(3);
278  valueStrides.push_back(3);
279 
280  weightValues.push_back(weights.getRawPtr());
281  weightValues.push_back(weights.getRawPtr() + numLocalIds);
282  weightStrides.push_back(1);
283  weightStrides.push_back(1);
284 
285  try{
287  numLocalIds, myIds.getRawPtr(), values, valueStrides,
288  weightValues, weightStrides));
289  }
290  catch (std::exception &e){
291  fail = 1;
292  }
293 
294  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 3", fail);
295 
296  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
297  myIds.getRawPtr(), xyz_values.getRawPtr(),
298  weights.getRawPtr(), ncoords, nweights);
299 
300  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 3", fail);
301 
302  // Try using default weight strides
303 
304  std::vector<int> emptyStrides;
305 
306  try{
308  numLocalIds, myIds.getRawPtr(), values, valueStrides,
309  weightValues, emptyStrides));
310  }
311  catch (std::exception &e){
312  fail = 1;
313  }
314 
315  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 4", fail);
316 
317  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
318  myIds.getRawPtr(), xyz_values.getRawPtr(),
319  weights.getRawPtr(), ncoords, nweights);
320 
321  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 4", fail);
322  }
323 
324  {
326  // 1-dimensional coordinates with stride one and two weights
327 
328  int ncoords = 1;
329  int nweights = 2;
330 
331  std::vector<const zscalar_t *> values, weightValues;
332  std::vector<int> valueStrides, weightStrides;
333 
334  values.push_back(x_values);
335  valueStrides.push_back(1);
336 
337  weightValues.push_back(weights.getRawPtr());
338  weightValues.push_back(weights.getRawPtr() + numLocalIds);
339  weightStrides.push_back(1);
340  weightStrides.push_back(1);
341 
342  try{
344  numLocalIds, myIds.getRawPtr(), values, valueStrides,
345  weightValues, weightStrides));
346  }
347  catch (std::exception &e){
348  fail = 1;
349  }
350 
351  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 4", fail);
352 
353  fail = checkBasicCoordinate(ia.getRawPtr(), numLocalIds, numGlobalIds,
354  myIds.getRawPtr(), xyz_values.getRawPtr(),
355  weights.getRawPtr(), ncoords, nweights);
356 
357  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check adapter 4", fail);
358  }
359 
360  if (rank == 0)
361  std::cout << "PASS" << std::endl;
362 
363  return fail;
364 }
365 
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > userTypes_t
size_t getLocalNumIDs() const
Returns the number of objects on this process.
static ArrayRCP< ArrayRCP< zscalar_t > > weights
int checkBasicCoordinate(Zoltan2::BasicVectorAdapter< userTypes_t > *ia, int len, int glen, zgno_t *ids, zscalar_t *xyz, zscalar_t *weights, int nCoords, int nWeights)
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
int main(int narg, char **arg)
Definition: coloring1.cpp:164
common code used by tests
int getNumEntriesPerID() const
Return the number of vectors.
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
Provide pointer to a weight array with stride.
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
Tpetra::Map::local_ordinal_type zlno_t
static const std::string fail
#define TEST_FAIL_AND_RETURN_VALUE(comm, ok, s, rc)
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...
float zscalar_t
void getEntriesView(const scalar_t *&entries, int &stride, int idx=0) const
Provide a pointer to the elements of the specified vector.
Defines the BasicVectorAdapter class.
Tpetra::Map::global_ordinal_type zgno_t