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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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