Zoltan2
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Macros Pages
Zoltan2_AlgZoltan.hpp
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 #ifndef _ZOLTAN2_ALGZOLTAN_HPP_
11 #define _ZOLTAN2_ALGZOLTAN_HPP_
12 
13 #include <Zoltan2_Standards.hpp>
14 #include <Zoltan2_Algorithm.hpp>
16 #include <Zoltan2_Util.hpp>
17 #include <Zoltan2_TPLTraits.hpp>
18 
19 #include <Zoltan2_Model.hpp>
20 
22 #include <zoltan_cpp.h>
23 #include <zoltan_partition_tree.h>
24 
28 //
29 // This first design templates Zoltan's callback functions on the
30 // input adapter. This approach has the advantage of simplicity and
31 // is most similar to current usage of Zoltan (where the callbacks define
32 // the model).
33 // A better approach might template them on a model,
34 // allowing Zoltan2 greater flexibility in creating models from the input.
35 // Alternatively, different callback implementations could be provided to
36 // represent different models to Zoltan.
38 
39 namespace Zoltan2 {
40 
41 template <typename Adapter>
42 class AlgZoltan : public Algorithm<Adapter>
43 {
44 
45 private:
46 
47  typedef typename Adapter::lno_t lno_t;
48  typedef typename Adapter::gno_t gno_t;
49  typedef typename Adapter::scalar_t scalar_t;
50  typedef typename Adapter::part_t part_t;
51  typedef typename Adapter::user_t user_t;
52  typedef typename Adapter::userCoord_t userCoord_t;
53 
54  const RCP<const Environment> env;
55  const RCP<const Comm<int> > problemComm;
56  const RCP<const typename Adapter::base_adapter_t> adapter;
57  RCP<const Model<Adapter> > model;
58  RCP<Zoltan> zz;
59 
60  MPI_Comm mpicomm;
61 
62  void setMPIComm(const RCP<const Comm<int> > &problemComm__) {
63 # ifdef HAVE_ZOLTAN2_MPI
64  mpicomm = Teuchos::getRawMpiComm(*problemComm__);
65 # else
66  mpicomm = MPI_COMM_WORLD; // taken from siMPI
67 # endif
68  }
69 
70  void zoltanInit() {
71  // call Zoltan_Initialize to make sure MPI_Init is called (in MPI or siMPI).
72  int argc = 0;
73  char **argv = NULL;
74  float ver;
75  Zoltan_Initialize(argc, argv, &ver);
76  }
77 
78  void setCallbacksIDs()
79  {
80  zz->Set_Num_Obj_Fn(zoltanNumObj<Adapter>, (void *) &(*adapter));
81  zz->Set_Obj_List_Fn(zoltanObjList<Adapter>, (void *) &(*adapter));
82 
83  const part_t *myparts;
84  adapter->getPartsView(myparts);
85  if (myparts != NULL)
86  zz->Set_Part_Multi_Fn(zoltanParts<Adapter>, (void *) &(*adapter));
87 
88  char tmp[4];
89  sprintf(tmp, "%d", TPL_Traits<ZOLTAN_ID_PTR, gno_t>::NUM_ID);
90  zz->Set_Param("NUM_GID_ENTRIES", tmp);
91  sprintf(tmp, "%d", TPL_Traits<ZOLTAN_ID_PTR, lno_t>::NUM_ID);
92  zz->Set_Param("NUM_LID_ENTRIES", tmp);
93  }
94 
95  template <typename AdapterWithCoords>
96  void setCallbacksGeom(const AdapterWithCoords *ia)
97  {
98  // Coordinates may be provided by the MeshAdapter or VectorAdapter.
99  // VectorAdapter may be provided directly by user or indirectly through
100  // GraphAdapter or MatrixAdapter. So separate template type is needed.
101  zz->Set_Num_Geom_Fn(zoltanNumGeom<AdapterWithCoords>, (void *) ia);
102  zz->Set_Geom_Multi_Fn(zoltanGeom<AdapterWithCoords>, (void *) ia);
103  }
104 
105  void setCallbacksGraph(
106  const RCP<const GraphAdapter<user_t,userCoord_t> > &/* adp */)
107  {
108  // std::cout << "NotReadyForGraphCallbacksYet" << std::endl;
109  // TODO
110  }
111 
112  void setCallbacksGraph(
113  const RCP<const MatrixAdapter<user_t,userCoord_t> > &adp)
114  {
115  // std::cout << "NotReadyForGraphCallbacksYet" << std::endl;
116  // TODO
117  }
118 
119  void setCallbacksGraph(
120  const RCP<const MeshAdapter<user_t> > &adp)
121  {
122  // std::cout << "NotReadyForGraphCallbacksYet" << std::endl;
123  // TODO
124  }
125 
126 #ifdef HAVE_ZOLTAN2_HYPERGRAPH
127  void setCallbacksHypergraph(
128  const RCP<const MatrixAdapter<user_t,userCoord_t> > &adp)
129  {
130  // TODO: If add parameter list to this function, can register
131  // TODO: different callbacks depending on the hypergraph model to use
132 
133  zz->Set_HG_Size_CS_Fn(zoltanHGSizeCS_withMatrixAdapter<Adapter>,
134  (void *) &(*adp));
135  zz->Set_HG_CS_Fn(zoltanHGCS_withMatrixAdapter<Adapter>,
136  (void *) &(*adp));
137 
138  // zz->Set_HG_Size_Edge_Wts_Fn(zoltanHGSizeEdgeWts_withMatrixAdapter<Adapter>,
139  // (void *) &(*adapter));
140  // zz->Set_HG_Edge_Wts_Fn(zoltanHGSizeEdgeWts_withMatrixAdapter<Adapter>,
141  // (void *) &(*adapter));
142  }
143 
144  void setCallbacksHypergraph(
145  const RCP<const GraphAdapter<user_t,userCoord_t> > &adp)
146  {
147  zz->Set_HG_Size_CS_Fn(zoltanHGSizeCS_withGraphAdapter<Adapter>,
148  (void *) &(*adp));
149  zz->Set_HG_CS_Fn(zoltanHGCS_withGraphAdapter<Adapter>,
150  (void *) &(*adp));
151 
152  if (adp->getNumWeightsPerEdge() != 0) {
153  if (adp->getNumWeightsPerEdge() > 1) {
154  std::cout << "Zoltan2 warning: getNumWeightsPerEdge() returned "
155  << adp->getNumWeightsPerEdge() << " but PHG supports only "
156  << " one weight per edge; only first weight will be used."
157  << std::endl;
158  }
159  zz->Set_HG_Size_Edge_Wts_Fn(zoltanHGSizeEdgeWts_withGraphAdapter<Adapter>,
160  (void *) &(*adapter));
161  zz->Set_HG_Edge_Wts_Fn(zoltanHGEdgeWts_withGraphAdapter<Adapter>,
162  (void *) &(*adapter));
163  }
164  }
165 
166  void setCallbacksHypergraph(const RCP<const MeshAdapter<user_t> > &adp)
167  {
168 
169  const Teuchos::ParameterList &pl = env->getParameters();
170 
171  const Teuchos::ParameterEntry *pe = pl.getEntryPtr("hypergraph_model_type");
172  std::string model_type("traditional");
173  if (pe){
174  model_type = pe->getValue<std::string>(&model_type);
175  }
176 
177  if (model_type=="ghosting" ||
178  !adp->areEntityIDsUnique(adp->getPrimaryEntityType())) {
179  Zoltan2::modelFlag_t flags;
181  problemComm, flags,
183  model = rcp(static_cast<const Model<Adapter>* >(mdl),true);
184 
185  zz->Set_Num_Obj_Fn(zoltanHGNumObj_withModel<Adapter>, (void *) &(*mdl));
186  zz->Set_Obj_List_Fn(zoltanHGObjList_withModel<Adapter>, (void *) &(*mdl));
187 
188  zz->Set_HG_Size_CS_Fn(zoltanHGSizeCS_withModel<Adapter>, (void *) &(*mdl));
189  zz->Set_HG_CS_Fn(zoltanHGCS_withModel<Adapter>, (void *) &(*mdl));
190  }
191  else {
192  //If entities are unique we dont need the extra cost of the model
193  zz->Set_HG_Size_CS_Fn(zoltanHGSizeCS_withMeshAdapter<Adapter>,
194  (void *) &(*adp));
195  zz->Set_HG_CS_Fn(zoltanHGCS_withMeshAdapter<Adapter>,
196  (void *) &(*adp));
197  }
198  // zz->Set_HG_Size_Edge_Wts_Fn(zoltanHGSizeEdgeWts_withMeshAdapter<Adapter>,
199  // (void *) &(*adp));
200  // zz->Set_HG_Edge_Wts_Fn(zoltanHGSizeEdgeWts_withMeshAdapter<Adapter>,
201  // (void *) &(*adp));
202  }
203  #endif
204 
206  virtual bool isPartitioningTreeBinary() const
207  {
208  return true;
209  }
210 
212  void rcb_recursive_partitionTree_calculations(
213  part_t arrayIndex,
214  part_t numParts,
215  std::vector<part_t> & splitRangeBeg,
216  std::vector<part_t> & splitRangeEnd) const
217  {
218  // Note the purpose of the recursive method is make sure the children of a
219  // node have updated their values for splitRangeBeg and splitRangeEnd
220  // Then we can set our own values simply based on the union
221  // first load the rcb data for the node
222  int parent = -1;
223  int left_leaf = -1;
224  int right_leaf = -1;
225  int err = Zoltan_RCB_Partition_Tree(zz->Get_C_Handle(),
226  arrayIndex - numParts + 1, // rcb starts as 1 but does not include terminals
227  &parent, &left_leaf, &right_leaf);
228  if(err != 0) {
229  throw std::logic_error( "Zoltan_RCB_Partition_Tree returned in error." );
230  }
231  // check that children both have their ranges set and if not, do those
232  // range first so we can build them to make our range
233  if(left_leaf > 0) { // neg is terminal and always already built
234  rcb_recursive_partitionTree_calculations(left_leaf+numParts-1, numParts,
235  splitRangeBeg, splitRangeEnd);
236  }
237  if(right_leaf > 0) { // neg is terminal and always already built
238  rcb_recursive_partitionTree_calculations(right_leaf+numParts-1, numParts,
239  splitRangeBeg, splitRangeEnd);
240  }
241  // now we can build our ranges from the children
242  // note this exploits the rcb conventions for right and left so we know
243  // that left_leaf will be our smaller indices
244  int leftIndex = (left_leaf > 0) ? (left_leaf-1+numParts) : (-left_leaf);
245  int rightIndex = (right_leaf > 0) ? (right_leaf-1+numParts) : (-right_leaf);
246  splitRangeBeg[arrayIndex] = splitRangeBeg[leftIndex];
247  splitRangeEnd[arrayIndex] = splitRangeEnd[rightIndex];
248  // for debugging sanity check verify left_leaf is a set of indices which
249  // goes continuously into the right_leaf
250  if(splitRangeBeg[rightIndex] != splitRangeEnd[leftIndex]) { // end is non-inclusive
251  throw std::logic_error( "RCB expected left_leaf indices and right leaf"
252  " indices to be continuous but it was not so." );
253  }
254  }
255 
257  void rcb_getPartitionTree(part_t numParts,
258  part_t & numTreeVerts,
259  std::vector<part_t> & permPartNums,
260  std::vector<part_t> & splitRangeBeg,
261  std::vector<part_t> & splitRangeEnd,
262  std::vector<part_t> & treeVertParents) const
263  {
264  // CALCULATE: numTreeVerts
265  // For rcb a tree node always takes 2 nodes and turns them into 1 node
266  // That means it takes numParts - 1 nodes to reduce a tree of numParts to
267  // a single root node - but we do 2 * numParts - 1 because we are currently
268  // treating all of the 'trivial' terminals as tree nodes - something we
269  // discussed we may change later
270  part_t numTreeNodes = 2 * numParts - 1;
271  numTreeVerts = numTreeNodes - 1; // by design convention root not included
272  // CALCULATE: permPartNums
273  permPartNums.resize(numParts);
274  for(part_t n = 0; n < numParts; ++n) {
275  permPartNums[n] = n; // for rcb we can assume this is trivial and in order
276  }
277  // CALCULATE: treeVertParents
278  treeVertParents.resize(numTreeNodes); // allocate space for numTreeNodes array
279  // scan all the non terminal nodes and set all children to have us as parent
280  // that will set all parents except for the root node which we will set to -1
281  // track the children of the root and final node for swapping later. Couple
282  // ways to do this - all seem a bit awkward but at least this is efficient.
283  part_t rootNode = 0; // track index of the root node for swapping
284  // a bit awkward but efficient - save the children of root and final node
285  // for swap at end to satisfy convention that root is highest index node
286  part_t saveRootNodeChildrenA = -1;
287  part_t saveRootNodeChildrenB = -1;
288  part_t saveFinalNodeChildrenA = -1;
289  part_t saveFinalNodeChildrenB = -1;
290  for(part_t n = numParts; n < numTreeNodes; ++n) { // scan and set all parents
291  int parent = -1;
292  int left_leaf = -1;
293  int right_leaf = -1;
294  int err = Zoltan_RCB_Partition_Tree(zz->Get_C_Handle(),
295  static_cast<int>(n - numParts) + 1, // rcb starts as 1 but does not include terminals
296  &parent, &left_leaf, &right_leaf);
297  if(err != 0) {
298  throw std::logic_error("Zoltan_RCB_Partition_Tree returned in error.");
299  }
300  part_t leftIndex = (left_leaf > 0) ? (left_leaf-1+numParts) : (-left_leaf);
301  part_t rightIndex = (right_leaf > 0) ? (right_leaf-1+numParts) : (-right_leaf);
302  treeVertParents[leftIndex] = n;
303  treeVertParents[rightIndex] = n;
304  // save root node for final swap
305  if(parent == 1 || parent == -1) { // is it the root?
306  rootNode = n; // remember I am the root
307  saveRootNodeChildrenA = leftIndex;
308  saveRootNodeChildrenB = rightIndex;
309  }
310  if(n == numTreeNodes-1) {
311  saveFinalNodeChildrenA = leftIndex;
312  saveFinalNodeChildrenB = rightIndex;
313  }
314  }
315  treeVertParents[rootNode] = -1; // convention parent is root -1
316  // splitRangeBeg and splitRangeEnd
317  splitRangeBeg.resize(numTreeNodes);
318  splitRangeEnd.resize(numTreeNodes);
319  // for terminal nodes this is trivial
320  for(part_t n = 0; n < numParts; ++n) {
321  splitRangeBeg[n] = n;
322  splitRangeEnd[n] = n + 1;
323  }
324  if(numParts > 1) { // not relevant for 1 part
325  // build the splitRangeBeg and splitRangeEnd recursively which forces the
326  // children of each node to be calculated before the parent so parent can
327  // just take the union of the two children
328  rcb_recursive_partitionTree_calculations(rootNode, numParts, splitRangeBeg,
329  splitRangeEnd);
330  // now as a final step handle the swap to root is the highest index node
331  // swap the parent of the two nodes
332  std::swap(treeVertParents[rootNode], treeVertParents[numTreeNodes-1]);
333  // get the children of the swapped nodes to have updated parents
334  treeVertParents[saveFinalNodeChildrenA] = rootNode;
335  treeVertParents[saveFinalNodeChildrenB] = rootNode;
336  // handle case where final node is child of the root
337  if(saveRootNodeChildrenA == numTreeNodes - 1) {
338  saveRootNodeChildrenA = rootNode;
339  }
340  if(saveRootNodeChildrenB == numTreeNodes - 1) {
341  saveRootNodeChildrenB = rootNode;
342  }
343  treeVertParents[saveRootNodeChildrenA] = numTreeNodes - 1;
344  treeVertParents[saveRootNodeChildrenB] = numTreeNodes - 1;
345  // update the beg and end indices - simply swap them
346  std::swap(splitRangeBeg[rootNode], splitRangeBeg[numTreeNodes-1]);
347  std::swap(splitRangeEnd[rootNode], splitRangeEnd[numTreeNodes-1]);
348  }
349  }
350 
352  void phg_getPartitionTree(part_t numParts,
353  part_t & numTreeVerts,
354  std::vector<part_t> & permPartNums,
355  std::vector<part_t> & splitRangeBeg,
356  std::vector<part_t> & splitRangeEnd,
357  std::vector<part_t> & treeVertParents) const
358  {
359  // First thing is to get the length of the tree from zoltan.
360  // The tree is a list of pairs (lo,hi) for each node.
361  // Here tree_array_size is the number of pairs.
362  // In phg indexing the first pairt (i=0) is empty garbage.
363  // The second pair (index 1) will be the root.
364  // Some nodes will be empty nodes, determined by hi = -1.
365  int tree_array_size = -1; // will be set by Zoltan_PHG_Partition_Tree_Size
366  int err = Zoltan_PHG_Partition_Tree_Size(
367  zz->Get_C_Handle(), &tree_array_size);
368  if(err != 0) {
369  throw std::logic_error("Zoltan_PHG_Partition_Tree_Size returned error.");
370  }
371  // Determine the number of valid nodes (PHG will have empty slots)
372  // We scan the list of pairs and count each node which does not have hi = -1
373  // During the loop we will also construct mapIndex which maps initial n
374  // to final n due to some conversions we apply to meet the design specs.
375  // The conversions implemented by mapIndex are:
376  // Move all terminals to the beginning (terminals have hi = lo)
377  // Resort the terminals in order (simply map to index lo works)
378  // Move non-terminals after the terminals (they start at index numParts)
379  // Map the first pair (root) to the be last to meet the design spec
380  part_t numTreeNodes = 0;
381  std::vector<part_t> mapIndex(tree_array_size); // maps n to final index
382  part_t trackNonTerminalIndex = numParts; // starts after terminals
383  for(part_t n = 0; n < static_cast<part_t>(tree_array_size); ++n) {
384  part_t phgIndex = n + 1; // phg indexing starts at 1
385  int lo_index = -1;
386  int hi_index = -1;
387  err = Zoltan_PHG_Partition_Tree(
388  zz->Get_C_Handle(), phgIndex, &lo_index, &hi_index);
389  if(hi_index != -1) { // hi -1 means it's an unused node
390  ++numTreeNodes; // increase the total count because this is a real node
391  if(n != 0) { // the root is mapped last but we don't know the length yet
392  mapIndex[n] = (lo_index == hi_index) ? // is it a terminal?
393  lo_index : // terminals map in sequence - lo_index is correct
394  (trackNonTerminalIndex++); // set then bump trackNonTerminalIndex +1
395  }
396  }
397  }
398  // now complete the map by setting root to the length-1 for the design specs
399  mapIndex[0] = numTreeNodes - 1;
400  // CALCULATE: numTreeVerts
401  numTreeVerts = numTreeNodes - 1; // this is the design - root not included
402  // CALCULATE: permPartNums
403  permPartNums.resize(numParts);
404  for(part_t n = 0; n < numParts; ++n) {
405  permPartNums[n] = n; // for phg we can assume this is trivial and in order
406  }
407  // CALCULATE: treeVertParents, splitRangeBeg, splitRangeEnd
408  // we will determine all of these in this second loop using mapIndex
409  // First set the arrays to have the proper length
410  treeVertParents.resize(numTreeNodes);
411  splitRangeBeg.resize(numTreeNodes);
412  splitRangeEnd.resize(numTreeNodes);
413  // Now loop a second time
414  for(part_t n = 0; n < tree_array_size; ++n) {
415  part_t phgIndex = n + 1; // phg indexing starts at 1
416  int lo_index = -1; // zoltan Zoltan_PHG_Partition_Tree will set this
417  int hi_index = -1; // zoltan Zoltan_PHG_Partition_Tree will set this
418  err = Zoltan_PHG_Partition_Tree( // access zoltan phg tree data
419  zz->Get_C_Handle(), phgIndex, &lo_index, &hi_index);
420  if(err != 0) {
421  throw std::logic_error("Zoltan_PHG_Partition_Tree returned in error.");
422  }
423  if(hi_index != -1) { // hi -1 means it's an unused node (a gap)
424  // get final index using mapIndex - so convert from phg to design plan
425  part_t finalIndex = mapIndex[n]; // get the index for the final output
426  // now determine the parent
427  // In the original phg indexing, the parent can be directly calculated
428  // from the pair index using the following rules:
429  // if phgIndex even, then parent is phgIndex/2
430  // here we determine even by ((phgIndex%2) == 0)
431  // if phgIndex odd, then parent is (phgIndex-1)/2
432  // but after getting parentPHGIndex we convert back to this array
433  // indexing by subtracting 1
434  part_t parentPHGIndex =
435  ((phgIndex%2) == 0) ? (phgIndex/2) : ((phgIndex-1)/2);
436  // map the parent phg index to the final parent index
437  // however, for the special case of the root (n=1), set the parent to -1
438  treeVertParents[finalIndex] = (n==0) ? -1 : mapIndex[parentPHGIndex-1];
439  // set begin (inclusive) and end (non-inclusive), so end is hi+1
440  splitRangeBeg[finalIndex] = static_cast<part_t>(lo_index);
441  splitRangeEnd[finalIndex] = static_cast<part_t>(hi_index+1);
442  }
443  }
444  }
445 
447  void getPartitionTree(part_t numParts,
448  part_t & numTreeVerts,
449  std::vector<part_t> & permPartNums,
450  std::vector<part_t> & splitRangeBeg,
451  std::vector<part_t> & splitRangeEnd,
452  std::vector<part_t> & treeVertParents) const
453  {
454  // first check that our parameters requested we keep the tree
455  const Teuchos::ParameterList &pl = env->getParameters();
456  bool keep_partition_tree = false;
457  const Teuchos::ParameterEntry * pe = pl.getEntryPtr("keep_partition_tree");
458  if(pe) {
459  keep_partition_tree = pe->getValue(&keep_partition_tree);
460  if(!keep_partition_tree) {
461  throw std::logic_error(
462  "Requested tree when param keep_partition_tree is off.");
463  }
464  }
465 
466  // now call the appropriate method based on LB_METHOD in the zoltan
467  // parameters list.
468  const Teuchos::ParameterList & zoltan_pl = pl.sublist("zoltan_parameters");
469  std::string lb_method;
470  pe = zoltan_pl.getEntryPtr("LB_METHOD");
471  if(pe) {
472  lb_method = pe->getValue(&lb_method);
473  }
474  if(lb_method == "phg") {
475  phg_getPartitionTree(numParts, numTreeVerts, permPartNums,
476  splitRangeBeg, splitRangeEnd, treeVertParents);
477  }
478  else if(lb_method == "rcb") {
479  rcb_getPartitionTree(numParts, numTreeVerts, permPartNums,
480  splitRangeBeg, splitRangeEnd, treeVertParents);
481  }
482  else {
483  throw std::logic_error("Did not recognize LB_METHOD: " + lb_method);
484  }
485  }
486 
487 public:
488 
494  AlgZoltan(const RCP<const Environment> &env__,
495  const RCP<const Comm<int> > &problemComm__,
496  const RCP<const IdentifierAdapter<user_t> > &adapter__):
497  env(env__), problemComm(problemComm__), adapter(adapter__)
498  {
499  setMPIComm(problemComm__);
500  zoltanInit();
501  zz = rcp(new Zoltan(mpicomm));
502  setCallbacksIDs();
503  }
504 
505  AlgZoltan(const RCP<const Environment> &env__,
506  const RCP<const Comm<int> > &problemComm__,
507  const RCP<const VectorAdapter<user_t> > &adapter__) :
508  env(env__), problemComm(problemComm__), adapter(adapter__)
509  {
510  setMPIComm(problemComm__);
511  zoltanInit();
512  zz = rcp(new Zoltan(mpicomm));
513  setCallbacksIDs();
514  setCallbacksGeom(&(*adapter));
515  }
516 
517  AlgZoltan(const RCP<const Environment> &env__,
518  const RCP<const Comm<int> > &problemComm__,
519  const RCP<const GraphAdapter<user_t,userCoord_t> > &adapter__) :
520  env(env__), problemComm(problemComm__), adapter(adapter__)
521  {
522  setMPIComm(problemComm__);
523  zoltanInit();
524  zz = rcp(new Zoltan(mpicomm));
525  setCallbacksIDs();
526  setCallbacksGraph(adapter);
527 #ifdef HAVE_ZOLTAN2_HYPERGRAPH
528  setCallbacksHypergraph(adapter);
529 #endif
530  if (adapter->coordinatesAvailable()) {
531  setCallbacksGeom(adapter->getCoordinateInput());
532  }
533  }
534 
535  AlgZoltan(const RCP<const Environment> &env__,
536  const RCP<const Comm<int> > &problemComm__,
537  const RCP<const MatrixAdapter<user_t,userCoord_t> > &adapter__) :
538  env(env__), problemComm(problemComm__), adapter(adapter__)
539  {
540  setMPIComm(problemComm__);
541  zoltanInit();
542  zz = rcp(new Zoltan(mpicomm));
543  setCallbacksIDs();
544  setCallbacksGraph(adapter);
545 #ifdef HAVE_ZOLTAN2_HYPERGRAPH
546  setCallbacksHypergraph(adapter);
547 #endif
548  if (adapter->coordinatesAvailable()) {
549  setCallbacksGeom(adapter->getCoordinateInput());
550  }
551  }
552 
553  AlgZoltan(const RCP<const Environment> &env__,
554  const RCP<const Comm<int> > &problemComm__,
555  const RCP<const MeshAdapter<user_t> > &adapter__) :
556  env(env__), problemComm(problemComm__), adapter(adapter__)
557  {
558  setMPIComm(problemComm__);
559  zoltanInit();
560  zz = rcp(new Zoltan(mpicomm));
561  setCallbacksIDs();
562  setCallbacksGraph(adapter);
563  //TODO:: check parameter list to see if hypergraph is needed. We dont want to build the model
564  // if we don't have to and we shouldn't as it can take a decent amount of time if the
565  // primary entity is copied
566 #ifdef HAVE_ZOLTAN2_HYPERGRAPH
567  setCallbacksHypergraph(adapter);
568 #endif
569  setCallbacksGeom(&(*adapter));
570  }
571 
572  void partition(const RCP<PartitioningSolution<Adapter> > &solution);
573  // void color(const RCP<ColoringSolution<Adapter> > &solution);
574 };
575 
577 template <typename Adapter>
579  const RCP<PartitioningSolution<Adapter> > &solution
580 )
581 {
582  HELLO;
583  char paramstr[128];
584 
585  size_t numGlobalParts = solution->getTargetGlobalNumberOfParts();
586 
587  sprintf(paramstr, "%lu", numGlobalParts);
588  zz->Set_Param("NUM_GLOBAL_PARTS", paramstr);
589 
590  int wdim = adapter->getNumWeightsPerID();
591  sprintf(paramstr, "%d", wdim);
592  zz->Set_Param("OBJ_WEIGHT_DIM", paramstr);
593 
594  const Teuchos::ParameterList &pl = env->getParameters();
595 
596  double tolerance;
597  const Teuchos::ParameterEntry *pe = pl.getEntryPtr("imbalance_tolerance");
598  if (pe){
599  char str[30];
600  tolerance = pe->getValue<double>(&tolerance);
601  sprintf(str, "%f", tolerance);
602  zz->Set_Param("IMBALANCE_TOL", str);
603  }
604 
605  pe = pl.getEntryPtr("partitioning_approach");
606  if (pe){
607  std::string approach;
608  approach = pe->getValue<std::string>(&approach);
609  if (approach == "partition")
610  zz->Set_Param("LB_APPROACH", "PARTITION");
611  else
612  zz->Set_Param("LB_APPROACH", "REPARTITION");
613  }
614 
615  pe = pl.getEntryPtr("partitioning_objective");
616  if (pe){
617  std::string strChoice = pe->getValue<std::string>(&strChoice);
618  if (strChoice == std::string("multicriteria_minimize_total_weight"))
619  zz->Set_Param("RCB_MULTICRITERIA_NORM", "1");
620  else if (strChoice == std::string("multicriteria_balance_total_maximum"))
621  zz->Set_Param("RCB_MULTICRITERIA_NORM", "2");
622  else if (strChoice == std::string("multicriteria_minimize_maximum_weight"))
623  zz->Set_Param("RCB_MULTICRITERIA_NORM", "3");
624  }
625 
626  // perhaps make this a bool stored in the AlgZoltan if we want to follow
627  // the pattern of multijagged mj_keep_part_boxes for example. However we can
628  // collect the error straight from Zoltan if we attempt to access the tree
629  // when we never stored it so that may not be necessary
630  bool keep_partition_tree = false;
631  pe = pl.getEntryPtr("keep_partition_tree");
632  if (pe) {
633  keep_partition_tree = pe->getValue(&keep_partition_tree);
634  if (keep_partition_tree) {
635  // need to resolve the organization of this
636  // when do we want to use the zoltan parameters directly versus
637  // using the zoltan2 parameters like this
638  zz->Set_Param("KEEP_CUTS", "1"); // rcb zoltan setting
639  zz->Set_Param("PHG_KEEP_TREE", "1"); // phg zoltan setting
640  }
641  }
642 
643  pe = pl.getEntryPtr("rectilinear");
644  if (pe) {
645  bool val = pe->getValue(&val);
646  if (val)
647  zz->Set_Param("RCB_RECTILINEAR_BLOCKS", "1");
648  }
649 
650  // Look for zoltan_parameters sublist; pass all zoltan parameters to Zoltan
651  try {
652  const Teuchos::ParameterList &zpl = pl.sublist("zoltan_parameters");
653  for (ParameterList::ConstIterator iter = zpl.begin();
654  iter != zpl.end(); iter++) {
655  const std::string &zname = pl.name(iter);
656  // Convert the value to a string to pass to Zoltan
657  std::string zval = pl.entry(iter).getValue(&zval);
658  zz->Set_Param(zname.c_str(), zval.c_str());
659  }
660  }
661  catch (std::exception &) {
662  // No zoltan_parameters sublist found; no Zoltan parameters to register
663  }
664 
665  // Get target part sizes
666  int pdim = (wdim > 1 ? wdim : 1);
667  for (int d = 0; d < pdim; d++) {
668  if (!solution->criteriaHasUniformPartSizes(d)) {
669  float *partsizes = new float[numGlobalParts];
670  int *partidx = new int[numGlobalParts];
671  int *wgtidx = new int[numGlobalParts];
672  for (size_t i=0; i<numGlobalParts; i++) partidx[i] = i;
673  for (size_t i=0; i<numGlobalParts; i++) wgtidx[i] = d;
674  for (size_t i=0; i<numGlobalParts; i++)
675  partsizes[i] = solution->getCriteriaPartSize(0, i);
676  zz->LB_Set_Part_Sizes(1, numGlobalParts, partidx, wgtidx, partsizes);
677  delete [] partsizes;
678  delete [] partidx;
679  delete [] wgtidx;
680  }
681  }
682 
683  // Make the call to LB_Partition
684  int changed = 0;
687 
688  int nDummy = -1; // Dummy vars to satisfy arglist
689  ZOLTAN_ID_PTR dGids = NULL, dLids = NULL;
690  int *dProcs = NULL, *dParts = NULL;
691  int nObj = -1; // Output vars with new part info
692  ZOLTAN_ID_PTR oGids = NULL, oLids = NULL;
693  int *oProcs = NULL, *oParts = NULL;
694 
695  zz->Set_Param("RETURN_LISTS", "PARTS"); // required format for Zoltan2;
696  // results in last five arguments
697 
698  int ierr = zz->LB_Partition(changed, nGidEnt, nLidEnt,
699  nDummy, dGids, dLids, dProcs, dParts,
700  nObj, oGids, oLids, oProcs, oParts);
701 
702  env->globalInputAssertion(__FILE__, __LINE__, "Zoltan LB_Partition",
703  (ierr==ZOLTAN_OK || ierr==ZOLTAN_WARN), BASIC_ASSERTION, problemComm);
704 
705  int numObjects=nObj;
706 #ifdef HAVE_ZOLTAN2_HYPERGRAPH
707  // The number of objects may be larger than zoltan knows due to copies that
708  // were removed by the hypergraph model
709  if (model!=RCP<const Model<Adapter> >() &&
710  dynamic_cast<const HyperGraphModel<Adapter>* >(&(*model)) &&
711  !(dynamic_cast<const HyperGraphModel<Adapter>* >(&(*model))->areVertexIDsUnique())) {
712  numObjects=model->getLocalNumObjects();
713  }
714 #endif
715 
716  // Load answer into the solution.
717  ArrayRCP<part_t> partList(new part_t[numObjects], 0, numObjects, true);
718  for (int i = 0; i < nObj; i++) {
719  lno_t tmp;
720  TPL_Traits<lno_t, ZOLTAN_ID_PTR>::ASSIGN(tmp, &(oLids[i*nLidEnt]));
721  partList[tmp] = oParts[i];
722  }
723 #ifdef HAVE_ZOLTAN2_HYPERGRAPH
724  if (model!=RCP<const Model<Adapter> >() &&
725  dynamic_cast<const HyperGraphModel<Adapter>* >(&(*model)) &&
726  !(dynamic_cast<const HyperGraphModel<Adapter>* >(&(*model))->areVertexIDsUnique())) {
727  // Setup the part ids for copied entities removed by ownership in
728  // hypergraph model.
729  const HyperGraphModel<Adapter>* mdl =
730  static_cast<const HyperGraphModel<Adapter>* >(&(*model));
731 
732  typedef typename HyperGraphModel<Adapter>::map_t map_t;
733  Teuchos::RCP<const map_t> mapWithCopies;
734  Teuchos::RCP<const map_t> oneToOneMap;
735  mdl->getVertexMaps(mapWithCopies,oneToOneMap);
736 
737  typedef Tpetra::Vector<scalar_t, lno_t, gno_t> vector_t;
738  vector_t vecWithCopies(mapWithCopies);
739  vector_t oneToOneVec(oneToOneMap);
740 
741  // Set values in oneToOneVec: each entry == rank
742  assert(nObj == lno_t(oneToOneMap->getLocalNumElements()));
743  for (lno_t i = 0; i < nObj; i++)
744  oneToOneVec.replaceLocalValue(i, oParts[i]);
745 
746  // Now import oneToOneVec's values back to vecWithCopies
747  Teuchos::RCP<const Tpetra::Import<lno_t, gno_t> > importer =
748  Tpetra::createImport<lno_t, gno_t>(oneToOneMap, mapWithCopies);
749  vecWithCopies.doImport(oneToOneVec, *importer, Tpetra::REPLACE);
750 
751  // Should see copied vector values when print VEC WITH COPIES
752  lno_t nlocal = lno_t(mapWithCopies->getLocalNumElements());
753  for (lno_t i = 0; i < nlocal; i++)
754  partList[i] = vecWithCopies.getData()[i];
755  }
756 #endif
757 
758  solution->setParts(partList);
759 
760  // Clean up
761  zz->LB_Free_Part(&oGids, &oLids, &oProcs, &oParts);
762 }
763 
764 } // namespace Zoltan2
765 
766 #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
map_t::global_ordinal_type gno_t
Definition: mapRemotes.cpp:27
Defines the PartitioningSolution class.
SparseMatrixAdapter_t::part_t part_t
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.
Tpetra::Map map_t
Definition: mapRemotes.cpp:25
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.
map_t::local_ordinal_type lno_t
Definition: mapRemotes.cpp:26
static void ASSIGN(first_t &a, second_t b)
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__)
Zoltan2::BasicUserTypes< zscalar_t, zlno_t, zgno_t > user_t
Definition: Metric.cpp:39