Actual source code: Partitioner.hh

  1: #ifndef included_ALE_Partitioner_hh
  2: #define included_ALE_Partitioner_hh

  4: #ifndef  included_ALE_Completion_hh
  5: #include <Completion.hh>
  6: #endif

  8: #ifdef PETSC_HAVE_ZOLTAN
  9: #include <zoltan.h>

 12:   // Inputs

 18:   int getNumVertices_Zoltan(void *, int *);

 20:   void getLocalElements_Zoltan(void *, int, int, ZOLTAN_ID_PTR, ZOLTAN_ID_PTR, int, float *, int *);

 22:   void getHgSizes_Zoltan(void *, int *, int *, int *, int *);

 24:   void getHg_Zoltan(void *, int, int, int, int, ZOLTAN_ID_PTR, int *, ZOLTAN_ID_PTR, int *);
 25: }

 27: #endif

 29: #ifdef PETSC_HAVE_CHACO
 30: /* Chaco does not have an include file */
 33:                        float *ewgts, float *x, float *y, float *z, char *outassignname,
 34:                        char *outfilename, short *assignment, int architecture, int ndims_tot,
 35:                        int mesh_dims[3], double *goal, int global_method, int local_method,
 36:                        int rqi_flag, int vmax, int ndims, double eigtol, long seed);

 39: }
 40: #endif
 41: #ifdef PETSC_HAVE_PARMETIS
 43:   #include <parmetis.h>
 45: }
 46: #endif
 47: #ifdef PETSC_HAVE_HMETIS
 50: }
 51: #endif

 53: namespace ALE {
 54:   namespace New {
 55:     template<typename Bundle_, typename Alloc_ = typename Bundle_::alloc_type>
 56:     class Partitioner {
 57:     public:
 58:       typedef Bundle_                          bundle_type;
 59:       typedef Alloc_                           alloc_type;
 60:       typedef typename bundle_type::sieve_type sieve_type;
 61:       typedef typename bundle_type::point_type point_type;
 62:     public:
 65:       // This creates a CSR representation of the adjacency matrix for cells
 66:       // - We allow an exception to contiguous numbering.
 67:       //   If the cell id > numElements, we assign a new number starting at
 68:       //     the top and going downward. I know these might not match up with
 69:       //     the iterator order, but we can fix it later.
 70:       static void buildDualCSR(const Obj<bundle_type>& bundle, const int dim, int **offsets, int **adjacency) {
 71:         ALE_LOG_EVENT_BEGIN;
 72:         typedef typename ALE::New::Completion<bundle_type, point_type, alloc_type> completion;
 73:         const Obj<sieve_type>&                           sieve        = bundle->getSieve();
 74:         const Obj<typename bundle_type::label_sequence>& elements     = bundle->heightStratum(0);
 75:         Obj<sieve_type>                                  overlapSieve = new sieve_type(bundle->comm(), bundle->debug());
 76:         std::map<point_type, point_type>                 newCells;
 77:         int  numElements = elements->size();
 78:         int  newCell     = numElements;
 79:         int *off         = new int[numElements+1];
 80:         int  offset      = 0;
 81:         int *adj;

 83:         completion::scatterSupports(sieve, overlapSieve, bundle->getSendOverlap(), bundle->getRecvOverlap(), bundle);
 84:         if (numElements == 0) {
 85:           *offsets   = NULL;
 86:           *adjacency = NULL;
 87:           ALE_LOG_EVENT_END;
 88:           return;
 89:         }
 90:         if (bundle->depth() == dim) {
 91:           int e = 1;

 93:           off[0] = 0;
 94:           for(typename bundle_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
 95:             const Obj<typename sieve_type::traits::coneSequence>& faces  = sieve->cone(*e_iter);
 96:             typename sieve_type::traits::coneSequence::iterator   fBegin = faces->begin();
 97:             typename sieve_type::traits::coneSequence::iterator   fEnd   = faces->end();

 99:             off[e] = off[e-1];
100:             for(typename sieve_type::traits::coneSequence::iterator f_iter = fBegin; f_iter != fEnd; ++f_iter) {
101:               if (sieve->support(*f_iter)->size() == 2) {
102:                 off[e]++;
103:               } else if ((sieve->support(*f_iter)->size() == 1) && (overlapSieve->support(*f_iter)->size() == 1)) {
104:                 off[e]++;
105:               }
106:             }
107:             e++;
108:           }
109:           adj = new int[off[numElements]];
110:           for(typename bundle_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
111:             const Obj<typename sieve_type::traits::coneSequence>& faces  = sieve->cone(*e_iter);
112:             typename sieve_type::traits::coneSequence::iterator   fBegin = faces->begin();
113:             typename sieve_type::traits::coneSequence::iterator   fEnd   = faces->end();

115:             for(typename sieve_type::traits::coneSequence::iterator f_iter = fBegin; f_iter != fEnd; ++f_iter) {
116:               const Obj<typename sieve_type::traits::supportSequence>& neighbors = sieve->support(*f_iter);
117:               typename sieve_type::traits::supportSequence::iterator   nBegin    = neighbors->begin();
118:               typename sieve_type::traits::supportSequence::iterator   nEnd      = neighbors->end();

120:               for(typename sieve_type::traits::supportSequence::iterator n_iter = nBegin; n_iter != nEnd; ++n_iter) {
121:                 if (*n_iter != *e_iter) adj[offset++] = *n_iter;
122:               }
123:               const Obj<typename sieve_type::traits::supportSequence>& oNeighbors = overlapSieve->support(*f_iter);
124:               typename sieve_type::traits::supportSequence::iterator   onBegin    = oNeighbors->begin();
125:               typename sieve_type::traits::supportSequence::iterator   onEnd      = oNeighbors->end();

127:               for(typename sieve_type::traits::supportSequence::iterator n_iter = onBegin; n_iter != onEnd; ++n_iter) {
128:                 adj[offset++] = *n_iter;
129:               }
130:             }
131:           }
132:         } else if (bundle->depth() == 1) {
133:           std::set<point_type> *neighborCells = new std::set<point_type>[numElements];
134:           int corners      = sieve->cone(*elements->begin())->size();
135:           int faceVertices = -1;

137:           if (corners == dim+1) {
138:             faceVertices = dim;
139:           } else if ((dim == 2) && (corners == 4)) {
140:             faceVertices = 2;
141:           } else if ((dim == 3) && (corners == 8)) {
142:             faceVertices = 4;
143:           } else {
144:             throw ALE::Exception("Could not determine number of face vertices");
145:           }
146:           for(typename bundle_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
147:             const Obj<typename sieve_type::traits::coneSequence>& vertices  = sieve->cone(*e_iter);
148:             typename sieve_type::traits::coneSequence::iterator vEnd = vertices->end();

150:             for(typename sieve_type::traits::coneSequence::iterator v_iter = vertices->begin(); v_iter != vEnd; ++v_iter) {
151:               const Obj<typename sieve_type::traits::supportSequence>& neighbors = sieve->support(*v_iter);
152:               typename sieve_type::traits::supportSequence::iterator nEnd = neighbors->end();

154:               for(typename sieve_type::traits::supportSequence::iterator n_iter = neighbors->begin(); n_iter != nEnd; ++n_iter) {
155:                 if (*e_iter == *n_iter) continue;
156:                 if ((int) sieve->nMeet(*e_iter, *n_iter, 1)->size() == faceVertices) {
157:                   if ((*e_iter < numElements) && (*n_iter < numElements)) {
158:                     neighborCells[*e_iter].insert(*n_iter);
159:                   } else {
160:                     point_type e = *e_iter, n = *n_iter;

162:                     if (*e_iter >= numElements) {
163:                       if (newCells.find(*e_iter) == newCells.end()) newCells[*e_iter] = --newCell;
164:                       e = newCells[*e_iter];
165:                     }
166:                     if (*n_iter >= numElements) {
167:                       if (newCells.find(*n_iter) == newCells.end()) newCells[*n_iter] = --newCell;
168:                       n = newCells[*n_iter];
169:                     }
170:                     neighborCells[e].insert(n);
171:                   }
172:                 }
173:               }
174:             }
175:           }
176:           off[0] = 0;
177:           for(int e = 1; e <= numElements; e++) {
178:             off[e] = neighborCells[e-1].size() + off[e-1];
179:           }
180:           adj = new int[off[numElements]];
181:           for(int e = 0; e < numElements; e++) {
182:             for(typename std::set<point_type>::iterator n_iter = neighborCells[e].begin(); n_iter != neighborCells[e].end(); ++n_iter) {
183:               adj[offset++] = *n_iter;
184:             }
185:           }
186:           delete [] neighborCells;
187:         } else {
188:           throw ALE::Exception("Dual creation not defined for partially interpolated meshes");
189:         }
190:         if (offset != off[numElements]) {
191:           ostringstream msg;
192:           msg << "ERROR: Total number of neighbors " << offset << " does not match the offset array " << off[numElements];
193:           throw ALE::Exception(msg.str().c_str());
194:         }
195:         //std::cout << "numElements: " << numElements << " newCell: " << newCell << std::endl;
196:         *offsets   = off;
197:         *adjacency = adj;
198:         ALE_LOG_EVENT_END;
199:       };
202:       // This creates a CSR representation of the adjacency hypergraph for faces
203:       static void buildFaceCSR(const Obj<bundle_type>& bundle, const int dim, const Obj<typename bundle_type::numbering_type>& fNumbering, int *numEdges, int **offsets, int **adjacency) {
204:         ALE_LOG_EVENT_BEGIN;
205:         const Obj<sieve_type>&                           sieve    = bundle->getSieve();
206:         const Obj<typename bundle_type::label_sequence>& elements = bundle->heightStratum(0);
207:         int  numElements = elements->size();
208:         int *off         = new int[numElements+1];
209:         int  e;

211:         if (bundle->depth() != dim) {
212:           throw ALE::Exception("Not yet implemented for non-interpolated meshes");
213:         }
214:         off[0] = 0;
215:         e      = 1;
216:         for(typename bundle_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
217:           off[e] = sieve->cone(*e_iter)->size() + off[e-1];
218:           e++;
219:         }
220:         int *adj    = new int[off[numElements]];
221:         int  offset = 0;
222:         for(typename bundle_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
223:           const Obj<typename sieve_type::traits::coneSequence>& faces = sieve->cone(*e_iter);
224:           typename sieve_type::traits::coneSequence::iterator   fEnd  = faces->end();

226:           for(typename sieve_type::traits::coneSequence::iterator f_iter = faces->begin(); f_iter != fEnd; ++f_iter) {
227:             adj[offset++] = fNumbering->getIndex(*f_iter);
228:           }
229:         }
230:         if (offset != off[numElements]) {
231:           ostringstream msg;
232:           msg << "ERROR: Total number of neighbors " << offset << " does not match the offset array " << off[numElements];
233:           throw ALE::Exception(msg.str().c_str());
234:         }
235:         *numEdges  = numElements;
236:         *offsets   = off;
237:         *adjacency = adj;
238:         ALE_LOG_EVENT_END;
239:       };
240:       template<typename PartitionType>
241:       static PartitionType *subordinatePartition(const Obj<bundle_type>& bundle, int levels, const Obj<bundle_type>& subBundle, const PartitionType assignment[]) {
242:         const Obj<typename bundle_type::numbering_type>& cNumbering = bundle->getFactory()->getLocalNumbering(bundle, bundle->depth());
243:         const Obj<typename bundle_type::label_sequence>& cells      = subBundle->heightStratum(0);
244:         const Obj<typename bundle_type::numbering_type>& sNumbering = bundle->getFactory()->getLocalNumbering(subBundle, subBundle->depth());
245:         const int        numCells      = cells->size();
246:         PartitionType   *subAssignment = new PartitionType[numCells];

248:         if (levels != 1) {
249:           throw ALE::Exception("Cannot calculate subordinate partition for any level separation other than 1");
250:         } else {
251:           const Obj<typename bundle_type::sieve_type>&   sieve    = bundle->getSieve();
252:           const Obj<typename bundle_type::sieve_type>&   subSieve = subBundle->getSieve();
253:           Obj<typename bundle_type::sieve_type::coneSet> tmpSet   = new typename bundle_type::sieve_type::coneSet();
254:           Obj<typename bundle_type::sieve_type::coneSet> tmpSet2  = new typename bundle_type::sieve_type::coneSet();

256:           for(typename bundle_type::label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
257:             const Obj<typename bundle_type::sieve_type::coneSequence>& cone = subSieve->cone(*c_iter);

259:             Obj<typename bundle_type::sieve_type::supportSet> cell = sieve->nJoin1(cone);
260:             if (cell->size() != 1) {
261:               std::cout << "Indeterminate subordinate partition for face " << *c_iter << std::endl;
262:               for(typename bundle_type::sieve_type::supportSet::iterator s_iter = cell->begin(); s_iter != cell->end(); ++s_iter) {
263:                 std::cout << "  cell " << *s_iter << std::endl;
264:               }
265:               // Could relax this to choosing the first one
266:               throw ALE::Exception("Indeterminate subordinate partition");
267:             }
268:             subAssignment[sNumbering->getIndex(*c_iter)] = assignment[cNumbering->getIndex(*cell->begin())];
269:             tmpSet->clear();
270:             tmpSet2->clear();
271:           }
272:         }
273:         return subAssignment;
274:       };
275:     };
276: #ifdef PETSC_HAVE_CHACO
277:     namespace Chaco {
278:       template<typename Bundle_>
279:       class Partitioner {
280:       public:
281:         typedef Bundle_                          bundle_type;
282:         typedef typename bundle_type::sieve_type sieve_type;
283:         typedef typename bundle_type::point_type point_type;
284:         typedef short int                        part_type;
285:       public:
288:         static part_type *partitionSieve(const Obj<bundle_type>& bundle, const int dim) {
289:           part_type *assignment = NULL; /* set number of each vtx (length n) */
290:           int       *start;             /* start of edge list for each vertex */
291:           int       *adjacency;         /* = adj -> j; edge list data  */

293:           ALE_LOG_EVENT_BEGIN;
294:           ALE::New::Partitioner<bundle_type>::buildDualCSR(bundle, dim, &start, &adjacency);
295:           if (bundle->commRank() == 0) {
296:             /* arguments for Chaco library */
297:             FREE_GRAPH = 0;                         /* Do not let Chaco free my memory */
298:             int nvtxs;                              /* number of vertices in full graph */
299:             int *vwgts = NULL;                      /* weights for all vertices */
300:             float *ewgts = NULL;                    /* weights for all edges */
301:             float *x = NULL, *y = NULL, *z = NULL;  /* coordinates for inertial method */
302:             char *outassignname = NULL;             /*  name of assignment output file */
303:             char *outfilename = NULL;               /* output file name */
304:             int architecture = dim;                 /* 0 => hypercube, d => d-dimensional mesh */
305:             int ndims_tot = 0;                      /* total number of cube dimensions to divide */
306:             int mesh_dims[3];                       /* dimensions of mesh of processors */
307:             double *goal = NULL;                    /* desired set sizes for each set */
308:             int global_method = 1;                  /* global partitioning algorithm */
309:             int local_method = 1;                   /* local partitioning algorithm */
310:             int rqi_flag = 0;                       /* should I use RQI/Symmlq eigensolver? */
311:             int vmax = 200;                         /* how many vertices to coarsen down to? */
312:             int ndims = 1;                          /* number of eigenvectors (2^d sets) */
313:             double eigtol = 0.001;                  /* tolerance on eigenvectors */
314:             long seed = 123636512;                  /* for random graph mutations */
315:             float *vCoords[3];

318:             PetscOptionsGetInt(PETSC_NULL, "-partitioner_chaco_global_method", &global_method, PETSC_NULL);CHKERROR(ierr, "Error in PetscOptionsGetInt");
319:             PetscOptionsGetInt(PETSC_NULL, "-partitioner_chaco_local_method",  &local_method,  PETSC_NULL);CHKERROR(ierr, "Error in PetscOptionsGetInt");
320:             if (global_method == 3) {
321:               // Inertial Partitioning
322:               PetscMalloc3(nvtxs,float,&x,nvtxs,float,&y,nvtxs,float,&z);CHKERROR(ierr, "Error in PetscMalloc");
323:               vCoords[0] = x; vCoords[1] = y; vCoords[2] = z;
324:               const Obj<typename bundle_type::label_sequence>&    cells       = bundle->heightStratum(0);
325:               const Obj<typename bundle_type::real_section_type>& coordinates = bundle->getRealSection("coordinates");
326:               const int corners = bundle->size(coordinates, *(cells->begin()))/dim;
327:               int       c       = 0;

329:               for(typename bundle_type::label_sequence::iterator c_iter = cells->begin(); c_iter !=cells->end(); ++c_iter, ++c) {
330:                 const double *coords = bundle->restrict(coordinates, *c_iter);

332:                 for(int d = 0; d < dim; ++d) {
333:                   vCoords[d][c] = 0.0;
334:                 }
335:                 for(int v = 0; v < corners; ++v) {
336:                   for(int d = 0; d < dim; ++d) {
337:                     vCoords[d][c] += coords[v*dim+d];
338:                   }
339:                 }
340:                 for(int d = 0; d < dim; ++d) {
341:                   vCoords[d][c] /= corners;
342:                 }
343:               }
344:             }

346:             nvtxs = bundle->heightStratum(0)->size();
347:             mesh_dims[0] = bundle->commSize(); mesh_dims[1] = 1; mesh_dims[2] = 1;
348:             for(int e = 0; e < start[nvtxs]; e++) {
349:               adjacency[e]++;
350:             }
351:             assignment = new part_type[nvtxs];
352:             PetscMemzero(assignment, nvtxs * sizeof(part_type));CHKERROR(ierr, "Error in PetscMemzero");

354:             /* redirect output to buffer: chaco -> msgLog */
355: #ifdef PETSC_HAVE_UNISTD_H
356:             char *msgLog;
357:             int fd_stdout, fd_pipe[2], count;

359:             fd_stdout = dup(1);
360:             pipe(fd_pipe);
361:             close(1);
362:             dup2(fd_pipe[1], 1);
363:             msgLog = new char[16284];
364: #endif

366:             interface(nvtxs, start, adjacency, vwgts, ewgts, x, y, z,
367:                              outassignname, outfilename, assignment, architecture, ndims_tot,
368:                              mesh_dims, goal, global_method, local_method, rqi_flag, vmax, ndims,
369:                              eigtol, seed);

371: #ifdef PETSC_HAVE_UNISTD_H
372:             int SIZE_LOG  = 10000;

374:             fflush(stdout);
375:             count = read(fd_pipe[0], msgLog, (SIZE_LOG - 1) * sizeof(char));
376:             if (count < 0) count = 0;
377:             msgLog[count] = 0;
378:             close(1);
379:             dup2(fd_stdout, 1);
380:             close(fd_stdout);
381:             close(fd_pipe[0]);
382:             close(fd_pipe[1]);
383:             if (bundle->debug()) {
384:               std::cout << msgLog << std::endl;
385:             }
386:             delete [] msgLog;
387: #endif
388:             if (global_method == 3) {
389:               // Inertial Partitioning
390:               PetscFree3(x, y, z);CHKERROR(ierr, "Error in PetscFree");
391:             }
392:           }
393:           if (adjacency) delete [] adjacency;
394:           if (start)     delete [] start;
395:           ALE_LOG_EVENT_END;
396:           return assignment;
397:         };
398:         static part_type *partitionSieveByFace(const Obj<bundle_type>& bundle, const int dim) {
399:           throw ALE::Exception("Chaco cannot partition a mesh by faces");
400:         };
401:       };
402:     };
403: #endif
404: #ifdef PETSC_HAVE_PARMETIS
405:     namespace ParMetis {
406:       template<typename Bundle_>
407:       class Partitioner {
408:       public:
409:         typedef Bundle_                          bundle_type;
410:         typedef typename bundle_type::sieve_type sieve_type;
411:         typedef typename bundle_type::point_type point_type;
412:         typedef int                              part_type;
413:       public:
416:         static part_type *partitionSieve(const Obj<bundle_type>& bundle, const int dim) {
417:           int    nvtxs      = 0;    // The number of vertices in full graph
418:           int   *vtxdist;           // Distribution of vertices across processes
419:           int   *xadj;              // Start of edge list for each vertex
420:           int   *adjncy;            // Edge lists for all vertices
421:           int   *vwgt       = NULL; // Vertex weights
422:           int   *adjwgt     = NULL; // Edge weights
423:           int    wgtflag    = 0;    // Indicates which weights are present
424:           int    numflag    = 0;    // Indicates initial offset (0 or 1)
425:           int    ncon       = 1;    // The number of weights per vertex
426:           int    nparts     = bundle->commSize(); // The number of partitions
427:           float *tpwgts;            // The fraction of vertex weights assigned to each partition
428:           float *ubvec;             // The balance intolerance for vertex weights
429:           int    options[5];        // Options
430:           // Outputs
431:           int    edgeCut;           // The number of edges cut by the partition
432:           int   *assignment = NULL; // The vertex partition

434:           options[0] = 0; // Use all defaults
435:           vtxdist    = new int[nparts+1];
436:           vtxdist[0] = 0;
437:           tpwgts     = new float[ncon*nparts];
438:           for(int p = 0; p < nparts; ++p) {
439:             tpwgts[p] = 1.0/nparts;
440:           }
441:           ubvec      = new float[ncon];
442:           ubvec[0]   = 1.05;
443:           nvtxs      = bundle->heightStratum(0)->size();
444:           assignment = new part_type[nvtxs];
445:           MPI_Allgather(&nvtxs, 1, MPI_INT, &vtxdist[1], 1, MPI_INT, bundle->comm());
446:           for(int p = 2; p <= nparts; ++p) {
447:             vtxdist[p] += vtxdist[p-1];
448:           }
449:           if (bundle->commSize() == 1) {
450:             PetscMemzero(assignment, nvtxs * sizeof(part_type));
451:           } else {
452:             ALE::New::Partitioner<bundle_type>::buildDualCSR(bundle, dim, &xadj, &adjncy);

454:             if (bundle->debug() && nvtxs) {
455:               for(int p = 0; p <= nvtxs; ++p) {
456:                 std::cout << "["<<bundle->commRank()<<"]xadj["<<p<<"] = " << xadj[p] << std::endl;
457:               }
458:               for(int i = 0; i < xadj[nvtxs]; ++i) {
459:                 std::cout << "["<<bundle->commRank()<<"]adjncy["<<i<<"] = " << adjncy[i] << std::endl;
460:               }
461:             }
462:             if (vtxdist[1] == vtxdist[nparts]) {
463:               if (bundle->commRank() == 0) {
464:                 METIS_PartGraphKway(&nvtxs, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &nparts, options, &edgeCut, assignment);
465:                 if (bundle->debug()) {std::cout << "Metis: edgecut is " << edgeCut << std::endl;}
466:               }
467:             } else {
468:               MPI_Comm comm = bundle->comm();

470:               ParMETIS_V3_PartKway(vtxdist, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgeCut, assignment, &comm);
471:               if (bundle->debug()) {std::cout << "ParMetis: edgecut is " << edgeCut << std::endl;}
472:             }
473:             if (xadj   != NULL) delete [] xadj;
474:             if (adjncy != NULL) delete [] adjncy;
475:           }
476:           delete [] vtxdist;
477:           delete [] tpwgts;
478:           delete [] ubvec;
479:           return assignment;
480:         };
483:         static part_type *partitionSieveByFace(const Obj<bundle_type>& bundle, const int dim) {
484: #ifdef PETSC_HAVE_HMETIS
485:           int   *assignment = NULL; // The vertex partition
486:           int    nvtxs;      // The number of vertices
487:           int    nhedges;    // The number of hyperedges
488:           int   *vwgts;      // The vertex weights
489:           int   *eptr;       // The offsets of each hyperedge
490:           int   *eind;       // The vertices in each hyperedge, indexed by eptr
491:           int   *hewgts;     // The hyperedge weights
492:           int    nparts;     // The number of partitions
493:           int    ubfactor;   // The allowed load imbalance (1-50)
494:           int    options[9]; // Options
495:           // Outputs
496:           int    edgeCut;    // The number of edges cut by the partition
497:           const Obj<ALE::Mesh::numbering_type>& fNumbering = bundle->getFactory()->getNumbering(bundle, bundle->depth()-1);

499:           if (topology->commRank() == 0) {
500:             nvtxs      = bundle->heightStratum(1)->size();
501:             vwgts      = NULL;
502:             hewgts     = NULL;
503:             nparts     = bundle->commSize();
504:             ubfactor   = 5;
505:             options[0] = 1;  // Use all defaults
506:             options[1] = 10; // Number of bisections tested
507:             options[2] = 1;  // Vertex grouping scheme
508:             options[3] = 1;  // Objective function
509:             options[4] = 1;  // V-cycle refinement
510:             options[5] = 0;
511:             options[6] = 0;
512:             options[7] = 1; // Random seed
513:             options[8] = 24; // Debugging level
514:             assignment = new part_type[nvtxs];

516:             if (bundle->commSize() == 1) {
517:               PetscMemzero(assignment, nvtxs * sizeof(part_type));
518:             } else {
519:               ALE::New::Partitioner<bundle_type>::buildFaceCSR(bundle, dim, fNumbering, &nhedges, &eptr, &eind);
520:               HMETIS_PartKway(nvtxs, nhedges, vwgts, eptr, eind, hewgts, nparts, ubfactor, options, assignment, &edgeCut);

522:               delete [] eptr;
523:               delete [] eind;
524:             }
525:             if (bundle->debug()) {for (int i = 0; i<nvtxs; i++) printf("[%d] %d\n", PetscGlobalRank, assignment[i]);}
526:           } else {
527:             assignment = NULL;
528:           }
529:           return assignment;
530: #else
531:           throw ALE::Exception("hmetis partitioner is not available.");
532: #endif
533:         };
534:       };
535:     };
536: #endif
537: #ifdef PETSC_HAVE_ZOLTAN
538:     namespace Zoltan {
539:       template<typename Bundle_>
540:       class Partitioner {
541:       public:
542:         typedef Bundle_                          bundle_type;
543:         typedef typename bundle_type::sieve_type sieve_type;
544:         typedef typename bundle_type::point_type point_type;
545:         typedef int                              part_type;
546:       public:
547:         static part_type *partitionSieve(const Obj<bundle_type>& bundle, const int dim) {
548:           throw ALE::Exception("Zoltan partition by cells not implemented");
549:         };
552:         static part_type *partitionSieveByFace(const Obj<bundle_type>& bundle, const int dim) {
553:           // Outputs
554:           float         version;           // The library version
555:           int           changed;           // Did the partition change?
556:           int           numGidEntries;     // Number of array entries for a single global ID (1)
557:           int           numLidEntries;     // Number of array entries for a single local ID (1)
558:           int           numImport;         // The number of imported points
559:           ZOLTAN_ID_PTR import_global_ids; // The imported points
560:           ZOLTAN_ID_PTR import_local_ids;  // The imported points
561:           int          *import_procs;      // The proc each point was imported from
562:           int          *import_to_part;    // The partition of each imported point
563:           int           numExport;         // The number of exported points
564:           ZOLTAN_ID_PTR export_global_ids; // The exported points
565:           ZOLTAN_ID_PTR export_local_ids;  // The exported points
566:           int          *export_procs;      // The proc each point was exported to
567:           int          *export_to_part;    // The partition assignment of all local points
568:           int          *assignment;        // The partition assignment of all local points
569:           const Obj<typename bundle_type::numbering_type>& fNumbering = bundle->getFactory()->getNumbering(bundle, bundle->depth()-1);

571:           if (bundle->commSize() == 1) {
572:             PetscMemzero(assignment, bundle->heightStratum(1)->size() * sizeof(part_type));
573:           } else {
574:             if (bundle->commRank() == 0) {
575:               nvtxs_Zoltan = bundle->heightStratum(1)->size();
576:               ALE::New::Partitioner<bundle_type>::buildFaceCSR(bundle, dim, fNumbering, &nhedges_Zoltan, &eptr_Zoltan, &eind_Zoltan);
577:               assignment = new int[nvtxs_Zoltan];
578:             } else {
579:               nvtxs_Zoltan   = bundle->heightStratum(1)->size();
580:               nhedges_Zoltan = 0;
581:               eptr_Zoltan    = new int[1];
582:               eind_Zoltan    = new int[1];
583:               eptr_Zoltan[0] = 0;
584:               assignment     = NULL;
585:             }

587:             int Zoltan_Initialize(0, NULL, &version);
588:             struct Zoltan_Struct *zz = Zoltan_Create(bundle->comm());
589:             // General parameters
590:             Zoltan_Set_Param(zz, "DEBUG_LEVEL", "2");
591:             Zoltan_Set_Param(zz, "LB_METHOD", "PHG");
592:             Zoltan_Set_Param(zz, "RETURN_LISTS", "PARTITION");
593:             // PHG parameters
594:             Zoltan_Set_Param(zz, "PHG_OUTPUT_LEVEL", "2");
595:             Zoltan_Set_Param(zz, "PHG_EDGE_SIZE_THRESHOLD", "1.0"); // Do not throw out dense edges
596:             // Call backs
597:             Zoltan_Set_Num_Obj_Fn(zz, getNumVertices_Zoltan, NULL);
598:             Zoltan_Set_Obj_List_Fn(zz, getLocalElements_Zoltan, NULL);
599:             Zoltan_Set_HG_Size_CS_Fn(zz, getHgSizes_Zoltan, NULL);
600:             Zoltan_Set_HG_CS_Fn(zz, getHg_Zoltan, NULL);
601:             // Debugging
602:             //Zoltan_Generate_Files(zz, "zoltan.debug", 1, 0, 0, 1); // if using hypergraph callbacks

604:             Zoltan_LB_Partition(zz, &changed, &numGidEntries, &numLidEntries,
605:                                        &numImport, &import_global_ids, &import_local_ids, &import_procs, &import_to_part,
606:                                        &numExport, &export_global_ids, &export_local_ids, &export_procs, &export_to_part);
607:             for(int v = 0; v < nvtxs_Zoltan; ++v) {
608:               assignment[v] = export_to_part[v];
609:             }
610:             Zoltan_LB_Free_Part(&import_global_ids, &import_local_ids, &import_procs, &import_to_part);
611:             Zoltan_LB_Free_Part(&export_global_ids, &export_local_ids, &export_procs, &export_to_part);
612:             Zoltan_Destroy(&zz);

614:             delete [] eptr_Zoltan;
615:             delete [] eind_Zoltan;
616:           }
617:           if (assignment) {for (int i=0; i<nvtxs_Zoltan; i++) printf("[%d] %d\n",PetscGlobalRank,assignment[i]);}
618:           return assignment;
619:         };
620:       };
621:     };
622: #endif
623:   }
624: }

626: namespace ALECompat {
627:   namespace New {
628:     template<typename Topology_>
629:     class Partitioner {
630:     public:
631:       typedef Topology_                          topology_type;
632:       typedef typename topology_type::sieve_type sieve_type;
633:       typedef typename topology_type::patch_type patch_type;
634:       typedef typename topology_type::point_type point_type;
635:     public:
638:       // This creates a CSR representation of the adjacency matrix for cells
639:       static void buildDualCSR(const Obj<topology_type>& topology, const int dim, const patch_type& patch, int **offsets, int **adjacency) {
640:         ALE_LOG_EVENT_BEGIN;
641:         typedef typename ALECompat::New::Completion<topology_type, typename Mesh::sieve_type::point_type> completion;
642:         const Obj<sieve_type>&                             sieve        = topology->getPatch(patch);
643:         const Obj<typename topology_type::label_sequence>& elements     = topology->heightStratum(patch, 0);
644:         Obj<sieve_type>                                    overlapSieve = new sieve_type(topology->comm(), topology->debug());
645:         int  numElements = elements->size();
646:         int *off         = new int[numElements+1];
647:         int  offset      = 0;
648:         int *adj;

650:         completion::scatterSupports(sieve, overlapSieve, topology->getSendOverlap(), topology->getRecvOverlap(), topology);
651:         if (numElements == 0) {
652:           *offsets   = NULL;
653:           *adjacency = NULL;
654:           ALE_LOG_EVENT_END;
655:           return;
656:         }
657:         if (topology->depth(patch) == dim) {
658:           int e = 1;

660:           off[0] = 0;
661:           for(typename topology_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
662:             const Obj<typename sieve_type::traits::coneSequence>& faces  = sieve->cone(*e_iter);
663:             typename sieve_type::traits::coneSequence::iterator   fBegin = faces->begin();
664:             typename sieve_type::traits::coneSequence::iterator   fEnd   = faces->end();

666:             off[e] = off[e-1];
667:             for(typename sieve_type::traits::coneSequence::iterator f_iter = fBegin; f_iter != fEnd; ++f_iter) {
668:               if (sieve->support(*f_iter)->size() == 2) {
669:                 off[e]++;
670:               } else if ((sieve->support(*f_iter)->size() == 1) && (overlapSieve->support(*f_iter)->size() == 1)) {
671:                 off[e]++;
672:               }
673:             }
674:             e++;
675:           }
676:           adj = new int[off[numElements]];
677:           for(typename topology_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
678:             const Obj<typename sieve_type::traits::coneSequence>& faces  = sieve->cone(*e_iter);
679:             typename sieve_type::traits::coneSequence::iterator   fBegin = faces->begin();
680:             typename sieve_type::traits::coneSequence::iterator   fEnd   = faces->end();

682:             for(typename sieve_type::traits::coneSequence::iterator f_iter = fBegin; f_iter != fEnd; ++f_iter) {
683:               const Obj<typename sieve_type::traits::supportSequence>& neighbors = sieve->support(*f_iter);
684:               typename sieve_type::traits::supportSequence::iterator   nBegin    = neighbors->begin();
685:               typename sieve_type::traits::supportSequence::iterator   nEnd      = neighbors->end();

687:               for(typename sieve_type::traits::supportSequence::iterator n_iter = nBegin; n_iter != nEnd; ++n_iter) {
688:                 if (*n_iter != *e_iter) adj[offset++] = *n_iter;
689:               }
690:               const Obj<typename sieve_type::traits::supportSequence>& oNeighbors = overlapSieve->support(*f_iter);
691:               typename sieve_type::traits::supportSequence::iterator   onBegin    = oNeighbors->begin();
692:               typename sieve_type::traits::supportSequence::iterator   onEnd      = oNeighbors->end();

694:               for(typename sieve_type::traits::supportSequence::iterator n_iter = onBegin; n_iter != onEnd; ++n_iter) {
695:                 adj[offset++] = *n_iter;
696:               }
697:             }
698:           }
699:         } else if (topology->depth(patch) == 1) {
700:           std::set<point_type> *neighborCells = new std::set<point_type>[numElements];
701:           int corners      = sieve->cone(*elements->begin())->size();
702:           int faceVertices = -1;

704:           if (corners == dim+1) {
705:             faceVertices = dim;
706:           } else if ((dim == 2) && (corners == 4)) {
707:             faceVertices = 2;
708:           } else if ((dim == 3) && (corners == 8)) {
709:             faceVertices = 4;
710:           } else {
711:             throw ALE::Exception("Could not determine number of face vertices");
712:           }
713:           for(typename topology_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
714:             const Obj<typename sieve_type::traits::coneSequence>& vertices  = sieve->cone(*e_iter);
715:             typename sieve_type::traits::coneSequence::iterator vEnd = vertices->end();

717:             for(typename sieve_type::traits::coneSequence::iterator v_iter = vertices->begin(); v_iter != vEnd; ++v_iter) {
718:               const Obj<typename sieve_type::traits::supportSequence>& neighbors = sieve->support(*v_iter);
719:               typename sieve_type::traits::supportSequence::iterator nEnd = neighbors->end();

721:               for(typename sieve_type::traits::supportSequence::iterator n_iter = neighbors->begin(); n_iter != nEnd; ++n_iter) {
722:                 if (*e_iter == *n_iter) continue;
723:                 if ((int) sieve->meet(*e_iter, *n_iter)->size() == faceVertices) {
724:                   neighborCells[*e_iter].insert(*n_iter);
725:                 }
726:               }
727:             }
728:           }
729:           off[0] = 0;
730:           for(int e = 1; e <= numElements; e++) {
731:             off[e] = neighborCells[e-1].size() + off[e-1];
732:           }
733:           adj = new int[off[numElements]];
734:           for(int e = 0; e < numElements; e++) {
735:             for(typename std::set<point_type>::iterator n_iter = neighborCells[e].begin(); n_iter != neighborCells[e].end(); ++n_iter) {
736:               adj[offset++] = *n_iter;
737:             }
738:           }
739:           delete [] neighborCells;
740:         } else {
741:           throw ALE::Exception("Dual creation not defined for partially interpolated meshes");
742:         }
743:         if (offset != off[numElements]) {
744:           ostringstream msg;
745:           msg << "ERROR: Total number of neighbors " << offset << " does not match the offset array " << off[numElements];
746:           throw ALE::Exception(msg.str().c_str());
747:         }
748:         *offsets   = off;
749:         *adjacency = adj;
750:         ALE_LOG_EVENT_END;
751:       };
752:       template<typename PartitionType>
753:       static PartitionType *subordinatePartition(const Obj<topology_type>& topology, int levels, const Obj<topology_type>& subTopology, const PartitionType assignment[]) {
754:         typedef ALECompat::New::NumberingFactory<topology_type> NumberingFactory;
755:         const patch_type patch = 0;
756:         const Obj<typename NumberingFactory::numbering_type>& cNumbering = NumberingFactory::singleton(topology->debug())->getLocalNumbering(topology, patch, topology->depth(patch));
757:         const Obj<typename topology_type::label_sequence>&    cells      = subTopology->heightStratum(patch, 0);
758:         const Obj<typename NumberingFactory::numbering_type>& sNumbering = NumberingFactory::singleton(subTopology->debug())->getLocalNumbering(subTopology, patch, subTopology->depth(patch));
759:         const int        numCells      = cells->size();
760:         PartitionType   *subAssignment = new PartitionType[numCells];

762:         if (levels != 1) {
763:           throw ALE::Exception("Cannot calculate subordinate partition for any level separation other than 1");
764:         } else {
765:           const Obj<typename topology_type::sieve_type>&   sieve    = topology->getPatch(patch);
766:           const Obj<typename topology_type::sieve_type>&   subSieve = subTopology->getPatch(patch);
767:           Obj<typename topology_type::sieve_type::coneSet> tmpSet   = new typename topology_type::sieve_type::coneSet();
768:           Obj<typename topology_type::sieve_type::coneSet> tmpSet2  = new typename topology_type::sieve_type::coneSet();

770:           for(typename topology_type::label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
771:             const Obj<typename topology_type::sieve_type::coneSequence>& cone = subSieve->cone(*c_iter);

773:             Obj<typename topology_type::sieve_type::supportSet> cell = sieve->nJoin1(cone);
774:             if (cell->size() != 1) {
775:               std::cout << "Indeterminate subordinate partition for face " << *c_iter << std::endl;
776:               for(typename topology_type::sieve_type::supportSet::iterator s_iter = cell->begin(); s_iter != cell->end(); ++s_iter) {
777:                 std::cout << "  cell " << *s_iter << std::endl;
778:               }
779:               // Could relax this to choosing the first one
780:               throw ALE::Exception("Indeterminate subordinate partition");
781:             }
782:             subAssignment[sNumbering->getIndex(*c_iter)] = assignment[cNumbering->getIndex(*cell->begin())];
783:             tmpSet->clear();
784:             tmpSet2->clear();
785:           }
786:         }
787:         return subAssignment;
788:       };
789:     };
790: #ifdef PETSC_HAVE_CHACO
791:     namespace Chaco {
792:       template<typename Topology_>
793:       class Partitioner {
794:       public:
795:         typedef Topology_                          topology_type;
796:         typedef typename topology_type::sieve_type sieve_type;
797:         typedef typename topology_type::patch_type patch_type;
798:         typedef typename topology_type::point_type point_type;
799:         typedef short int                          part_type;
800:       public:
803:         static part_type *partitionSieve(const Obj<topology_type>& topology, const int dim) {
804:           part_type *assignment = NULL; /* set number of each vtx (length n) */
805:           int       *start;             /* start of edge list for each vertex */
806:           int       *adjacency;         /* = adj -> j; edge list data  */
807:           typename topology_type::patch_type patch = 0;

809:           ALE_LOG_EVENT_BEGIN;
810:           ALECompat::New::Partitioner<topology_type>::buildDualCSR(topology, dim, patch, &start, &adjacency);
811:           if (topology->commRank() == 0) {
812:             /* arguments for Chaco library */
813:             FREE_GRAPH = 0;                         /* Do not let Chaco free my memory */
814:             int nvtxs;                              /* number of vertices in full graph */
815:             int *vwgts = NULL;                      /* weights for all vertices */
816:             float *ewgts = NULL;                    /* weights for all edges */
817:             float *x = NULL, *y = NULL, *z = NULL;  /* coordinates for inertial method */
818:             char *outassignname = NULL;             /*  name of assignment output file */
819:             char *outfilename = NULL;               /* output file name */
820:             int architecture = 1;                   /* 0 => hypercube, d => d-dimensional mesh */
821:             int ndims_tot = 0;                      /* total number of cube dimensions to divide */
822:             int mesh_dims[3];                       /* dimensions of mesh of processors */
823:             double *goal = NULL;                    /* desired set sizes for each set */
824:             int global_method = 1;                  /* global partitioning algorithm */
825:             int local_method = 1;                   /* local partitioning algorithm */
826:             int rqi_flag = 0;                       /* should I use RQI/Symmlq eigensolver? */
827:             int vmax = 200;                         /* how many vertices to coarsen down to? */
828:             int ndims = 1;                          /* number of eigenvectors (2^d sets) */
829:             double eigtol = 0.001;                  /* tolerance on eigenvectors */
830:             long seed = 123636512;                  /* for random graph mutations */

833:             nvtxs = topology->heightStratum(patch, 0)->size();
834:             mesh_dims[0] = topology->commSize(); mesh_dims[1] = 1; mesh_dims[2] = 1;
835:             for(int e = 0; e < start[nvtxs]; e++) {
836:               adjacency[e]++;
837:             }
838:             assignment = new part_type[nvtxs];
839:             PetscMemzero(assignment, nvtxs * sizeof(part_type));

841:             /* redirect output to buffer: chaco -> msgLog */
842: #ifdef PETSC_HAVE_UNISTD_H
843:             char *msgLog;
844:             int fd_stdout, fd_pipe[2], count;

846:             fd_stdout = dup(1);
847:             pipe(fd_pipe);
848:             close(1);
849:             dup2(fd_pipe[1], 1);
850:             msgLog = new char[16284];
851: #endif

853:             interface(nvtxs, start, adjacency, vwgts, ewgts, x, y, z,
854:                              outassignname, outfilename, assignment, architecture, ndims_tot,
855:                              mesh_dims, goal, global_method, local_method, rqi_flag, vmax, ndims,
856:                              eigtol, seed);

858: #ifdef PETSC_HAVE_UNISTD_H
859:             int SIZE_LOG  = 10000;

861:             fflush(stdout);
862:             count = read(fd_pipe[0], msgLog, (SIZE_LOG - 1) * sizeof(char));
863:             if (count < 0) count = 0;
864:             msgLog[count] = 0;
865:             close(1);
866:             dup2(fd_stdout, 1);
867:             close(fd_stdout);
868:             close(fd_pipe[0]);
869:             close(fd_pipe[1]);
870:             if (topology->debug()) {
871:               std::cout << msgLog << std::endl;
872:             }
873:             delete [] msgLog;
874: #endif
875:           }
876:           if (adjacency) delete [] adjacency;
877:           if (start)     delete [] start;
878:           ALE_LOG_EVENT_END;
879:           return assignment;
880:         };
881:       };
882:     };
883: #endif
884: #ifdef PETSC_HAVE_PARMETIS
885:     namespace ParMetis {
886:       template<typename Topology_>
887:       class Partitioner {
888:       public:
889:         typedef Topology_                          topology_type;
890:         typedef typename topology_type::sieve_type sieve_type;
891:         typedef typename topology_type::patch_type patch_type;
892:         typedef typename topology_type::point_type point_type;
893:         typedef int                                part_type;
894:       public:
897:         static part_type *partitionSieve(const Obj<topology_type>& topology, const int dim) {
898:           int    nvtxs      = 0;    // The number of vertices in full graph
899:           int   *vtxdist;           // Distribution of vertices across processes
900:           int   *xadj;              // Start of edge list for each vertex
901:           int   *adjncy;            // Edge lists for all vertices
902:           int   *vwgt       = NULL; // Vertex weights
903:           int   *adjwgt     = NULL; // Edge weights
904:           int    wgtflag    = 0;    // Indicates which weights are present
905:           int    numflag    = 0;    // Indicates initial offset (0 or 1)
906:           int    ncon       = 1;    // The number of weights per vertex
907:           int    nparts     = topology->commSize(); // The number of partitions
908:           float *tpwgts;            // The fraction of vertex weights assigned to each partition
909:           float *ubvec;             // The balance intolerance for vertex weights
910:           int    options[5];        // Options
911:           // Outputs
912:           int    edgeCut;           // The number of edges cut by the partition
913:           int   *assignment = NULL; // The vertex partition
914:           const typename topology_type::patch_type patch = 0;

916:           options[0] = 0; // Use all defaults
917:           vtxdist    = new int[nparts+1];
918:           vtxdist[0] = 0;
919:           tpwgts     = new float[ncon*nparts];
920:           for(int p = 0; p < nparts; ++p) {
921:             tpwgts[p] = 1.0/nparts;
922:           }
923:           ubvec      = new float[ncon];
924:           ubvec[0]   = 1.05;
925:           if (topology->hasPatch(patch)) {
926:             nvtxs      = topology->heightStratum(patch, 0)->size();
927:             assignment = new part_type[nvtxs];
928:           }
929:           MPI_Allgather(&nvtxs, 1, MPI_INT, &vtxdist[1], 1, MPI_INT, topology->comm());
930:           for(int p = 2; p <= nparts; ++p) {
931:             vtxdist[p] += vtxdist[p-1];
932:           }
933:           if (topology->commSize() == 1) {
934:             PetscMemzero(assignment, nvtxs * sizeof(part_type));
935:           } else {
936:             ALECompat::New::Partitioner<topology_type>::buildDualCSR(topology, dim, patch, &xadj, &adjncy);

938:             if (topology->debug() && nvtxs) {
939:               for(int p = 0; p <= nvtxs; ++p) {
940:                 std::cout << "["<<topology->commRank()<<"]xadj["<<p<<"] = " << xadj[p] << std::endl;
941:               }
942:               for(int i = 0; i < xadj[nvtxs]; ++i) {
943:                 std::cout << "["<<topology->commRank()<<"]adjncy["<<i<<"] = " << adjncy[i] << std::endl;
944:               }
945:             }
946:             if (vtxdist[1] == vtxdist[nparts]) {
947:               if (topology->commRank() == 0) {
948:                 METIS_PartGraphKway(&nvtxs, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &nparts, options, &edgeCut, assignment);
949:                 if (topology->debug()) {std::cout << "Metis: edgecut is " << edgeCut << std::endl;}
950:               }
951:             } else {
952:               MPI_Comm comm = topology->comm();

954:               ParMETIS_V3_PartKway(vtxdist, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgeCut, assignment, &comm);
955:               if (topology->debug()) {std::cout << "ParMetis: edgecut is " << edgeCut << std::endl;}
956:             }
957:             if (xadj   != NULL) delete [] xadj;
958:             if (adjncy != NULL) delete [] adjncy;
959:           }
960:           delete [] vtxdist;
961:           delete [] tpwgts;
962:           delete [] ubvec;
963:           return assignment;
964:         };
965:       };
966:     };
967: #endif
968:   }
969: }
970: #endif