Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
BasicVectorInput.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
47 
49 #include <Zoltan2_TestHelpers.hpp>
50 
51 #include <Teuchos_DefaultComm.hpp>
52 #include <Teuchos_RCP.hpp>
53 #include <Teuchos_CommHelpers.hpp>
54 
55 
57 
59  Zoltan2::BasicVectorAdapter<userTypes_t> *ia, int len, int glen,
60  zgno_t *ids, int mvdim, const zscalar_t **values, int *valueStrides,
61  int wdim, const zscalar_t **weights, int *weightStrides)
62 {
63  int fail = 0;
64  bool strideOne = false;
65 
66  if (valueStrides == NULL) strideOne = true;
67 
68  if (ia->getNumEntriesPerID() != mvdim)
69  fail = 100;
70 
71  if (!fail && ia->getNumWeightsPerID() != wdim)
72  fail = 101;
73 
74  if (!fail && ia->getLocalNumIDs() != size_t(len))
75  fail = 102;
76 
77  const zgno_t *idList;
78  ia->getIDsView(idList);
79  for (int i=0; !fail && i < len; i++)
80  if (!fail && idList[i] != ids[i])
81  fail = 107;
82 
83  for (int v=0; !fail && v < mvdim; v++){
84  const zscalar_t *vals;
85  int correctStride = (strideOne ? 1 : valueStrides[v]);
86  int stride;
87 
88  ia->getEntriesView(vals, stride, v);
89 
90  if (!fail && stride != correctStride)
91  fail = 105;
92 
93  for (int i=0; !fail && i < len; i++){
94 // TODO fix values check
95 // if (vals[stride*i] != values[v][correctStride*i])
96 // fail = 106;
97  }
98  }
99 
100  for (int w=0; !fail && w < wdim; w++){
101  const zscalar_t *wgts;
102  int stride;
103 
104  ia->getWeightsView(wgts, stride, w);
105 
106  if (!fail && stride != weightStrides[w])
107  fail = 109;
108 
109  for (int i=0; !fail && i < len; i++){
110  if (wgts[stride*i] != weights[w][weightStrides[w]*i])
111  fail = 110;
112  }
113  }
114 
115  return fail;
116 }
117 
118 
119 int main(int narg, char *arg[])
120 {
121  Tpetra::ScopeGuard tscope(&narg, &arg);
122  Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
123 
124  int rank = comm->getRank();
125  int nprocs = comm->getSize();
126  int fail = 0;
127 
128  // Create a single vector and a strided multi-vector with
129  // strided multi-weights.
130 
131  zlno_t numLocalIds = 10;
132  zgno_t *myIds = new zgno_t[numLocalIds];
133  zgno_t myFirstId = rank * numLocalIds;
134 
135  int wdim = 2;
136  zscalar_t *weights = new zscalar_t [numLocalIds*wdim];
137  int *weightStrides = new int [wdim];
138  const zscalar_t **weightPtrs = new const zscalar_t * [wdim];
139 
140  int mvdim = 3;
141  zscalar_t *v_values= new zscalar_t [numLocalIds];
142  zscalar_t *mv_values= new zscalar_t [mvdim * numLocalIds];
143  int *valueStrides = new int [mvdim];
144  const zscalar_t **valuePtrs = new const zscalar_t * [mvdim];
145 
146  for (zlno_t i=0; i < numLocalIds; i++){
147  myIds[i] = myFirstId+i;
148 
149  for (int w=0; w < wdim; w++)
150  weights[w*numLocalIds + i] = w + 1 + nprocs - rank;
151 
152  v_values[i] = numLocalIds-i;
153 
154  for (int v=0; v < mvdim; v++)
155  mv_values[i*mvdim + v] = (v+1) * (nprocs-rank) / (i+1);
156  }
157 
158  for (int w=0; w < wdim; w++){
159  weightStrides[w] = 1;
160  weightPtrs[w] = weights + numLocalIds*w;
161  }
162 
163  for (int v=0; v < mvdim; v++){
164  valueStrides[v] = mvdim;
165  valuePtrs[v] = mv_values + v;
166  }
167 
169 
170  {
171  // A Zoltan2::BasicVectorAdapter object with one vector and no weights
172 
173  std::vector<const zscalar_t *> weightValues;
174  std::vector<int> strides;
175 
176  try{
177  ia = new Zoltan2::BasicVectorAdapter<userTypes_t>(numLocalIds, myIds,
178  v_values, 1);
179  }
180  catch (std::exception &e){
181  fail = 1;
182  }
183 
184  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 1", fail);
185 
186  fail = checkBasicVector(ia, numLocalIds, numLocalIds*nprocs,
187  myIds, 1, valuePtrs, NULL, 0, NULL, NULL);
188 
189  delete ia;
190 
191  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check vector 1", fail);
192  }
193 
194  {
195  // A Zoltan2::BasicVectorAdapter object with one vector and weights
196 
197  std::vector<const zscalar_t *> weightValues;
198  std::vector<int> strides;
199 
200  weightValues.push_back(weightPtrs[0]);
201  weightValues.push_back(weightPtrs[1]);
202  strides.push_back(1);
203  strides.push_back(1);
204 
205  try{
206  ia = new Zoltan2::BasicVectorAdapter<userTypes_t>(numLocalIds, myIds,
207  v_values, 1, true, weightPtrs[0], 1);
208  }
209  catch (std::exception &e){
210  fail = 1;
211  }
212 
213  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 2", fail);
214 
215  fail = checkBasicVector(ia, numLocalIds, numLocalIds*nprocs,
216  myIds, 1, valuePtrs, NULL, 1, weightPtrs, weightStrides);
217 
218  delete ia;
219 
220  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check vector 2", fail);
221  }
222 
223  {
224  // A Zoltan2::BasicVectorAdapter object with a multivector and no weights
225 
226  std::vector<const zscalar_t *> weightValues, values;
227  std::vector<int> wstrides, vstrides;
228 
229  for (int dim=0; dim < mvdim; dim++){
230  values.push_back(valuePtrs[dim]);
231  vstrides.push_back(mvdim);
232  }
233 
234 
235  try{
237  numLocalIds, myIds, values, vstrides, weightValues, wstrides);
238  }
239  catch (std::exception &e){
240  fail = 1;
241  }
242 
243  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 3", fail);
244 
245  fail = checkBasicVector(ia, numLocalIds, numLocalIds*nprocs,
246  myIds, mvdim, valuePtrs, valueStrides, 0, NULL, NULL);
247 
248  delete ia;
249 
250  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check vector 3", fail);
251  }
252 
253  {
254  // A Zoltan2::BasicVectorAdapter object with a multivector with weights
255 
256  std::vector<const zscalar_t *> weightValues, values;
257  std::vector<int> wstrides, vstrides;
258 
259  for (int dim=0; dim < wdim; dim++){
260  weightValues.push_back(weightPtrs[dim]);
261  wstrides.push_back(1);
262  }
263 
264  for (int dim=0; dim < mvdim; dim++){
265  values.push_back(valuePtrs[dim]);
266  vstrides.push_back(mvdim);
267  }
268 
269  try{
271  numLocalIds, myIds, values, vstrides, weightValues, wstrides);
272 
273  }
274  catch (std::exception &e){
275  fail = 1;
276  }
277 
278  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 4", fail);
279 
280  fail = checkBasicVector(ia, numLocalIds, numLocalIds*nprocs,
281  myIds, mvdim, valuePtrs, valueStrides,
282  wdim, weightPtrs, weightStrides);
283 
284  delete ia;
285 
286  TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check vector 4", fail);
287  }
288 
289  if (rank == 0)
290  std::cout << "PASS" << std::endl;
291 
292  delete [] myIds;
293  delete [] weights;
294  delete [] weightStrides;
295  delete [] weightPtrs;
296  delete [] v_values;
297  delete [] mv_values;
298  delete [] valueStrides;
299  delete [] valuePtrs;
300 
301  return fail;
302 }
303 
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 main(int narg, char **arg)
Definition: coloring1.cpp:199
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.
int checkBasicVector(Zoltan2::BasicVectorAdapter< userTypes_t > *ia, int len, int glen, zgno_t *ids, int mvdim, const zscalar_t **values, int *valueStrides, int wdim, const zscalar_t **weights, int *weightStrides)
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