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: #ifndef included_ALE_ZOLTAN_hh
9: #include <src/dm/mesh/meshzoltan.h>
10: #endif
12: #ifdef PETSC_HAVE_CHACO
13: /* Chaco does not have an include file */
16: float *ewgts, float *x, float *y, float *z, char *outassignname,
17: char *outfilename, short *assignment, int architecture, int ndims_tot,
18: int mesh_dims[3], double *goal, int global_method, int local_method,
19: int rqi_flag, int vmax, int ndims, double eigtol, long seed);
22: }
23: #endif
24: #ifdef PETSC_HAVE_PARMETIS
26: #include <parmetis.h>
28: }
29: #endif
30: #ifdef PETSC_HAVE_HMETIS
33: }
34: #endif
36: namespace ALE {
37: namespace New {
38: template<typename Bundle_>
39: class Partitioner {
40: public:
41: typedef Bundle_ bundle_type;
42: typedef typename bundle_type::sieve_type sieve_type;
43: typedef typename bundle_type::point_type point_type;
44: public:
47: // This creates a CSR representation of the adjacency matrix for cells
48: static void buildDualCSR(const Obj<bundle_type>& bundle, const int dim, int **offsets, int **adjacency) {
49: ALE_LOG_EVENT_BEGIN;
50: typedef typename ALE::New::Completion<bundle_type, point_type> completion;
51: const Obj<sieve_type>& sieve = bundle->getSieve();
52: const Obj<typename bundle_type::label_sequence>& elements = bundle->heightStratum(0);
53: Obj<sieve_type> overlapSieve = new sieve_type(bundle->comm(), bundle->debug());
54: int numElements = elements->size();
55: int *off = new int[numElements+1];
56: int offset = 0;
57: int *adj;
59: completion::scatterSupports(sieve, overlapSieve, bundle->getSendOverlap(), bundle->getRecvOverlap(), bundle);
60: if (numElements == 0) {
61: *offsets = NULL;
62: *adjacency = NULL;
63: ALE_LOG_EVENT_END;
64: return;
65: }
66: if (bundle->depth() == dim) {
67: int e = 1;
69: off[0] = 0;
70: for(typename bundle_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
71: const Obj<typename sieve_type::traits::coneSequence>& faces = sieve->cone(*e_iter);
72: typename sieve_type::traits::coneSequence::iterator fBegin = faces->begin();
73: typename sieve_type::traits::coneSequence::iterator fEnd = faces->end();
75: off[e] = off[e-1];
76: for(typename sieve_type::traits::coneSequence::iterator f_iter = fBegin; f_iter != fEnd; ++f_iter) {
77: if (sieve->support(*f_iter)->size() == 2) {
78: off[e]++;
79: } else if ((sieve->support(*f_iter)->size() == 1) && (overlapSieve->support(*f_iter)->size() == 1)) {
80: off[e]++;
81: }
82: }
83: e++;
84: }
85: adj = new int[off[numElements]];
86: for(typename bundle_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
87: const Obj<typename sieve_type::traits::coneSequence>& faces = sieve->cone(*e_iter);
88: typename sieve_type::traits::coneSequence::iterator fBegin = faces->begin();
89: typename sieve_type::traits::coneSequence::iterator fEnd = faces->end();
91: for(typename sieve_type::traits::coneSequence::iterator f_iter = fBegin; f_iter != fEnd; ++f_iter) {
92: const Obj<typename sieve_type::traits::supportSequence>& neighbors = sieve->support(*f_iter);
93: typename sieve_type::traits::supportSequence::iterator nBegin = neighbors->begin();
94: typename sieve_type::traits::supportSequence::iterator nEnd = neighbors->end();
96: for(typename sieve_type::traits::supportSequence::iterator n_iter = nBegin; n_iter != nEnd; ++n_iter) {
97: if (*n_iter != *e_iter) adj[offset++] = *n_iter;
98: }
99: const Obj<typename sieve_type::traits::supportSequence>& oNeighbors = overlapSieve->support(*f_iter);
100: typename sieve_type::traits::supportSequence::iterator onBegin = oNeighbors->begin();
101: typename sieve_type::traits::supportSequence::iterator onEnd = oNeighbors->end();
103: for(typename sieve_type::traits::supportSequence::iterator n_iter = onBegin; n_iter != onEnd; ++n_iter) {
104: adj[offset++] = *n_iter;
105: }
106: }
107: }
108: } else if (bundle->depth() == 1) {
109: std::set<point_type> *neighborCells = new std::set<point_type>[numElements];
110: int corners = sieve->cone(*elements->begin())->size();
111: int faceVertices = -1;
113: if (corners == dim+1) {
114: faceVertices = dim;
115: } else if ((dim == 2) && (corners == 4)) {
116: faceVertices = 2;
117: } else if ((dim == 3) && (corners == 8)) {
118: faceVertices = 4;
119: } else {
120: throw ALE::Exception("Could not determine number of face vertices");
121: }
122: for(typename bundle_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
123: const Obj<typename sieve_type::traits::coneSequence>& vertices = sieve->cone(*e_iter);
124: typename sieve_type::traits::coneSequence::iterator vEnd = vertices->end();
126: for(typename sieve_type::traits::coneSequence::iterator v_iter = vertices->begin(); v_iter != vEnd; ++v_iter) {
127: const Obj<typename sieve_type::traits::supportSequence>& neighbors = sieve->support(*v_iter);
128: typename sieve_type::traits::supportSequence::iterator nEnd = neighbors->end();
130: for(typename sieve_type::traits::supportSequence::iterator n_iter = neighbors->begin(); n_iter != nEnd; ++n_iter) {
131: if (*e_iter == *n_iter) continue;
132: if ((int) sieve->meet(*e_iter, *n_iter)->size() == faceVertices) {
133: neighborCells[*e_iter].insert(*n_iter);
134: }
135: }
136: }
137: }
138: off[0] = 0;
139: for(int e = 1; e <= numElements; e++) {
140: off[e] = neighborCells[e-1].size() + off[e-1];
141: }
142: adj = new int[off[numElements]];
143: for(int e = 0; e < numElements; e++) {
144: for(typename std::set<point_type>::iterator n_iter = neighborCells[e].begin(); n_iter != neighborCells[e].end(); ++n_iter) {
145: adj[offset++] = *n_iter;
146: }
147: }
148: delete [] neighborCells;
149: } else {
150: throw ALE::Exception("Dual creation not defined for partially interpolated meshes");
151: }
152: if (offset != off[numElements]) {
153: ostringstream msg;
154: msg << "ERROR: Total number of neighbors " << offset << " does not match the offset array " << off[numElements];
155: throw ALE::Exception(msg.str().c_str());
156: }
157: *offsets = off;
158: *adjacency = adj;
159: ALE_LOG_EVENT_END;
160: };
163: // This creates a CSR representation of the adjacency hypergraph for faces
164: 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) {
165: ALE_LOG_EVENT_BEGIN;
166: const Obj<sieve_type>& sieve = bundle->getSieve();
167: const Obj<typename bundle_type::label_sequence>& elements = bundle->heightStratum(0);
168: int numElements = elements->size();
169: int *off = new int[numElements+1];
170: int e;
172: if (bundle->depth() != dim) {
173: throw ALE::Exception("Not yet implemented for non-interpolated meshes");
174: }
175: off[0] = 0;
176: e = 1;
177: for(typename bundle_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
178: off[e] = sieve->cone(*e_iter)->size() + off[e-1];
179: e++;
180: }
181: int *adj = new int[off[numElements]];
182: int offset = 0;
183: for(typename bundle_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
184: const Obj<typename sieve_type::traits::coneSequence>& faces = sieve->cone(*e_iter);
185: typename sieve_type::traits::coneSequence::iterator fEnd = faces->end();
187: for(typename sieve_type::traits::coneSequence::iterator f_iter = faces->begin(); f_iter != fEnd; ++f_iter) {
188: adj[offset++] = fNumbering->getIndex(*f_iter);
189: }
190: }
191: if (offset != off[numElements]) {
192: ostringstream msg;
193: msg << "ERROR: Total number of neighbors " << offset << " does not match the offset array " << off[numElements];
194: throw ALE::Exception(msg.str().c_str());
195: }
196: *numEdges = numElements;
197: *offsets = off;
198: *adjacency = adj;
199: ALE_LOG_EVENT_END;
200: };
201: template<typename PartitionType>
202: static PartitionType *subordinatePartition(const Obj<bundle_type>& bundle, int levels, const Obj<bundle_type>& subBundle, const PartitionType assignment[]) {
203: const Obj<typename bundle_type::numbering_type>& cNumbering = bundle->getFactory()->getLocalNumbering(bundle, bundle->depth());
204: const Obj<typename bundle_type::label_sequence>& cells = subBundle->heightStratum(0);
205: const Obj<typename bundle_type::numbering_type>& sNumbering = bundle->getFactory()->getLocalNumbering(subBundle, subBundle->depth());
206: const int numCells = cells->size();
207: PartitionType *subAssignment = new PartitionType[numCells];
209: if (levels != 1) {
210: throw ALE::Exception("Cannot calculate subordinate partition for any level separation other than 1");
211: } else {
212: const Obj<typename bundle_type::sieve_type>& sieve = bundle->getSieve();
213: const Obj<typename bundle_type::sieve_type>& subSieve = subBundle->getSieve();
214: Obj<typename bundle_type::sieve_type::coneSet> tmpSet = new typename bundle_type::sieve_type::coneSet();
215: Obj<typename bundle_type::sieve_type::coneSet> tmpSet2 = new typename bundle_type::sieve_type::coneSet();
217: for(typename bundle_type::label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
218: const Obj<typename bundle_type::sieve_type::coneSequence>& cone = subSieve->cone(*c_iter);
220: Obj<typename bundle_type::sieve_type::supportSet> cell = sieve->nJoin1(cone);
221: if (cell->size() != 1) {
222: std::cout << "Indeterminate subordinate partition for face " << *c_iter << std::endl;
223: for(typename bundle_type::sieve_type::supportSet::iterator s_iter = cell->begin(); s_iter != cell->end(); ++s_iter) {
224: std::cout << " cell " << *s_iter << std::endl;
225: }
226: // Could relax this to choosing the first one
227: throw ALE::Exception("Indeterminate subordinate partition");
228: }
229: subAssignment[sNumbering->getIndex(*c_iter)] = assignment[cNumbering->getIndex(*cell->begin())];
230: tmpSet->clear();
231: tmpSet2->clear();
232: }
233: }
234: return subAssignment;
235: };
236: };
237: #ifdef PETSC_HAVE_CHACO
238: namespace Chaco {
239: template<typename Bundle_>
240: class Partitioner {
241: public:
242: typedef Bundle_ bundle_type;
243: typedef typename bundle_type::sieve_type sieve_type;
244: typedef typename bundle_type::point_type point_type;
245: typedef short int part_type;
246: public:
249: static part_type *partitionSieve(const Obj<bundle_type>& bundle, const int dim) {
250: part_type *assignment = NULL; /* set number of each vtx (length n) */
251: int *start; /* start of edge list for each vertex */
252: int *adjacency; /* = adj -> j; edge list data */
254: ALE_LOG_EVENT_BEGIN;
255: ALE::New::Partitioner<bundle_type>::buildDualCSR(bundle, dim, &start, &adjacency);
256: if (bundle->commRank() == 0) {
257: /* arguments for Chaco library */
258: FREE_GRAPH = 0; /* Do not let Chaco free my memory */
259: int nvtxs; /* number of vertices in full graph */
260: int *vwgts = NULL; /* weights for all vertices */
261: float *ewgts = NULL; /* weights for all edges */
262: float *x = NULL, *y = NULL, *z = NULL; /* coordinates for inertial method */
263: char *outassignname = NULL; /* name of assignment output file */
264: char *outfilename = NULL; /* output file name */
265: int architecture = 1; /* 0 => hypercube, d => d-dimensional mesh */
266: int ndims_tot = 0; /* total number of cube dimensions to divide */
267: int mesh_dims[3]; /* dimensions of mesh of processors */
268: double *goal = NULL; /* desired set sizes for each set */
269: int global_method = 1; /* global partitioning algorithm */
270: int local_method = 1; /* local partitioning algorithm */
271: int rqi_flag = 0; /* should I use RQI/Symmlq eigensolver? */
272: int vmax = 200; /* how many vertices to coarsen down to? */
273: int ndims = 1; /* number of eigenvectors (2^d sets) */
274: double eigtol = 0.001; /* tolerance on eigenvectors */
275: long seed = 123636512; /* for random graph mutations */
278: nvtxs = bundle->heightStratum(0)->size();
279: mesh_dims[0] = bundle->commSize(); mesh_dims[1] = 1; mesh_dims[2] = 1;
280: for(int e = 0; e < start[nvtxs]; e++) {
281: adjacency[e]++;
282: }
283: assignment = new part_type[nvtxs];
284: PetscMemzero(assignment, nvtxs * sizeof(part_type));CHKERROR(ierr, "Error in PetscMemzero");
286: /* redirect output to buffer: chaco -> msgLog */
287: #ifdef PETSC_HAVE_UNISTD_H
288: char *msgLog;
289: int fd_stdout, fd_pipe[2], count;
291: fd_stdout = dup(1);
292: pipe(fd_pipe);
293: close(1);
294: dup2(fd_pipe[1], 1);
295: msgLog = new char[16284];
296: #endif
298: interface(nvtxs, start, adjacency, vwgts, ewgts, x, y, z,
299: outassignname, outfilename, assignment, architecture, ndims_tot,
300: mesh_dims, goal, global_method, local_method, rqi_flag, vmax, ndims,
301: eigtol, seed);
303: #ifdef PETSC_HAVE_UNISTD_H
304: int SIZE_LOG = 10000;
306: fflush(stdout);
307: count = read(fd_pipe[0], msgLog, (SIZE_LOG - 1) * sizeof(char));
308: if (count < 0) count = 0;
309: msgLog[count] = 0;
310: close(1);
311: dup2(fd_stdout, 1);
312: close(fd_stdout);
313: close(fd_pipe[0]);
314: close(fd_pipe[1]);
315: if (bundle->debug()) {
316: std::cout << msgLog << std::endl;
317: }
318: delete [] msgLog;
319: #endif
320: }
321: if (adjacency) delete [] adjacency;
322: if (start) delete [] start;
323: ALE_LOG_EVENT_END;
324: return assignment;
325: };
326: static part_type *partitionSieveByFace(const Obj<bundle_type>& bundle, const int dim) {
327: throw ALE::Exception("Chaco cannot partition a mesh by faces");
328: };
329: };
330: };
331: #endif
332: #ifdef PETSC_HAVE_PARMETIS
333: namespace ParMetis {
334: template<typename Bundle_>
335: class Partitioner {
336: public:
337: typedef Bundle_ bundle_type;
338: typedef typename bundle_type::sieve_type sieve_type;
339: typedef typename bundle_type::point_type point_type;
340: typedef int part_type;
341: public:
344: static part_type *partitionSieve(const Obj<bundle_type>& bundle, const int dim) {
345: int nvtxs = 0; // The number of vertices in full graph
346: int *vtxdist; // Distribution of vertices across processes
347: int *xadj; // Start of edge list for each vertex
348: int *adjncy; // Edge lists for all vertices
349: int *vwgt = NULL; // Vertex weights
350: int *adjwgt = NULL; // Edge weights
351: int wgtflag = 0; // Indicates which weights are present
352: int numflag = 0; // Indicates initial offset (0 or 1)
353: int ncon = 1; // The number of weights per vertex
354: int nparts = bundle->commSize(); // The number of partitions
355: float *tpwgts; // The fraction of vertex weights assigned to each partition
356: float *ubvec; // The balance intolerance for vertex weights
357: int options[5]; // Options
358: // Outputs
359: int edgeCut; // The number of edges cut by the partition
360: int *assignment = NULL; // The vertex partition
362: options[0] = 0; // Use all defaults
363: vtxdist = new int[nparts+1];
364: vtxdist[0] = 0;
365: tpwgts = new float[ncon*nparts];
366: for(int p = 0; p < nparts; ++p) {
367: tpwgts[p] = 1.0/nparts;
368: }
369: ubvec = new float[ncon];
370: ubvec[0] = 1.05;
371: nvtxs = bundle->heightStratum(0)->size();
372: assignment = new part_type[nvtxs];
373: MPI_Allgather(&nvtxs, 1, MPI_INT, &vtxdist[1], 1, MPI_INT, bundle->comm());
374: for(int p = 2; p <= nparts; ++p) {
375: vtxdist[p] += vtxdist[p-1];
376: }
377: if (bundle->commSize() == 1) {
378: PetscMemzero(assignment, nvtxs * sizeof(part_type));
379: } else {
380: ALE::New::Partitioner<bundle_type>::buildDualCSR(bundle, dim, &xadj, &adjncy);
382: if (bundle->debug() && nvtxs) {
383: for(int p = 0; p <= nvtxs; ++p) {
384: std::cout << "["<<bundle->commRank()<<"]xadj["<<p<<"] = " << xadj[p] << std::endl;
385: }
386: for(int i = 0; i < xadj[nvtxs]; ++i) {
387: std::cout << "["<<bundle->commRank()<<"]adjncy["<<i<<"] = " << adjncy[i] << std::endl;
388: }
389: }
390: if (vtxdist[1] == vtxdist[nparts]) {
391: if (bundle->commRank() == 0) {
392: METIS_PartGraphKway(&nvtxs, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &nparts, options, &edgeCut, assignment);
393: if (bundle->debug()) {std::cout << "Metis: edgecut is " << edgeCut << std::endl;}
394: }
395: } else {
396: MPI_Comm comm = bundle->comm();
398: ParMETIS_V3_PartKway(vtxdist, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgeCut, assignment, &comm);
399: if (bundle->debug()) {std::cout << "ParMetis: edgecut is " << edgeCut << std::endl;}
400: }
401: if (xadj != NULL) delete [] xadj;
402: if (adjncy != NULL) delete [] adjncy;
403: }
404: delete [] vtxdist;
405: delete [] tpwgts;
406: delete [] ubvec;
407: return assignment;
408: };
411: static part_type *partitionSieveByFace(const Obj<bundle_type>& bundle, const int dim) {
412: #ifdef PETSC_HAVE_HMETIS
413: int *assignment = NULL; // The vertex partition
414: int nvtxs; // The number of vertices
415: int nhedges; // The number of hyperedges
416: int *vwgts; // The vertex weights
417: int *eptr; // The offsets of each hyperedge
418: int *eind; // The vertices in each hyperedge, indexed by eptr
419: int *hewgts; // The hyperedge weights
420: int nparts; // The number of partitions
421: int ubfactor; // The allowed load imbalance (1-50)
422: int options[9]; // Options
423: // Outputs
424: int edgeCut; // The number of edges cut by the partition
425: const Obj<ALE::Mesh::numbering_type>& fNumbering = bundle->getFactory()->getNumbering(bundle, bundle->depth()-1);
427: if (topology->commRank() == 0) {
428: nvtxs = bundle->heightStratum(1)->size();
429: vwgts = NULL;
430: hewgts = NULL;
431: nparts = bundle->commSize();
432: ubfactor = 5;
433: options[0] = 1; // Use all defaults
434: options[1] = 10; // Number of bisections tested
435: options[2] = 1; // Vertex grouping scheme
436: options[3] = 1; // Objective function
437: options[4] = 1; // V-cycle refinement
438: options[5] = 0;
439: options[6] = 0;
440: options[7] = 1; // Random seed
441: options[8] = 24; // Debugging level
442: assignment = new part_type[nvtxs];
444: if (bundle->commSize() == 1) {
445: PetscMemzero(assignment, nvtxs * sizeof(part_type));
446: } else {
447: ALE::New::Partitioner<bundle_type>::buildFaceCSR(bundle, dim, fNumbering, &nhedges, &eptr, &eind);
448: HMETIS_PartKway(nvtxs, nhedges, vwgts, eptr, eind, hewgts, nparts, ubfactor, options, assignment, &edgeCut);
450: delete [] eptr;
451: delete [] eind;
452: }
453: if (bundle->debug()) {for (int i = 0; i<nvtxs; i++) printf("[%d] %d\n", PetscGlobalRank, assignment[i]);}
454: } else {
455: assignment = NULL;
456: }
457: return assignment;
458: #else
459: throw ALE::Exception("hmetis partitioner is not available.");
460: #endif
461: };
462: };
463: };
464: #endif
465: #ifdef PETSC_HAVE_ZOLTAN
466: namespace Zoltan {
467: template<typename Bundle_>
468: class Partitioner {
469: public:
470: typedef Bundle_ bundle_type;
471: typedef typename bundle_type::sieve_type sieve_type;
472: typedef typename bundle_type::point_type point_type;
473: typedef int part_type;
474: public:
475: static part_type *partitionSieve(const Obj<bundle_type>& bundle, const int dim) {
476: throw ALE::Exception("Zoltan partition by cells not implemented");
477: };
480: static part_type *partitionSieveByFace(const Obj<bundle_type>& bundle, const int dim) {
481: // Outputs
482: float version; // The library version
483: int changed; // Did the partition change?
484: int numGidEntries; // Number of array entries for a single global ID (1)
485: int numLidEntries; // Number of array entries for a single local ID (1)
486: int numImport; // The number of imported points
487: ZOLTAN_ID_PTR import_global_ids; // The imported points
488: ZOLTAN_ID_PTR import_local_ids; // The imported points
489: int *import_procs; // The proc each point was imported from
490: int *import_to_part; // The partition of each imported point
491: int numExport; // The number of exported points
492: ZOLTAN_ID_PTR export_global_ids; // The exported points
493: ZOLTAN_ID_PTR export_local_ids; // The exported points
494: int *export_procs; // The proc each point was exported to
495: int *export_to_part; // The partition assignment of all local points
496: int *assignment; // The partition assignment of all local points
497: const Obj<typename bundle_type::numbering_type>& fNumbering = bundle->getFactory()->getNumbering(bundle, bundle->depth()-1);
499: if (bundle->commSize() == 1) {
500: PetscMemzero(assignment, bundle->heightStratum(1)->size() * sizeof(part_type));
501: } else {
502: if (bundle->commRank() == 0) {
503: nvtxs_Zoltan = bundle->heightStratum(1)->size();
504: ALE::New::Partitioner<bundle_type>::buildFaceCSR(bundle, dim, fNumbering, &nhedges_Zoltan, &eptr_Zoltan, &eind_Zoltan);
505: assignment = new int[nvtxs_Zoltan];
506: } else {
507: nvtxs_Zoltan = bundle->heightStratum(1)->size();
508: nhedges_Zoltan = 0;
509: eptr_Zoltan = new int[1];
510: eind_Zoltan = new int[1];
511: eptr_Zoltan[0] = 0;
512: assignment = NULL;
513: }
515: int Zoltan_Initialize(0, NULL, &version);
516: struct Zoltan_Struct *zz = Zoltan_Create(bundle->comm());
517: // General parameters
518: Zoltan_Set_Param(zz, "DEBUG_LEVEL", "2");
519: Zoltan_Set_Param(zz, "LB_METHOD", "PHG");
520: Zoltan_Set_Param(zz, "RETURN_LISTS", "PARTITION");
521: // PHG parameters
522: Zoltan_Set_Param(zz, "PHG_OUTPUT_LEVEL", "2");
523: Zoltan_Set_Param(zz, "PHG_EDGE_SIZE_THRESHOLD", "1.0"); // Do not throw out dense edges
524: // Call backs
525: Zoltan_Set_Num_Obj_Fn(zz, getNumVertices_Zoltan, NULL);
526: Zoltan_Set_Obj_List_Fn(zz, getLocalElements_Zoltan, NULL);
527: Zoltan_Set_HG_Size_CS_Fn(zz, getHgSizes_Zoltan, NULL);
528: Zoltan_Set_HG_CS_Fn(zz, getHg_Zoltan, NULL);
529: // Debugging
530: //Zoltan_Generate_Files(zz, "zoltan.debug", 1, 0, 0, 1); // if using hypergraph callbacks
532: Zoltan_LB_Partition(zz, &changed, &numGidEntries, &numLidEntries,
533: &numImport, &import_global_ids, &import_local_ids, &import_procs, &import_to_part,
534: &numExport, &export_global_ids, &export_local_ids, &export_procs, &export_to_part);
535: for(int v = 0; v < nvtxs_Zoltan; ++v) {
536: assignment[v] = export_to_part[v];
537: }
538: Zoltan_LB_Free_Part(&import_global_ids, &import_local_ids, &import_procs, &import_to_part);
539: Zoltan_LB_Free_Part(&export_global_ids, &export_local_ids, &export_procs, &export_to_part);
540: Zoltan_Destroy(&zz);
542: delete [] eptr_Zoltan;
543: delete [] eind_Zoltan;
544: }
545: if (assignment) {for (int i=0; i<nvtxs_Zoltan; i++) printf("[%d] %d\n",PetscGlobalRank,assignment[i]);}
546: return assignment;
547: };
548: };
549: };
550: #endif
551: }
552: }
554: namespace ALECompat {
555: namespace New {
556: template<typename Topology_>
557: class Partitioner {
558: public:
559: typedef Topology_ topology_type;
560: typedef typename topology_type::sieve_type sieve_type;
561: typedef typename topology_type::patch_type patch_type;
562: typedef typename topology_type::point_type point_type;
563: public:
566: // This creates a CSR representation of the adjacency matrix for cells
567: static void buildDualCSR(const Obj<topology_type>& topology, const int dim, const patch_type& patch, int **offsets, int **adjacency) {
568: ALE_LOG_EVENT_BEGIN;
569: typedef typename ALECompat::New::Completion<topology_type, typename Mesh::sieve_type::point_type> completion;
570: const Obj<sieve_type>& sieve = topology->getPatch(patch);
571: const Obj<typename topology_type::label_sequence>& elements = topology->heightStratum(patch, 0);
572: Obj<sieve_type> overlapSieve = new sieve_type(topology->comm(), topology->debug());
573: int numElements = elements->size();
574: int *off = new int[numElements+1];
575: int offset = 0;
576: int *adj;
578: completion::scatterSupports(sieve, overlapSieve, topology->getSendOverlap(), topology->getRecvOverlap(), topology);
579: if (numElements == 0) {
580: *offsets = NULL;
581: *adjacency = NULL;
582: ALE_LOG_EVENT_END;
583: return;
584: }
585: if (topology->depth(patch) == dim) {
586: int e = 1;
588: off[0] = 0;
589: for(typename topology_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
590: const Obj<typename sieve_type::traits::coneSequence>& faces = sieve->cone(*e_iter);
591: typename sieve_type::traits::coneSequence::iterator fBegin = faces->begin();
592: typename sieve_type::traits::coneSequence::iterator fEnd = faces->end();
594: off[e] = off[e-1];
595: for(typename sieve_type::traits::coneSequence::iterator f_iter = fBegin; f_iter != fEnd; ++f_iter) {
596: if (sieve->support(*f_iter)->size() == 2) {
597: off[e]++;
598: } else if ((sieve->support(*f_iter)->size() == 1) && (overlapSieve->support(*f_iter)->size() == 1)) {
599: off[e]++;
600: }
601: }
602: e++;
603: }
604: adj = new int[off[numElements]];
605: for(typename topology_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
606: const Obj<typename sieve_type::traits::coneSequence>& faces = sieve->cone(*e_iter);
607: typename sieve_type::traits::coneSequence::iterator fBegin = faces->begin();
608: typename sieve_type::traits::coneSequence::iterator fEnd = faces->end();
610: for(typename sieve_type::traits::coneSequence::iterator f_iter = fBegin; f_iter != fEnd; ++f_iter) {
611: const Obj<typename sieve_type::traits::supportSequence>& neighbors = sieve->support(*f_iter);
612: typename sieve_type::traits::supportSequence::iterator nBegin = neighbors->begin();
613: typename sieve_type::traits::supportSequence::iterator nEnd = neighbors->end();
615: for(typename sieve_type::traits::supportSequence::iterator n_iter = nBegin; n_iter != nEnd; ++n_iter) {
616: if (*n_iter != *e_iter) adj[offset++] = *n_iter;
617: }
618: const Obj<typename sieve_type::traits::supportSequence>& oNeighbors = overlapSieve->support(*f_iter);
619: typename sieve_type::traits::supportSequence::iterator onBegin = oNeighbors->begin();
620: typename sieve_type::traits::supportSequence::iterator onEnd = oNeighbors->end();
622: for(typename sieve_type::traits::supportSequence::iterator n_iter = onBegin; n_iter != onEnd; ++n_iter) {
623: adj[offset++] = *n_iter;
624: }
625: }
626: }
627: } else if (topology->depth(patch) == 1) {
628: std::set<point_type> *neighborCells = new std::set<point_type>[numElements];
629: int corners = sieve->cone(*elements->begin())->size();
630: int faceVertices = -1;
632: if (corners == dim+1) {
633: faceVertices = dim;
634: } else if ((dim == 2) && (corners == 4)) {
635: faceVertices = 2;
636: } else if ((dim == 3) && (corners == 8)) {
637: faceVertices = 4;
638: } else {
639: throw ALE::Exception("Could not determine number of face vertices");
640: }
641: for(typename topology_type::label_sequence::iterator e_iter = elements->begin(); e_iter != elements->end(); ++e_iter) {
642: const Obj<typename sieve_type::traits::coneSequence>& vertices = sieve->cone(*e_iter);
643: typename sieve_type::traits::coneSequence::iterator vEnd = vertices->end();
645: for(typename sieve_type::traits::coneSequence::iterator v_iter = vertices->begin(); v_iter != vEnd; ++v_iter) {
646: const Obj<typename sieve_type::traits::supportSequence>& neighbors = sieve->support(*v_iter);
647: typename sieve_type::traits::supportSequence::iterator nEnd = neighbors->end();
649: for(typename sieve_type::traits::supportSequence::iterator n_iter = neighbors->begin(); n_iter != nEnd; ++n_iter) {
650: if (*e_iter == *n_iter) continue;
651: if ((int) sieve->meet(*e_iter, *n_iter)->size() == faceVertices) {
652: neighborCells[*e_iter].insert(*n_iter);
653: }
654: }
655: }
656: }
657: off[0] = 0;
658: for(int e = 1; e <= numElements; e++) {
659: off[e] = neighborCells[e-1].size() + off[e-1];
660: }
661: adj = new int[off[numElements]];
662: for(int e = 0; e < numElements; e++) {
663: for(typename std::set<point_type>::iterator n_iter = neighborCells[e].begin(); n_iter != neighborCells[e].end(); ++n_iter) {
664: adj[offset++] = *n_iter;
665: }
666: }
667: delete [] neighborCells;
668: } else {
669: throw ALE::Exception("Dual creation not defined for partially interpolated meshes");
670: }
671: if (offset != off[numElements]) {
672: ostringstream msg;
673: msg << "ERROR: Total number of neighbors " << offset << " does not match the offset array " << off[numElements];
674: throw ALE::Exception(msg.str().c_str());
675: }
676: *offsets = off;
677: *adjacency = adj;
678: ALE_LOG_EVENT_END;
679: };
680: template<typename PartitionType>
681: static PartitionType *subordinatePartition(const Obj<topology_type>& topology, int levels, const Obj<topology_type>& subTopology, const PartitionType assignment[]) {
682: typedef ALECompat::New::NumberingFactory<topology_type> NumberingFactory;
683: const patch_type patch = 0;
684: const Obj<typename NumberingFactory::numbering_type>& cNumbering = NumberingFactory::singleton(topology->debug())->getLocalNumbering(topology, patch, topology->depth(patch));
685: const Obj<typename topology_type::label_sequence>& cells = subTopology->heightStratum(patch, 0);
686: const Obj<typename NumberingFactory::numbering_type>& sNumbering = NumberingFactory::singleton(subTopology->debug())->getLocalNumbering(subTopology, patch, subTopology->depth(patch));
687: const int numCells = cells->size();
688: PartitionType *subAssignment = new PartitionType[numCells];
690: if (levels != 1) {
691: throw ALE::Exception("Cannot calculate subordinate partition for any level separation other than 1");
692: } else {
693: const Obj<typename topology_type::sieve_type>& sieve = topology->getPatch(patch);
694: const Obj<typename topology_type::sieve_type>& subSieve = subTopology->getPatch(patch);
695: Obj<typename topology_type::sieve_type::coneSet> tmpSet = new typename topology_type::sieve_type::coneSet();
696: Obj<typename topology_type::sieve_type::coneSet> tmpSet2 = new typename topology_type::sieve_type::coneSet();
698: for(typename topology_type::label_sequence::iterator c_iter = cells->begin(); c_iter != cells->end(); ++c_iter) {
699: const Obj<typename topology_type::sieve_type::coneSequence>& cone = subSieve->cone(*c_iter);
701: Obj<typename topology_type::sieve_type::supportSet> cell = sieve->nJoin1(cone);
702: if (cell->size() != 1) {
703: std::cout << "Indeterminate subordinate partition for face " << *c_iter << std::endl;
704: for(typename topology_type::sieve_type::supportSet::iterator s_iter = cell->begin(); s_iter != cell->end(); ++s_iter) {
705: std::cout << " cell " << *s_iter << std::endl;
706: }
707: // Could relax this to choosing the first one
708: throw ALE::Exception("Indeterminate subordinate partition");
709: }
710: subAssignment[sNumbering->getIndex(*c_iter)] = assignment[cNumbering->getIndex(*cell->begin())];
711: tmpSet->clear();
712: tmpSet2->clear();
713: }
714: }
715: return subAssignment;
716: };
717: };
718: #ifdef PETSC_HAVE_CHACO
719: namespace Chaco {
720: template<typename Topology_>
721: class Partitioner {
722: public:
723: typedef Topology_ topology_type;
724: typedef typename topology_type::sieve_type sieve_type;
725: typedef typename topology_type::patch_type patch_type;
726: typedef typename topology_type::point_type point_type;
727: typedef short int part_type;
728: public:
731: static part_type *partitionSieve(const Obj<topology_type>& topology, const int dim) {
732: part_type *assignment = NULL; /* set number of each vtx (length n) */
733: int *start; /* start of edge list for each vertex */
734: int *adjacency; /* = adj -> j; edge list data */
735: typename topology_type::patch_type patch = 0;
737: ALE_LOG_EVENT_BEGIN;
738: ALECompat::New::Partitioner<topology_type>::buildDualCSR(topology, dim, patch, &start, &adjacency);
739: if (topology->commRank() == 0) {
740: /* arguments for Chaco library */
741: FREE_GRAPH = 0; /* Do not let Chaco free my memory */
742: int nvtxs; /* number of vertices in full graph */
743: int *vwgts = NULL; /* weights for all vertices */
744: float *ewgts = NULL; /* weights for all edges */
745: float *x = NULL, *y = NULL, *z = NULL; /* coordinates for inertial method */
746: char *outassignname = NULL; /* name of assignment output file */
747: char *outfilename = NULL; /* output file name */
748: int architecture = 1; /* 0 => hypercube, d => d-dimensional mesh */
749: int ndims_tot = 0; /* total number of cube dimensions to divide */
750: int mesh_dims[3]; /* dimensions of mesh of processors */
751: double *goal = NULL; /* desired set sizes for each set */
752: int global_method = 1; /* global partitioning algorithm */
753: int local_method = 1; /* local partitioning algorithm */
754: int rqi_flag = 0; /* should I use RQI/Symmlq eigensolver? */
755: int vmax = 200; /* how many vertices to coarsen down to? */
756: int ndims = 1; /* number of eigenvectors (2^d sets) */
757: double eigtol = 0.001; /* tolerance on eigenvectors */
758: long seed = 123636512; /* for random graph mutations */
761: nvtxs = topology->heightStratum(patch, 0)->size();
762: mesh_dims[0] = topology->commSize(); mesh_dims[1] = 1; mesh_dims[2] = 1;
763: for(int e = 0; e < start[nvtxs]; e++) {
764: adjacency[e]++;
765: }
766: assignment = new part_type[nvtxs];
767: PetscMemzero(assignment, nvtxs * sizeof(part_type));
769: /* redirect output to buffer: chaco -> msgLog */
770: #ifdef PETSC_HAVE_UNISTD_H
771: char *msgLog;
772: int fd_stdout, fd_pipe[2], count;
774: fd_stdout = dup(1);
775: pipe(fd_pipe);
776: close(1);
777: dup2(fd_pipe[1], 1);
778: msgLog = new char[16284];
779: #endif
781: interface(nvtxs, start, adjacency, vwgts, ewgts, x, y, z,
782: outassignname, outfilename, assignment, architecture, ndims_tot,
783: mesh_dims, goal, global_method, local_method, rqi_flag, vmax, ndims,
784: eigtol, seed);
786: #ifdef PETSC_HAVE_UNISTD_H
787: int SIZE_LOG = 10000;
789: fflush(stdout);
790: count = read(fd_pipe[0], msgLog, (SIZE_LOG - 1) * sizeof(char));
791: if (count < 0) count = 0;
792: msgLog[count] = 0;
793: close(1);
794: dup2(fd_stdout, 1);
795: close(fd_stdout);
796: close(fd_pipe[0]);
797: close(fd_pipe[1]);
798: if (topology->debug()) {
799: std::cout << msgLog << std::endl;
800: }
801: delete [] msgLog;
802: #endif
803: }
804: if (adjacency) delete [] adjacency;
805: if (start) delete [] start;
806: ALE_LOG_EVENT_END;
807: return assignment;
808: };
809: };
810: };
811: #endif
812: #ifdef PETSC_HAVE_PARMETIS
813: namespace ParMetis {
814: template<typename Topology_>
815: class Partitioner {
816: public:
817: typedef Topology_ topology_type;
818: typedef typename topology_type::sieve_type sieve_type;
819: typedef typename topology_type::patch_type patch_type;
820: typedef typename topology_type::point_type point_type;
821: typedef int part_type;
822: public:
825: static part_type *partitionSieve(const Obj<topology_type>& topology, const int dim) {
826: int nvtxs = 0; // The number of vertices in full graph
827: int *vtxdist; // Distribution of vertices across processes
828: int *xadj; // Start of edge list for each vertex
829: int *adjncy; // Edge lists for all vertices
830: int *vwgt = NULL; // Vertex weights
831: int *adjwgt = NULL; // Edge weights
832: int wgtflag = 0; // Indicates which weights are present
833: int numflag = 0; // Indicates initial offset (0 or 1)
834: int ncon = 1; // The number of weights per vertex
835: int nparts = topology->commSize(); // The number of partitions
836: float *tpwgts; // The fraction of vertex weights assigned to each partition
837: float *ubvec; // The balance intolerance for vertex weights
838: int options[5]; // Options
839: // Outputs
840: int edgeCut; // The number of edges cut by the partition
841: int *assignment = NULL; // The vertex partition
842: const typename topology_type::patch_type patch = 0;
844: options[0] = 0; // Use all defaults
845: vtxdist = new int[nparts+1];
846: vtxdist[0] = 0;
847: tpwgts = new float[ncon*nparts];
848: for(int p = 0; p < nparts; ++p) {
849: tpwgts[p] = 1.0/nparts;
850: }
851: ubvec = new float[ncon];
852: ubvec[0] = 1.05;
853: if (topology->hasPatch(patch)) {
854: nvtxs = topology->heightStratum(patch, 0)->size();
855: assignment = new part_type[nvtxs];
856: }
857: MPI_Allgather(&nvtxs, 1, MPI_INT, &vtxdist[1], 1, MPI_INT, topology->comm());
858: for(int p = 2; p <= nparts; ++p) {
859: vtxdist[p] += vtxdist[p-1];
860: }
861: if (topology->commSize() == 1) {
862: PetscMemzero(assignment, nvtxs * sizeof(part_type));
863: } else {
864: ALECompat::New::Partitioner<topology_type>::buildDualCSR(topology, dim, patch, &xadj, &adjncy);
866: if (topology->debug() && nvtxs) {
867: for(int p = 0; p <= nvtxs; ++p) {
868: std::cout << "["<<topology->commRank()<<"]xadj["<<p<<"] = " << xadj[p] << std::endl;
869: }
870: for(int i = 0; i < xadj[nvtxs]; ++i) {
871: std::cout << "["<<topology->commRank()<<"]adjncy["<<i<<"] = " << adjncy[i] << std::endl;
872: }
873: }
874: if (vtxdist[1] == vtxdist[nparts]) {
875: if (topology->commRank() == 0) {
876: METIS_PartGraphKway(&nvtxs, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &nparts, options, &edgeCut, assignment);
877: if (topology->debug()) {std::cout << "Metis: edgecut is " << edgeCut << std::endl;}
878: }
879: } else {
880: MPI_Comm comm = topology->comm();
882: ParMETIS_V3_PartKway(vtxdist, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgeCut, assignment, &comm);
883: if (topology->debug()) {std::cout << "ParMetis: edgecut is " << edgeCut << std::endl;}
884: }
885: if (xadj != NULL) delete [] xadj;
886: if (adjncy != NULL) delete [] adjncy;
887: }
888: delete [] vtxdist;
889: delete [] tpwgts;
890: delete [] ubvec;
891: return assignment;
892: };
893: };
894: };
895: #endif
896: }
897: }
898: #endif