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