Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
fix4785.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 
55 #include <Zoltan2_InputTraits.hpp>
56 #include <Tpetra_Map.hpp>
57 #include <vector>
58 #include <cstdlib>
59 
60 typedef Tpetra::Map<> myMap_t;
61 typedef myMap_t::local_ordinal_type myLocalId_t;
62 typedef myMap_t::global_ordinal_type myGlobalId_t;
63 typedef double myScalar_t;
64 
66 template <typename A, typename B>
67 void copyArrays(size_t n, A *&a, B *b)
68 {
69  a = new A[n];
70  for (size_t i = 0; i < n; i++) a[i] = b[i];
71 }
72 
73 template <typename A, typename B>
74 void copyArrays(size_t n, A *&a, A *b)
75 {
76  a = b;
77 }
78 
80 template <typename A, typename B>
81 void freeArrays(A *&a, B *b)
82 {
83  delete [] a;
84 }
85 
86 template <typename A, typename B>
87 void freeArrays(A *&a, A *b)
88 {
89  // no delete needed since only copied pointer
90 }
91 
94  const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
95  Teuchos::ParameterList &params,
96  size_t localCount,
97  myGlobalId_t *globalIds,
98  myScalar_t *coords,
99  int &nFail
100 )
101 {
103  typedef Zoltan2::BasicVectorAdapter<myTypes> inputAdapter_t;
105 
106  myScalar_t *sx = coords;
107  myScalar_t *sy = sx + localCount;
108  myScalar_t *sz = sy + localCount;
109 
110  inputAdapter_t *ia = new inputAdapter_t(localCount,globalIds,sx,sy,sz,1,1,1);
111 
114 
115  problem->solve();
116 
117  quality_t *metricObject = new quality_t(ia, &params, comm,
118  &problem->getSolution());
119  if (comm->getRank() == 0){
120 
121  metricObject->printMetrics(std::cout);
122 
123  double imb = metricObject->getObjectCountImbalance();
124  if (imb <= 1.01) // Should get perfect balance
125  std::cout << "no weights -- balance satisfied: " << imb << std::endl;
126  else {
127  std::cout << "no weights -- balance failure: " << imb << std::endl;
128  nFail++;
129  }
130  std::cout << std::endl;
131  }
132 
133  delete metricObject;
134  delete problem;
135  delete ia;
136 }
137 
140  const Teuchos::RCP<const Teuchos::Comm<int> > &comm,
141  Teuchos::ParameterList &params,
142  size_t localCount,
143  myGlobalId_t *globalIds,
144  myScalar_t *coords,
146  int &nFail
147 )
148 {
150  typedef Zoltan2::BasicVectorAdapter<myTypes> inputAdapter_t;
152 
153  std::vector<const myScalar_t *> coordVec(3);
154  std::vector<int> coordStrides(3);
155 
156  coordVec[0] = coords; coordStrides[0] = 1;
157  coordVec[1] = coords + localCount; coordStrides[1] = 1;
158  coordVec[2] = coords + localCount + localCount; coordStrides[2] = 1;
159 
160  std::vector<const myScalar_t *> weightVec(1);
161  std::vector<int> weightStrides(1);
162 
163  weightVec[0] = weights; weightStrides[0] = 1;
164 
165  inputAdapter_t *ia=new inputAdapter_t(localCount, globalIds, coordVec,
166  coordStrides,weightVec,weightStrides);
167 
170 
171  problem->solve();
172 
173  quality_t *metricObject = new quality_t(ia, &params, comm,
174  &problem->getSolution());
175  if (comm->getRank() == 0){
176 
177  metricObject->printMetrics(std::cout);
178 
179  double imb = metricObject->getWeightImbalance(0);
180  if (imb <= 1.01)
181  std::cout << "weighted -- balance satisfied " << imb << std::endl;
182  else {
183  std::cout << "weighted -- balance failed " << imb << std::endl;
184  nFail++;
185  }
186  std::cout << std::endl;
187  }
188 
189  delete metricObject;
190  delete problem;
191  delete ia;
192 }
193 
195 int main(int narg, char *arg[])
196 {
197  Tpetra::ScopeGuard scope(&narg, &arg);
198  const Teuchos::RCP<const Teuchos::Comm<int> > comm = Tpetra::getDefaultComm();
199  int me = comm->getRank();
200  int np = comm->getSize();
201  int nFail = 0;
202 
204  // Create parameters for an MJ problem
205 
206  Teuchos::ParameterList params("test params");
207  params.set("debug_level", "basic_status");
208  params.set("error_check_level", "debug_mode_assertions");
209 
210  params.set("algorithm", "multijagged");
211  params.set("num_global_parts", 64);
212 
214  // Create input data.
215 
216  const size_t minN = 9000000;
217  const size_t maxN = 16000000;
218  const size_t incN = 2000000;
219  size_t maxLocalCount = maxN / np; // biggest test we'll run
220 
221  // Create coordinates that range from 0 to 999
222  const int dim = 3;
223  myScalar_t *coords = new myScalar_t[dim * maxLocalCount];
224 
225  srand(me);
226  for (size_t i=0; i < maxLocalCount*dim; i++)
227  coords[i] = myScalar_t(0);
228 
229  // Create weights for the coordinates
230  myScalar_t *weights = new myScalar_t[maxLocalCount];
231  for (size_t i=0; i < maxLocalCount; i++)
232  weights[i] = myScalar_t(1+me);
233 
234  // Allocate space for global IDs; they will be generated below
235  myGlobalId_t *globalIds = new myGlobalId_t[maxLocalCount];
236 
237  for (size_t N = minN; N < maxN; N += incN) {
238 
239  size_t localCount = N / np;
240 
241  // Create consecutive global ids for the coordinates
242  myGlobalId_t offset = me * localCount;
243  for (size_t i=0; i < localCount; i++)
244  globalIds[i] = offset++;
245 
246  if (me == 0) {
247  std::cout << "---------------------------------------------" << std::endl;
248  std::cout << "myGlobalId_t = " << typeid(offset).name()
249  << " " << sizeof(offset)
250  << "; localCount = " << localCount
251  << "; globalCount = " << np * localCount << std::endl;
252  }
253 
255  // Test one: No weights
256  if (me == 0) std::cout << "Test: no weights, scalar = double" << std::endl;
257  test_no_weights(comm, params, localCount, globalIds, coords, nFail);
258 
260  // Test two: weighted
261  if (me == 0) std::cout << "Test: weights, scalar = double" << std::endl;
262  test_weights(comm, params, localCount, globalIds, coords, weights, nFail);
263 
264  // Early exit when failure is detected
265  int gnFail;
266  Teuchos::reduceAll<int, int>(*comm, Teuchos::REDUCE_MAX, 1, &nFail, &gnFail);
267  if (gnFail > 0) break;
268 
269  }
270 
271  delete [] weights;
272  delete [] coords;
273  delete [] globalIds;
274 
275  if (me == 0) {
276  if (nFail == 0) std::cout << "PASS" << std::endl;
277  else std::cout << "FAIL: " << nFail << " tests failed" << std::endl;
278  }
279 
280  return 0;
281 }
282 
myMap_t::local_ordinal_type myLocalId_t
Definition: fix4785.cpp:61
Tpetra::Map myMap_t
Definition: fix4785.cpp:60
int main(int narg, char *arg[])
A simple class that can be the User template argument for an InputAdapter.
static ArrayRCP< ArrayRCP< zscalar_t > > weights
Defines the PartitioningSolution class.
void test_no_weights(const Teuchos::RCP< const Teuchos::Comm< int > > &comm, Teuchos::ParameterList &params, size_t localCount, myGlobalId_t *globalIds, myScalar_t *coords, int &nFail)
Definition: fix4785.cpp:93
void test_weights(const Teuchos::RCP< const Teuchos::Comm< int > > &comm, Teuchos::ParameterList &params, size_t localCount, myGlobalId_t *globalIds, myScalar_t *coords, myScalar_t *weights, int &nFail)
Definition: fix4785.cpp:139
void copyArrays(size_t n, A *&a, B *b)
Definition: fix4785.cpp:67
myMap_t::global_ordinal_type myGlobalId_t
Definition: fix4785.cpp:62
double myScalar_t
Definition: fix4785.cpp:63
BasicVectorAdapter represents a vector (plus optional weights) supplied by the user as pointers to st...
Traits for application input objects.
void freeArrays(A *&a, B *b)
Definition: fix4785.cpp:81
const PartitioningSolution< Adapter > & getSolution()
Get the solution to the problem.
PartitioningProblem sets up partitioning problems for the user.
Defines the PartitioningProblem class.
A class that computes and returns quality metrics.
Defines the BasicVectorAdapter class.
void solve(bool updateInputData=true)
Direct the problem to create a solution.