Zoltan2
Loading...
Searching...
No Matches
BasicVectorInput.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
12
15
16#include <Teuchos_DefaultComm.hpp>
17#include <Teuchos_RCP.hpp>
18#include <Teuchos_CommHelpers.hpp>
19
20
22
24 Zoltan2::BasicVectorAdapter<userTypes_t> *ia, int len, int glen,
25 zgno_t *ids, int mvdim, const zscalar_t **values, int *valueStrides,
26 int wdim, const zscalar_t **weights, int *weightStrides)
27{
28 int fail = 0;
29 bool strideOne = false;
30
31 if (valueStrides == NULL) strideOne = true;
32
33 if (ia->getNumEntriesPerID() != mvdim)
34 fail = 100;
35
36 if (!fail && ia->getNumWeightsPerID() != wdim)
37 fail = 101;
38
39 if (!fail && ia->getLocalNumIDs() != size_t(len))
40 fail = 102;
41
42 const zgno_t *idList;
43 ia->getIDsView(idList);
44 for (int i=0; !fail && i < len; i++)
45 if (!fail && idList[i] != ids[i])
46 fail = 107;
47
48 for (int v=0; !fail && v < mvdim; v++){
49 const zscalar_t *vals;
50 int correctStride = (strideOne ? 1 : valueStrides[v]);
51 int stride;
52
53 ia->getEntriesView(vals, stride, v);
54
55 if (!fail && stride != correctStride)
56 fail = 105;
57
58 for (int i=0; !fail && i < len; i++){
59// TODO fix values check
60// if (vals[stride*i] != values[v][correctStride*i])
61// fail = 106;
62 }
63 }
64
65 for (int w=0; !fail && w < wdim; w++){
66 const zscalar_t *wgts;
67 int stride;
68
69 ia->getWeightsView(wgts, stride, w);
70
71 if (!fail && stride != weightStrides[w])
72 fail = 109;
73
74 for (int i=0; !fail && i < len; i++){
75 if (wgts[stride*i] != weights[w][weightStrides[w]*i])
76 fail = 110;
77 }
78 }
79
80 return fail;
81}
82
83
84int main(int narg, char *arg[])
85{
86 Tpetra::ScopeGuard tscope(&narg, &arg);
87 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
88
89 int rank = comm->getRank();
90 int nprocs = comm->getSize();
91 int fail = 0;
92
93 // Create a single vector and a strided multi-vector with
94 // strided multi-weights.
95
96 zlno_t numLocalIds = 10;
97 zgno_t *myIds = new zgno_t[numLocalIds];
98 zgno_t myFirstId = rank * numLocalIds;
99
100 int wdim = 2;
101 zscalar_t *weights = new zscalar_t [numLocalIds*wdim];
102 int *weightStrides = new int [wdim];
103 const zscalar_t **weightPtrs = new const zscalar_t * [wdim];
104
105 int mvdim = 3;
106 zscalar_t *v_values= new zscalar_t [numLocalIds];
107 zscalar_t *mv_values= new zscalar_t [mvdim * numLocalIds];
108 int *valueStrides = new int [mvdim];
109 const zscalar_t **valuePtrs = new const zscalar_t * [mvdim];
110
111 for (zlno_t i=0; i < numLocalIds; i++){
112 myIds[i] = myFirstId+i;
113
114 for (int w=0; w < wdim; w++)
115 weights[w*numLocalIds + i] = w + 1 + nprocs - rank;
116
117 v_values[i] = numLocalIds-i;
118
119 for (int v=0; v < mvdim; v++)
120 mv_values[i*mvdim + v] = (v+1) * (nprocs-rank) / (i+1);
121 }
122
123 for (int w=0; w < wdim; w++){
124 weightStrides[w] = 1;
125 weightPtrs[w] = weights + numLocalIds*w;
126 }
127
128 for (int v=0; v < mvdim; v++){
129 valueStrides[v] = mvdim;
130 valuePtrs[v] = mv_values + v;
131 }
132
134
135 {
136 // A Zoltan2::BasicVectorAdapter object with one vector and no weights
137
138 std::vector<const zscalar_t *> weightValues;
139 std::vector<int> strides;
140
141 try{
142 ia = new Zoltan2::BasicVectorAdapter<userTypes_t>(numLocalIds, myIds,
143 v_values, 1);
144 }
145 catch (std::exception &e){
146 fail = 1;
147 }
148
149 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 1", fail);
150
151 fail = checkBasicVector(ia, numLocalIds, numLocalIds*nprocs,
152 myIds, 1, valuePtrs, NULL, 0, NULL, NULL);
153
154 delete ia;
155
156 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check vector 1", fail);
157 }
158
159 {
160 // A Zoltan2::BasicVectorAdapter object with one vector and weights
161
162 std::vector<const zscalar_t *> weightValues;
163 std::vector<int> strides;
164
165 weightValues.push_back(weightPtrs[0]);
166 weightValues.push_back(weightPtrs[1]);
167 strides.push_back(1);
168 strides.push_back(1);
169
170 try{
171 ia = new Zoltan2::BasicVectorAdapter<userTypes_t>(numLocalIds, myIds,
172 v_values, 1, true, weightPtrs[0], 1);
173 }
174 catch (std::exception &e){
175 fail = 1;
176 }
177
178 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 2", fail);
179
180 fail = checkBasicVector(ia, numLocalIds, numLocalIds*nprocs,
181 myIds, 1, valuePtrs, NULL, 1, weightPtrs, weightStrides);
182
183 delete ia;
184
185 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check vector 2", fail);
186 }
187
188 {
189 // A Zoltan2::BasicVectorAdapter object with a multivector and no weights
190
191 std::vector<const zscalar_t *> weightValues, values;
192 std::vector<int> wstrides, vstrides;
193
194 for (int dim=0; dim < mvdim; dim++){
195 values.push_back(valuePtrs[dim]);
196 vstrides.push_back(mvdim);
197 }
198
199
200 try{
202 numLocalIds, myIds, values, vstrides, weightValues, wstrides);
203 }
204 catch (std::exception &e){
205 fail = 1;
206 }
207
208 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 3", fail);
209
210 fail = checkBasicVector(ia, numLocalIds, numLocalIds*nprocs,
211 myIds, mvdim, valuePtrs, valueStrides, 0, NULL, NULL);
212
213 delete ia;
214
215 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check vector 3", fail);
216 }
217
218 {
219 // A Zoltan2::BasicVectorAdapter object with a multivector with weights
220
221 std::vector<const zscalar_t *> weightValues, values;
222 std::vector<int> wstrides, vstrides;
223
224 for (int dim=0; dim < wdim; dim++){
225 weightValues.push_back(weightPtrs[dim]);
226 wstrides.push_back(1);
227 }
228
229 for (int dim=0; dim < mvdim; dim++){
230 values.push_back(valuePtrs[dim]);
231 vstrides.push_back(mvdim);
232 }
233
234 try{
236 numLocalIds, myIds, values, vstrides, weightValues, wstrides);
237
238 }
239 catch (std::exception &e){
240 fail = 1;
241 }
242
243 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "constructor 4", fail);
244
245 fail = checkBasicVector(ia, numLocalIds, numLocalIds*nprocs,
246 myIds, mvdim, valuePtrs, valueStrides,
247 wdim, weightPtrs, weightStrides);
248
249 delete ia;
250
251 TEST_FAIL_AND_RETURN_VALUE(*comm, fail==0, "check vector 4", fail);
252 }
253
254 if (rank == 0)
255 std::cout << "PASS" << std::endl;
256
257 delete [] myIds;
258 delete [] weights;
259 delete [] weightStrides;
260 delete [] weightPtrs;
261 delete [] v_values;
262 delete [] mv_values;
263 delete [] valueStrides;
264 delete [] valuePtrs;
265
266 return fail;
267}
268
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)
#define TEST_FAIL_AND_RETURN_VALUE(comm, ok, s, rc)
Defines the BasicVectorAdapter class.
common code used by tests
float zscalar_t
Tpetra::Map ::local_ordinal_type zlno_t
Tpetra::Map ::global_ordinal_type zgno_t
int main()
A simple class that can be the User template argument for an InputAdapter.
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
size_t getLocalNumIDs() const
Returns the number of objects on this process.
void getEntriesView(const scalar_t *&entries, int &stride, int idx=0) const
Provide a pointer to the elements of the specified vector.
int getNumWeightsPerID() const
Returns the number of weights per object. Number of weights per object should be zero or greater....
void getWeightsView(const scalar_t *&weights, int &stride, int idx) const
void getIDsView(const gno_t *&ids) const
int getNumEntriesPerID() const
Return the number of vectors.
static const std::string fail
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > userTypes_t
static ArrayRCP< ArrayRCP< zscalar_t > > weights