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