Zoltan2
Loading...
Searching...
No Matches
CoordinateModel.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// Testing of CoordinateModel
12//
13
14#include "Kokkos_UnorderedMap.hpp"
18
19#include <set>
20#include <bitset>
21
22#include <Teuchos_Comm.hpp>
23#include <Teuchos_CommHelpers.hpp>
24#include <Teuchos_DefaultComm.hpp>
25#include <Teuchos_ArrayView.hpp>
26#include <Teuchos_OrdinalTraits.hpp>
27
28#include <Tpetra_CrsMatrix.hpp>
29
30using Teuchos::RCP;
31using Teuchos::Comm;
32
33void testCoordinateModel(std::string &fname, int nWeights,
34 const RCP<const Comm<int> > &comm,
35 bool nodeZeroHasAll, bool printInfo)
36{
37 int fail = 0, gfail = 0;
38
39 if (printInfo){
40 std::cout << "Test: " << fname << std::endl;
41 std::cout << "Num Weights: " << nWeights;
42 std::cout << " proc 0 has all: " << nodeZeroHasAll;
43 std::cout << std::endl;
44 }
45
47 // Input data
49
50 typedef Tpetra::MultiVector<zscalar_t, zlno_t, zgno_t, znode_t> mv_t;
51
52 RCP<UserInputForTests> uinput;
53
54 try{
55 uinput = rcp(new UserInputForTests(testDataFilePath, fname, comm, true));
56 }
57 catch(std::exception &e){
58 fail=1;
59 }
60
61 TEST_FAIL_AND_EXIT(*comm, !fail, "input constructor", 1);
62
63 RCP<mv_t> coords;
64
65 try{
66 coords = uinput->getUICoordinates();
67 }
68 catch(std::exception &e){
69 fail=2;
70 }
71
72 TEST_FAIL_AND_EXIT(*comm, !fail, "getting coordinates", 1);
73
74 int coordDim = coords->getNumVectors();
75
76 TEST_FAIL_AND_EXIT(*comm, coordDim <= 3, "dim 3 at most", 1);
77
78 const zscalar_t *x=NULL, *y=NULL, *z=NULL;
79
80 x = coords->getData(0).getRawPtr();
81 if (coordDim > 1){
82 y = coords->getData(1).getRawPtr();
83 if (coordDim > 2)
84 z = coords->getData(2).getRawPtr();
85 }
86
87 // Are these coordinates correct
88
89 int nLocalIds = coords->getLocalLength();
90 ArrayView<const zgno_t> idList = coords->getMap()->getLocalElementList();
91
92 int nGlobalIds = 0;
93 if (nodeZeroHasAll){
94 if (comm->getRank() > 0){
95 x = y = z = NULL;
96 nLocalIds = 0;
97 }
98 else{
99 nGlobalIds = nLocalIds;
100 }
101 Teuchos::broadcast<int, int>(*comm, 0, &nGlobalIds);
102 }
103 else{
104 nGlobalIds = coords->getGlobalLength();
105 }
106
107 Array<ArrayRCP<const zscalar_t> > coordWeights(nWeights);
108
109 if (nLocalIds > 0){
110 for (int wdim=0; wdim < nWeights; wdim++){
111 zscalar_t *w = new zscalar_t [nLocalIds];
112 for (int i=0; i < nLocalIds; i++){
113 w[i] = ((i%2) + 1) + wdim;
114 }
115 coordWeights[wdim] = Teuchos::arcp(w, 0, nLocalIds);
116 }
117 }
118
119
121 // Create a BasicVectorAdapter adapter object.
123
125 typedef Zoltan2::VectorAdapter<mv_t> base_ia_t;
126
127 RCP<ia_t> ia;
128
129 if (nWeights == 0){ // use the simpler constructor
130 try{
131 ia = rcp(new ia_t(nLocalIds, idList.getRawPtr(), x, y, z));
132 }
133 catch(std::exception &e){
134 fail=3;
135 }
136 }
137 else{
138 std::vector<const zscalar_t *> values, weights;
139 std::vector<int> valueStrides, weightStrides; // default is 1
140 values.push_back(x);
141 if (y) {
142 values.push_back(y);
143 if (z)
144 values.push_back(z);
145 }
146 for (int wdim=0; wdim < nWeights; wdim++){
147 weights.push_back(coordWeights[wdim].getRawPtr());
148 }
149
150 try{
151 ia = rcp(new ia_t(nLocalIds, idList.getRawPtr(),
152 values, valueStrides, weights, weightStrides));
153 }
154 catch(std::exception &e){
155 fail=4;
156 }
157 }
158
159 RCP<const base_ia_t> base_ia = Teuchos::rcp_dynamic_cast<const base_ia_t>(ia);
160
161 TEST_FAIL_AND_EXIT(*comm, !fail, "making input adapter", 1);
162
164 // Create an CoordinateModel with this input
166
168 typedef std::bitset<Zoltan2::NUM_MODEL_FLAGS> modelFlags_t;
170 modelFlags_t modelFlags;
171
172 RCP<const Zoltan2::Environment> env = rcp(new Zoltan2::Environment(comm));
173 RCP<model_t> model;
174
175
176 try{
177 model = rcp(new model_t(base_ia, env, comm, modelFlags));
178 }
179 catch (std::exception &e){
180 fail = 5;
181 }
182
183 TEST_FAIL_AND_EXIT(*comm, !fail, "making model", 1);
184
185 // Test the CoordinateModel interface
186
187 if (model->getCoordinateDim() != coordDim)
188 fail = 6;
189
190 if (!fail && model->getLocalNumCoordinates() != size_t(nLocalIds))
191 fail = 7;
192
193 if (!fail && model->getGlobalNumCoordinates() != size_t(nGlobalIds))
194 fail = 8;
195
196 if (!fail && model->getNumWeightsPerCoordinate() != nWeights)
197 fail = 9;
198
199 gfail = globalFail(*comm, fail);
200
201 if (gfail)
202 printFailureCode(*comm, fail);
203
204 ArrayView<const zgno_t> gids;
205 ArrayView<input_t> xyz;
206 ArrayView<input_t> wgts;
207
208 model->getCoordinates(gids, xyz, wgts);
209
210 if (!fail && gids.size() != nLocalIds)
211 fail = 10;
212
213 for (int i=0; !fail && i < nLocalIds; i++){
214 if (gids[i] != idList[i])
215 fail = 11;
216 }
217
218 if (!fail && wgts.size() != nWeights)
219 fail = 12;
220
221 const zscalar_t *vals[3] = {x, y, z};
222
223 for (int dim=0; !fail && dim < coordDim; dim++){
224 for (int i=0; !fail && i < nLocalIds; i++){
225 if (xyz[dim][i] != vals[dim][i])
226 fail = 13;
227 }
228 }
229
230 for (int wdim=0; !fail && wdim < nWeights; wdim++){
231 for (int i=0; !fail && i < nLocalIds; i++){
232 if (wgts[wdim][i] != coordWeights[wdim][i])
233 fail = 14;
234 }
235 }
236
237 gfail = globalFail(*comm, fail);
238
239 if (gfail)
240 printFailureCode(*comm, fail);
241
242
246 Kokkos::View<const zgno_t *, typename znode_t::device_type> gidsKokkos;
247
248 Kokkos::View<zscalar_t **, Kokkos::LayoutLeft, typename znode_t::device_type> xyzKokkos;
249 Kokkos::View<zscalar_t **, typename znode_t::device_type> wgtsKokkos;
250
251 model->getCoordinatesKokkos(gidsKokkos, xyzKokkos, wgtsKokkos);
252
253 if (!fail && gidsKokkos.extent(0) != static_cast<size_t>(nLocalIds))
254 fail = 10;
255
256 auto gidsKokkos_host = Kokkos::create_mirror_view(gidsKokkos);
257 Kokkos::deep_copy(gidsKokkos_host, gidsKokkos);
258
259 for (int i=0; !fail && i < nLocalIds; i++){
260 if (gidsKokkos_host(i) != idList[i])
261 fail = 11;
262 }
263
264 if (!fail && wgtsKokkos.extent(1) != static_cast<size_t>(nWeights))
265 fail = 12;
266
267 auto xyzKokkos_host = Kokkos::create_mirror_view(xyzKokkos);
268 Kokkos::deep_copy(xyzKokkos_host, xyzKokkos);
269
270 for (int dim=0; !fail && dim < coordDim; dim++){
271 for (int i=0; !fail && i < nLocalIds; i++){
272 if (xyzKokkos_host(i, dim) != vals[dim][i])
273 fail = 13;
274 }
275 }
276
277 auto wgtsKokkos_host = Kokkos::create_mirror_view(wgtsKokkos);
278 Kokkos::deep_copy(wgtsKokkos_host, wgtsKokkos);
279
280 for (int wdim=0; !fail && wdim < nWeights; wdim++){
281 for (int i=0; !fail && i < nLocalIds; i++){
282 if (wgtsKokkos_host(i, wdim) != coordWeights[wdim][i])
283 fail = 14;
284 }
285 }
286
287 gfail = globalFail(*comm, fail);
288
289 if (gfail)
290 printFailureCode(*comm, fail);
291}
292
293int main(int narg, char *arg[])
294{
295 Tpetra::ScopeGuard tscope(&narg, &arg);
296 Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
297
298 int rank = comm->getRank();
299 string fname("simple"); // reader will seek coord file
300
301 testCoordinateModel(fname, 0, comm, false, rank==0);
302
303 testCoordinateModel(fname, 1, comm, false, rank==0);
304
305 testCoordinateModel(fname, 2, comm, false, rank==0);
306
307 testCoordinateModel(fname, 0, comm, true, rank==0);
308
309 testCoordinateModel(fname, 1, comm, true, rank==0);
310
311 testCoordinateModel(fname, 2, comm, true, rank==0);
312
313 if (rank==0) std::cout << "PASS" << std::endl;
314
315 return 0;
316}
void testCoordinateModel(std::string &fname, int nWeights, const RCP< const Comm< int > > &comm, bool nodeZeroHasAll, bool printInfo)
int globalFail(const Comm< int > &comm, int fail)
void printFailureCode(const Comm< int > &comm, int fail)
#define TEST_FAIL_AND_EXIT(comm, ok, s, code)
Defines the BasicVectorAdapter class.
Defines the CoordinateModel classes.
common code used by tests
float zscalar_t
std::string testDataFilePath(".")
int main()
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
This class provides geometric coordinates with optional weights to the Zoltan2 algorithm.
The user parameters, debug, timing and memory profiling output objects, and error checking methods.
The StridedData class manages lists of weights or coordinates.
VectorAdapter defines the interface for vector input.
static const std::string fail
static ArrayRCP< ArrayRCP< zscalar_t > > weights