Zoltan2
Zoltan2_AlgZoltan.hpp
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 #ifndef _ZOLTAN2_ALGZOLTAN_HPP_
46 #define _ZOLTAN2_ALGZOLTAN_HPP_
47 
48 #include <Zoltan2_Standards.hpp>
49 #include <Zoltan2_Algorithm.hpp>
51 #include <Zoltan2_Util.hpp>
52 #include <Zoltan2_TPLTraits.hpp>
53 
54 #include <Zoltan2_Model.hpp>
55 
57 #include <zoltan_cpp.h>
58 
62 //
63 // This first design templates Zoltan's callback functions on the
64 // input adapter. This approach has the advantage of simplicity and
65 // is most similar to current usage of Zoltan (where the callbacks define
66 // the model).
67 // A better approach might template them on a model,
68 // allowing Zoltan2 greater flexibility in creating models from the input.
69 // Alternatively, different callback implementations could be provided to
70 // represent different models to Zoltan.
72 
73 namespace Zoltan2 {
74 
75 template <typename Adapter>
76 class AlgZoltan : public Algorithm<Adapter>
77 {
78 
79 private:
80 
81  typedef typename Adapter::lno_t lno_t;
82  typedef typename Adapter::gno_t gno_t;
83  typedef typename Adapter::scalar_t scalar_t;
84  typedef typename Adapter::part_t part_t;
85  typedef typename Adapter::user_t user_t;
86  typedef typename Adapter::userCoord_t userCoord_t;
87 
88  const RCP<const Environment> env;
89  const RCP<const Comm<int> > problemComm;
90  const RCP<const typename Adapter::base_adapter_t> adapter;
91  RCP<const Model<Adapter> > model;
92  RCP<Zoltan> zz;
93 
94  MPI_Comm mpicomm;
95 
96  void setMPIComm(const RCP<const Comm<int> > &problemComm__) {
97 # ifdef HAVE_ZOLTAN2_MPI
98  mpicomm = TeuchosConst2MPI(problemComm__);
99 # else
100  mpicomm = MPI_COMM_WORLD; // taken from siMPI
101 # endif
102  }
103 
104  void zoltanInit() {
105  // call Zoltan_Initialize to make sure MPI_Init is called (in MPI or siMPI).
106  int argc = 0;
107  char **argv = NULL;
108  float ver;
109  Zoltan_Initialize(argc, argv, &ver);
110  }
111 
112  void setCallbacksIDs()
113  {
114  zz->Set_Num_Obj_Fn(zoltanNumObj<Adapter>, (void *) &(*adapter));
115  zz->Set_Obj_List_Fn(zoltanObjList<Adapter>, (void *) &(*adapter));
116 
117  const part_t *myparts;
118  adapter->getPartsView(myparts);
119  if (myparts != NULL)
120  zz->Set_Part_Multi_Fn(zoltanParts<Adapter>, (void *) &(*adapter));
121  }
122 
123  template <typename AdapterWithCoords>
124  void setCallbacksGeom(const AdapterWithCoords *ia)
125  {
126  // Coordinates may be provided by the MeshAdapter or VectorAdapter.
127  // VectorAdapter may be provided directly by user or indirectly through
128  // GraphAdapter or MatrixAdapter. So separate template type is needed.
129  zz->Set_Num_Geom_Fn(zoltanNumGeom<AdapterWithCoords>, (void *) ia);
130  zz->Set_Geom_Multi_Fn(zoltanGeom<AdapterWithCoords>, (void *) ia);
131  }
132 
133  void setCallbacksGraph(
134  const RCP<const GraphAdapter<user_t,userCoord_t> > &adp)
135  {
136  std::cout << "NotReadyForGraphYet" << std::endl;
137  // TODO
138  }
139 
140  void setCallbacksGraph(
141  const RCP<const MatrixAdapter<user_t,userCoord_t> > &adp)
142  {
143  std::cout << "NotReadyForGraphYet" << std::endl;
144  // TODO
145  }
146 
147  void setCallbacksGraph(
148  const RCP<const MeshAdapter<user_t> > &adp)
149  {
150  std::cout << "NotReadyForGraphYet" << std::endl;
151  // TODO
152  }
153 
154  void setCallbacksHypergraph(
155  const RCP<const MatrixAdapter<user_t,userCoord_t> > &adp)
156  {
157  // TODO: If add parameter list to this function, can register
158  // TODO: different callbacks depending on the hypergraph model to use
159 
160  zz->Set_HG_Size_CS_Fn(zoltanHGSizeCSForMatrixAdapter<Adapter>,
161  (void *) &(*adp));
162  zz->Set_HG_CS_Fn(zoltanHGCSForMatrixAdapter<Adapter>,
163  (void *) &(*adp));
164 
165  // zz->Set_HG_Size_Edge_Wts_Fn(zoltanHGSizeEdgeWtsForMatrixAdapter<Adapter>,
166  // (void *) &(*adapter));
167  // zz->Set_HG_Edge_Wts_Fn(zoltanHGSizeEdgeWtsForMatrixAdapter<Adapter>,
168  // (void *) &(*adapter));
169  }
170 
171  void setCallbacksHypergraph(const RCP<const MeshAdapter<user_t> > &adp)
172  {
173 
174  const Teuchos::ParameterList &pl = env->getParameters();
175 
176  const Teuchos::ParameterEntry *pe = pl.getEntryPtr("hypergraph_model_type");
177  std::string model_type("traditional");
178  if (pe){
179  model_type = pe->getValue<std::string>(&model_type);
180  }
181 
182  if (model_type=="ghosting" ||
183  !adp->areEntityIDsUnique(adp->getPrimaryEntityType())) {
184  Zoltan2::modelFlag_t flags;
186  problemComm, flags,
188  model = rcp(static_cast<const Model<Adapter>* >(mdl),true);
189 
190  zz->Set_Num_Obj_Fn(zoltanHGModelNumObj<Adapter>, (void *) &(*mdl));
191  zz->Set_Obj_List_Fn(zoltanHGModelObjList<Adapter>, (void *) &(*mdl));
192 
193  zz->Set_HG_Size_CS_Fn(zoltanHGModelSizeCSForMeshAdapter<Adapter>,
194  (void *) &(*mdl));
195  zz->Set_HG_CS_Fn(zoltanHGModelCSForMeshAdapter<Adapter>,
196  (void *) &(*mdl));
197  }
198  else {
199  //If entities are unique we dont need the extra cost of the model
200  zz->Set_HG_Size_CS_Fn(zoltanHGSizeCSForMeshAdapter<Adapter>,
201  (void *) &(*adp));
202  zz->Set_HG_CS_Fn(zoltanHGCSForMeshAdapter<Adapter>,
203  (void *) &(*adp));
204  }
205  // zz->Set_HG_Size_Edge_Wts_Fn(zoltanHGSizeEdgeWtsForMeshAdapter<Adapter>,
206  // (void *) &(*adp));
207  // zz->Set_HG_Edge_Wts_Fn(zoltanHGSizeEdgeWtsForMeshAdapter<Adapter>,
208  // (void *) &(*adp));
209  }
210 
211 public:
212 
218  AlgZoltan(const RCP<const Environment> &env__,
219  const RCP<const Comm<int> > &problemComm__,
220  const RCP<const IdentifierAdapter<user_t> > &adapter__):
221  env(env__), problemComm(problemComm__), adapter(adapter__)
222  {
223  setMPIComm(problemComm__);
224  zoltanInit();
225  zz = rcp(new Zoltan(mpicomm));
226  setCallbacksIDs();
227  }
228 
229  AlgZoltan(const RCP<const Environment> &env__,
230  const RCP<const Comm<int> > &problemComm__,
231  const RCP<const VectorAdapter<user_t> > &adapter__) :
232  env(env__), problemComm(problemComm__), adapter(adapter__)
233  {
234  setMPIComm(problemComm__);
235  zoltanInit();
236  zz = rcp(new Zoltan(mpicomm));
237  setCallbacksIDs();
238  setCallbacksGeom(&(*adapter));
239  }
240 
241  AlgZoltan(const RCP<const Environment> &env__,
242  const RCP<const Comm<int> > &problemComm__,
243  const RCP<const GraphAdapter<user_t,userCoord_t> > &adapter__) :
244  env(env__), problemComm(problemComm__), adapter(adapter__)
245  {
246  setMPIComm(problemComm__);
247  zoltanInit();
248  zz = rcp(new Zoltan(mpicomm));
249  setCallbacksIDs();
250  setCallbacksGraph(adapter);
251  if (adapter->coordinatesAvailable()) {
252  setCallbacksGeom(adapter->getCoordinateInput());
253  }
254  }
255 
256  AlgZoltan(const RCP<const Environment> &env__,
257  const RCP<const Comm<int> > &problemComm__,
258  const RCP<const MatrixAdapter<user_t,userCoord_t> > &adapter__) :
259  env(env__), problemComm(problemComm__), adapter(adapter__)
260  {
261  setMPIComm(problemComm__);
262  zoltanInit();
263  zz = rcp(new Zoltan(mpicomm));
264  setCallbacksIDs();
265  setCallbacksGraph(adapter);
266  setCallbacksHypergraph(adapter);
267  if (adapter->coordinatesAvailable()) {
268  setCallbacksGeom(adapter->getCoordinateInput());
269  }
270  }
271 
272  AlgZoltan(const RCP<const Environment> &env__,
273  const RCP<const Comm<int> > &problemComm__,
274  const RCP<const MeshAdapter<user_t> > &adapter__) :
275  env(env__), problemComm(problemComm__), adapter(adapter__)
276  {
277  setMPIComm(problemComm__);
278  zoltanInit();
279  zz = rcp(new Zoltan(mpicomm));
280  setCallbacksIDs();
281  setCallbacksGraph(adapter);
282  //TODO:: check parameter list to see if hypergraph is needed. We dont want to build the model
283  // if we don't have to and we shouldn't as it can take a decent amount of time if the
284  // primary entity is copied
285  setCallbacksHypergraph(adapter);
286  setCallbacksGeom(&(*adapter));
287  }
288 
289  void partition(const RCP<PartitioningSolution<Adapter> > &solution);
290  // void color(const RCP<ColoringSolution<Adapter> > &solution);
291 
292 };
293 
295 template <typename Adapter>
297  const RCP<PartitioningSolution<Adapter> > &solution
298 )
299 {
300  HELLO;
301  char paramstr[128];
302 
303  size_t numGlobalParts = solution->getTargetGlobalNumberOfParts();
304 
305  sprintf(paramstr, "%lu", numGlobalParts);
306  zz->Set_Param("NUM_GLOBAL_PARTS", paramstr);
307 
308  int wdim = adapter->getNumWeightsPerID();
309  sprintf(paramstr, "%d", wdim);
310  zz->Set_Param("OBJ_WEIGHT_DIM", paramstr);
311 
312  const Teuchos::ParameterList &pl = env->getParameters();
313 
314  double tolerance;
315  const Teuchos::ParameterEntry *pe = pl.getEntryPtr("imbalance_tolerance");
316  if (pe){
317  char str[30];
318  tolerance = pe->getValue<double>(&tolerance);
319  sprintf(str, "%f", tolerance);
320  zz->Set_Param("IMBALANCE_TOL", str);
321  }
322 
323 
324  // Look for zoltan_parameters sublist; pass all zoltan parameters to Zoltan
325  try {
326  const Teuchos::ParameterList &zpl = pl.sublist("zoltan_parameters");
327  for (ParameterList::ConstIterator iter = zpl.begin();
328  iter != zpl.end(); iter++) {
329  const std::string &zname = pl.name(iter);
330  // Convert the value to a string to pass to Zoltan
331  std::string zval = pl.entry(iter).getValue(&zval);
332  zz->Set_Param(zname.c_str(), zval.c_str());
333  }
334  }
335  catch (std::exception &e) {
336  // No zoltan_parameters sublist found; no Zoltan parameters to register
337  }
338 
339  // Get target part sizes
340  int pdim = (wdim > 1 ? wdim : 1);
341  for (int d = 0; d < pdim; d++) {
342  if (!solution->criteriaHasUniformPartSizes(d)) {
343  float *partsizes = new float[numGlobalParts];
344  int *partidx = new int[numGlobalParts];
345  int *wgtidx = new int[numGlobalParts];
346  for (size_t i=0; i<numGlobalParts; i++) partidx[i] = i;
347  for (size_t i=0; i<numGlobalParts; i++) wgtidx[i] = d;
348  for (size_t i=0; i<numGlobalParts; i++)
349  partsizes[i] = solution->getCriteriaPartSize(0, i);
350  zz->LB_Set_Part_Sizes(1, numGlobalParts, partidx, wgtidx, partsizes);
351  delete [] partsizes;
352  delete [] partidx;
353  delete [] wgtidx;
354  }
355  }
356 
357  // Make the call to LB_Partition
358  int changed = 0;
359  int nGidEnt = 1, nLidEnt = 1;
360  int nDummy = -1; // Dummy vars to satisfy arglist
361  ZOLTAN_ID_PTR dGids = NULL, dLids = NULL;
362  int *dProcs = NULL, *dParts = NULL;
363  int nObj = -1; // Output vars with new part info
364  ZOLTAN_ID_PTR oGids = NULL, oLids = NULL;
365  int *oProcs = NULL, *oParts = NULL;
366 
367  zz->Set_Param("RETURN_LISTS", "PARTS"); // required format for Zoltan2;
368  // results in last five arguments
369  int ierr = zz->LB_Partition(changed, nGidEnt, nLidEnt,
370  nDummy, dGids, dLids, dProcs, dParts,
371  nObj, oGids, oLids, oProcs, oParts);
372 
373  env->globalInputAssertion(__FILE__, __LINE__, "Zoltan LB_Partition",
374  (ierr==ZOLTAN_OK || ierr==ZOLTAN_WARN), BASIC_ASSERTION, problemComm);
375 
376  int numObjects=nObj;
377  //The number of objects may be larger than zoltan knows due to copies that were removed by the hypergraph model
378  if (model!=RCP<const Model<Adapter> >() &&
379  dynamic_cast<const HyperGraphModel<Adapter>* >(&(*model)) &&
380  !dynamic_cast<const HyperGraphModel<Adapter>* >(&(*model))->areVertexIDsUnique()) {
381  numObjects=model->getLocalNumObjects();
382  }
383 
384  // Load answer into the solution.
385  ArrayRCP<part_t> partList(new part_t[numObjects], 0, numObjects, true);
386  for (int i = 0; i < nObj; i++) partList[oLids[i]] = oParts[i];
387 
388  if (model!=RCP<const Model<Adapter> >() &&
389  dynamic_cast<const HyperGraphModel<Adapter>* >(&(*model)) &&
390  !dynamic_cast<const HyperGraphModel<Adapter>* >(&(*model))->areVertexIDsUnique()) {
391  //Setup the part ids for copied entities removed by ownership in hypergraph model.
392  const HyperGraphModel<Adapter>* mdl = static_cast<const HyperGraphModel<Adapter>* >(&(*model));
393 
394  typedef typename HyperGraphModel<Adapter>::map_t map_t;
395  Teuchos::RCP<const map_t> mapWithCopies;
396  Teuchos::RCP<const map_t> oneToOneMap;
397  mdl->getVertexMaps(mapWithCopies,oneToOneMap);
398 
399  typedef Tpetra::Vector<scalar_t, lno_t, gno_t> vector_t;
400  vector_t vecWithCopies(mapWithCopies);
401  vector_t oneToOneVec(oneToOneMap);
402 
403  // Set values in oneToOneVec: each entry == rank
404  assert(nObj == lno_t(oneToOneMap->getNodeNumElements()));
405  for (lno_t i = 0; i < nObj; i++)
406  oneToOneVec.replaceLocalValue(i, oParts[i]);
407 
408  // Now import oneToOneVec's values back to vecWithCopies
409  Teuchos::RCP<const Tpetra::Import<lno_t, gno_t> > importer =
410  Tpetra::createImport<lno_t, gno_t>(oneToOneMap, mapWithCopies);
411  vecWithCopies.doImport(oneToOneVec, *importer, Tpetra::REPLACE);
412 
413  // Should see copied vector values when print VEC WITH COPIES
414  lno_t nlocal = lno_t(mapWithCopies->getNodeNumElements());
415  for (lno_t i = 0; i < nlocal; i++)
416  partList[i] = vecWithCopies.getData()[i];
417 
418  }
419 
420  solution->setParts(partList);
421 
422  // Clean up
423  zz->LB_Free_Part(&oGids, &oLids, &oProcs, &oParts);
424 }
425 
426 } // namespace Zoltan2
427 
428 #endif
AlgZoltan(const RCP< const Environment > &env__, const RCP< const Comm< int > > &problemComm__, const RCP< const VectorAdapter< user_t > > &adapter__)
#define HELLO
IdentifierAdapter defines the interface for identifiers.
fast typical checks for valid arguments
MatrixAdapter defines the adapter interface for matrices.
Defines the Model interface.
AlgZoltan(const RCP< const Environment > &env__, const RCP< const Comm< int > > &problemComm__, const RCP< const MatrixAdapter< user_t, userCoord_t > > &adapter__)
GraphAdapter defines the interface for graph-based user data.
MeshAdapter defines the interface for mesh input.
std::bitset< NUM_MODEL_FLAGS > modelFlag_t
Defines the PartitioningSolution class.
void partition(const RCP< PartitioningSolution< Adapter > > &solution)
Partitioning method.
callback functions for the Zoltan package (templated on Adapter)
A PartitioningSolution is a solution to a partitioning problem.
VectorAdapter defines the interface for vector input.
AlgZoltan(const RCP< const Environment > &env__, const RCP< const Comm< int > > &problemComm__, const RCP< const MeshAdapter< user_t > > &adapter__)
Algorithm defines the base class for all algorithms.
Traits class to handle conversions between gno_t/lno_t and TPL data types (e.g., ParMETIS&#39;s idx_t...
Gathering definitions used in software development.
The base class for all model classes.
AlgZoltan(const RCP< const Environment > &env__, const RCP< const Comm< int > > &problemComm__, const RCP< const GraphAdapter< user_t, userCoord_t > > &adapter__)
A gathering of useful namespace methods.
HyperGraphModel defines the interface required for hyper graph models.
void getVertexMaps(Teuchos::RCP< const map_t > &copiesMap, Teuchos::RCP< const map_t > &onetooneMap) const
Sets pointers to the vertex map with copies and the vertex map without copies Note: the pointers will...
AlgZoltan(const RCP< const Environment > &env__, const RCP< const Comm< int > > &problemComm__, const RCP< const IdentifierAdapter< user_t > > &adapter__)