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