xref: /petsc/src/dm/impls/plex/plexpartition.c (revision 9852e1235e6de1345fdd93ba94fdd27144d16762)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
270034214SMatthew G. Knepley 
377623264SMatthew G. Knepley PetscClassId PETSCPARTITIONER_CLASSID = 0;
477623264SMatthew G. Knepley 
577623264SMatthew G. Knepley PetscFunctionList PetscPartitionerList              = NULL;
677623264SMatthew G. Knepley PetscBool         PetscPartitionerRegisterAllCalled = PETSC_FALSE;
777623264SMatthew G. Knepley 
877623264SMatthew G. Knepley PetscBool ChacoPartitionercite = PETSC_FALSE;
977623264SMatthew G. Knepley const char ChacoPartitionerCitation[] = "@inproceedings{Chaco95,\n"
1077623264SMatthew G. Knepley                                "  author    = {Bruce Hendrickson and Robert Leland},\n"
1177623264SMatthew G. Knepley                                "  title     = {A multilevel algorithm for partitioning graphs},\n"
1277623264SMatthew G. Knepley                                "  booktitle = {Supercomputing '95: Proceedings of the 1995 ACM/IEEE Conference on Supercomputing (CDROM)},"
1377623264SMatthew G. Knepley                                "  isbn      = {0-89791-816-9},\n"
1477623264SMatthew G. Knepley                                "  pages     = {28},\n"
1577623264SMatthew G. Knepley                                "  doi       = {http://doi.acm.org/10.1145/224170.224228},\n"
1677623264SMatthew G. Knepley                                "  publisher = {ACM Press},\n"
1777623264SMatthew G. Knepley                                "  address   = {New York},\n"
1877623264SMatthew G. Knepley                                "  year      = {1995}\n}\n";
1977623264SMatthew G. Knepley 
2077623264SMatthew G. Knepley PetscBool ParMetisPartitionercite = PETSC_FALSE;
2177623264SMatthew G. Knepley const char ParMetisPartitionerCitation[] = "@article{KarypisKumar98,\n"
2277623264SMatthew G. Knepley                                "  author  = {George Karypis and Vipin Kumar},\n"
2377623264SMatthew G. Knepley                                "  title   = {A Parallel Algorithm for Multilevel Graph Partitioning and Sparse Matrix Ordering},\n"
2477623264SMatthew G. Knepley                                "  journal = {Journal of Parallel and Distributed Computing},\n"
2577623264SMatthew G. Knepley                                "  volume  = {48},\n"
2677623264SMatthew G. Knepley                                "  pages   = {71--85},\n"
2777623264SMatthew G. Knepley                                "  year    = {1998}\n}\n";
2877623264SMatthew G. Knepley 
29532c4e7dSMichael Lange /*@C
30532c4e7dSMichael Lange   DMPlexCreatePartitionerGraph - Create a CSR graph of point connections for the partitioner
31532c4e7dSMichael Lange 
32532c4e7dSMichael Lange   Input Parameters:
33532c4e7dSMichael Lange + dm      - The mesh DM dm
34532c4e7dSMichael Lange - height  - Height of the strata from which to construct the graph
35532c4e7dSMichael Lange 
36532c4e7dSMichael Lange   Output Parameter:
37532c4e7dSMichael Lange + numVertices     - Number of vertices in the graph
383fa7752dSToby Isaac . offsets         - Point offsets in the graph
393fa7752dSToby Isaac . adjacency       - Point connectivity in the graph
403fa7752dSToby Isaac - globalNumbering - A map from the local cell numbering to the global numbering used in "adjacency".  Negative indicates that the cell is a duplicate from another process.
41532c4e7dSMichael Lange 
42532c4e7dSMichael Lange   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
43532c4e7dSMichael Lange   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
44532c4e7dSMichael Lange   representation on the mesh.
45532c4e7dSMichael Lange 
46532c4e7dSMichael Lange   Level: developer
47532c4e7dSMichael Lange 
48532c4e7dSMichael Lange .seealso: PetscPartitionerGetType(), PetscPartitionerCreate(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
49532c4e7dSMichael Lange @*/
503fa7752dSToby Isaac PetscErrorCode DMPlexCreatePartitionerGraph(DM dm, PetscInt height, PetscInt *numVertices, PetscInt **offsets, PetscInt **adjacency, IS *globalNumbering)
51532c4e7dSMichael Lange {
5243f7d02bSMichael Lange   PetscInt       p, pStart, pEnd, a, adjSize, idx, size, nroots;
53389e55d8SMichael Lange   PetscInt      *adj = NULL, *vOffsets = NULL, *graph = NULL;
548cfe4c1fSMichael Lange   IS             cellNumbering;
558cfe4c1fSMichael Lange   const PetscInt *cellNum;
56532c4e7dSMichael Lange   PetscBool      useCone, useClosure;
57532c4e7dSMichael Lange   PetscSection   section;
58532c4e7dSMichael Lange   PetscSegBuffer adjBuffer;
598cfe4c1fSMichael Lange   PetscSF        sfPoint;
60532c4e7dSMichael Lange   PetscMPIInt    rank;
61532c4e7dSMichael Lange   PetscErrorCode ierr;
62532c4e7dSMichael Lange 
63532c4e7dSMichael Lange   PetscFunctionBegin;
64532c4e7dSMichael Lange   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
65532c4e7dSMichael Lange   ierr = DMPlexGetHeightStratum(dm, height, &pStart, &pEnd);CHKERRQ(ierr);
668cfe4c1fSMichael Lange   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
678cfe4c1fSMichael Lange   ierr = PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
68532c4e7dSMichael Lange   /* Build adjacency graph via a section/segbuffer */
69532c4e7dSMichael Lange   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &section);CHKERRQ(ierr);
70532c4e7dSMichael Lange   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
71532c4e7dSMichael Lange   ierr = PetscSegBufferCreate(sizeof(PetscInt),1000,&adjBuffer);CHKERRQ(ierr);
72532c4e7dSMichael Lange   /* Always use FVM adjacency to create partitioner graph */
73532c4e7dSMichael Lange   ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
74532c4e7dSMichael Lange   ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr);
75532c4e7dSMichael Lange   ierr = DMPlexSetAdjacencyUseCone(dm, PETSC_TRUE);CHKERRQ(ierr);
76532c4e7dSMichael Lange   ierr = DMPlexSetAdjacencyUseClosure(dm, PETSC_FALSE);CHKERRQ(ierr);
77f0927f4eSMatthew G. Knepley   ierr = DMPlexCreateCellNumbering_Internal(dm, PETSC_TRUE, &cellNumbering);CHKERRQ(ierr);
783fa7752dSToby Isaac   if (globalNumbering) {
793fa7752dSToby Isaac     ierr = PetscObjectReference((PetscObject)cellNumbering);CHKERRQ(ierr);
803fa7752dSToby Isaac     *globalNumbering = cellNumbering;
813fa7752dSToby Isaac   }
828cfe4c1fSMichael Lange   ierr = ISGetIndices(cellNumbering, &cellNum);CHKERRQ(ierr);
838cfe4c1fSMichael Lange   for (*numVertices = 0, p = pStart; p < pEnd; p++) {
848cfe4c1fSMichael Lange     /* Skip non-owned cells in parallel (ParMetis expects no overlap) */
858cfe4c1fSMichael Lange     if (nroots > 0) {if (cellNum[p] < 0) continue;}
86532c4e7dSMichael Lange     adjSize = PETSC_DETERMINE;
87532c4e7dSMichael Lange     ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
88532c4e7dSMichael Lange     for (a = 0; a < adjSize; ++a) {
89532c4e7dSMichael Lange       const PetscInt point = adj[a];
90532c4e7dSMichael Lange       if (point != p && pStart <= point && point < pEnd) {
91532c4e7dSMichael Lange         PetscInt *PETSC_RESTRICT pBuf;
92532c4e7dSMichael Lange         ierr = PetscSectionAddDof(section, p, 1);CHKERRQ(ierr);
93532c4e7dSMichael Lange         ierr = PetscSegBufferGetInts(adjBuffer, 1, &pBuf);CHKERRQ(ierr);
94532c4e7dSMichael Lange         *pBuf = point;
95532c4e7dSMichael Lange       }
96532c4e7dSMichael Lange     }
978cfe4c1fSMichael Lange     (*numVertices)++;
98532c4e7dSMichael Lange   }
99532c4e7dSMichael Lange   ierr = DMPlexSetAdjacencyUseCone(dm, useCone);CHKERRQ(ierr);
100532c4e7dSMichael Lange   ierr = DMPlexSetAdjacencyUseClosure(dm, useClosure);CHKERRQ(ierr);
101532c4e7dSMichael Lange   /* Derive CSR graph from section/segbuffer */
102532c4e7dSMichael Lange   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
103532c4e7dSMichael Lange   ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
104389e55d8SMichael Lange   ierr = PetscMalloc1(*numVertices+1, &vOffsets);CHKERRQ(ierr);
10543f7d02bSMichael Lange   for (idx = 0, p = pStart; p < pEnd; p++) {
10643f7d02bSMichael Lange     if (nroots > 0) {if (cellNum[p] < 0) continue;}
10743f7d02bSMichael Lange     ierr = PetscSectionGetOffset(section, p, &(vOffsets[idx++]));CHKERRQ(ierr);
10843f7d02bSMichael Lange   }
109532c4e7dSMichael Lange   vOffsets[*numVertices] = size;
110532c4e7dSMichael Lange   if (offsets) *offsets = vOffsets;
111389e55d8SMichael Lange   ierr = PetscSegBufferExtractAlloc(adjBuffer, &graph);CHKERRQ(ierr);
112bf4602e4SToby Isaac   {
1138cfe4c1fSMichael Lange     ISLocalToGlobalMapping ltogCells;
1148cfe4c1fSMichael Lange     PetscInt n, size, *cells_arr;
1158cfe4c1fSMichael Lange     /* In parallel, apply a global cell numbering to the graph */
1168cfe4c1fSMichael Lange     ierr = ISRestoreIndices(cellNumbering, &cellNum);CHKERRQ(ierr);
1178cfe4c1fSMichael Lange     ierr = ISLocalToGlobalMappingCreateIS(cellNumbering, &ltogCells);CHKERRQ(ierr);
1188cfe4c1fSMichael Lange     ierr = ISLocalToGlobalMappingGetSize(ltogCells, &size);CHKERRQ(ierr);
1198cfe4c1fSMichael Lange     ierr = ISLocalToGlobalMappingGetIndices(ltogCells, (const PetscInt**)&cells_arr);CHKERRQ(ierr);
1208cfe4c1fSMichael Lange     /* Convert to positive global cell numbers */
1218cfe4c1fSMichael Lange     for (n=0; n<size; n++) {if (cells_arr[n] < 0) cells_arr[n] = -(cells_arr[n]+1);}
1228cfe4c1fSMichael Lange     ierr = ISLocalToGlobalMappingRestoreIndices(ltogCells, (const PetscInt**)&cells_arr);CHKERRQ(ierr);
123389e55d8SMichael Lange     ierr = ISLocalToGlobalMappingApplyBlock(ltogCells, vOffsets[*numVertices], graph, graph);CHKERRQ(ierr);
1248cfe4c1fSMichael Lange     ierr = ISLocalToGlobalMappingDestroy(&ltogCells);CHKERRQ(ierr);
125f0927f4eSMatthew G. Knepley     ierr = ISDestroy(&cellNumbering);CHKERRQ(ierr);
1268cfe4c1fSMichael Lange   }
127389e55d8SMichael Lange   if (adjacency) *adjacency = graph;
128532c4e7dSMichael Lange   /* Clean up */
129532c4e7dSMichael Lange   ierr = PetscSegBufferDestroy(&adjBuffer);CHKERRQ(ierr);
130532c4e7dSMichael Lange   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
131532c4e7dSMichael Lange   ierr = PetscFree(adj);CHKERRQ(ierr);
132532c4e7dSMichael Lange   PetscFunctionReturn(0);
133532c4e7dSMichael Lange }
134532c4e7dSMichael Lange 
13570034214SMatthew G. Knepley PetscErrorCode DMPlexCreateNeighborCSR(DM dm, PetscInt cellHeight, PetscInt *numVertices, PetscInt **offsets, PetscInt **adjacency)
13670034214SMatthew G. Knepley {
13770034214SMatthew G. Knepley   const PetscInt maxFaceCases = 30;
13870034214SMatthew G. Knepley   PetscInt       numFaceCases = 0;
13970034214SMatthew G. Knepley   PetscInt       numFaceVertices[30]; /* maxFaceCases, C89 sucks sucks sucks */
14070034214SMatthew G. Knepley   PetscInt      *off, *adj;
14170034214SMatthew G. Knepley   PetscInt      *neighborCells = NULL;
14270034214SMatthew G. Knepley   PetscInt       dim, cellDim, depth = 0, faceDepth, cStart, cEnd, c, numCells, cell;
14370034214SMatthew G. Knepley   PetscErrorCode ierr;
14470034214SMatthew G. Knepley 
14570034214SMatthew G. Knepley   PetscFunctionBegin;
14670034214SMatthew G. Knepley   /* For parallel partitioning, I think you have to communicate supports */
147c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
14870034214SMatthew G. Knepley   cellDim = dim - cellHeight;
14970034214SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
15070034214SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
15170034214SMatthew G. Knepley   if (cEnd - cStart == 0) {
15270034214SMatthew G. Knepley     if (numVertices) *numVertices = 0;
15370034214SMatthew G. Knepley     if (offsets)   *offsets   = NULL;
15470034214SMatthew G. Knepley     if (adjacency) *adjacency = NULL;
15570034214SMatthew G. Knepley     PetscFunctionReturn(0);
15670034214SMatthew G. Knepley   }
15770034214SMatthew G. Knepley   numCells  = cEnd - cStart;
15870034214SMatthew G. Knepley   faceDepth = depth - cellHeight;
15970034214SMatthew G. Knepley   if (dim == depth) {
16070034214SMatthew G. Knepley     PetscInt f, fStart, fEnd;
16170034214SMatthew G. Knepley 
16270034214SMatthew G. Knepley     ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr);
16370034214SMatthew G. Knepley     /* Count neighboring cells */
16470034214SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, cellHeight+1, &fStart, &fEnd);CHKERRQ(ierr);
16570034214SMatthew G. Knepley     for (f = fStart; f < fEnd; ++f) {
16670034214SMatthew G. Knepley       const PetscInt *support;
16770034214SMatthew G. Knepley       PetscInt        supportSize;
16870034214SMatthew G. Knepley       ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr);
16970034214SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr);
17070034214SMatthew G. Knepley       if (supportSize == 2) {
1719ffc88e4SToby Isaac         PetscInt numChildren;
1729ffc88e4SToby Isaac 
1739ffc88e4SToby Isaac         ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr);
1749ffc88e4SToby Isaac         if (!numChildren) {
17570034214SMatthew G. Knepley           ++off[support[0]-cStart+1];
17670034214SMatthew G. Knepley           ++off[support[1]-cStart+1];
17770034214SMatthew G. Knepley         }
17870034214SMatthew G. Knepley       }
1799ffc88e4SToby Isaac     }
18070034214SMatthew G. Knepley     /* Prefix sum */
18170034214SMatthew G. Knepley     for (c = 1; c <= numCells; ++c) off[c] += off[c-1];
18270034214SMatthew G. Knepley     if (adjacency) {
18370034214SMatthew G. Knepley       PetscInt *tmp;
18470034214SMatthew G. Knepley 
18570034214SMatthew G. Knepley       ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr);
186854ce69bSBarry Smith       ierr = PetscMalloc1(numCells+1, &tmp);CHKERRQ(ierr);
18770034214SMatthew G. Knepley       ierr = PetscMemcpy(tmp, off, (numCells+1) * sizeof(PetscInt));CHKERRQ(ierr);
18870034214SMatthew G. Knepley       /* Get neighboring cells */
18970034214SMatthew G. Knepley       for (f = fStart; f < fEnd; ++f) {
19070034214SMatthew G. Knepley         const PetscInt *support;
19170034214SMatthew G. Knepley         PetscInt        supportSize;
19270034214SMatthew G. Knepley         ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr);
19370034214SMatthew G. Knepley         ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr);
19470034214SMatthew G. Knepley         if (supportSize == 2) {
1959ffc88e4SToby Isaac           PetscInt numChildren;
1969ffc88e4SToby Isaac 
1979ffc88e4SToby Isaac           ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr);
1989ffc88e4SToby Isaac           if (!numChildren) {
19970034214SMatthew G. Knepley             adj[tmp[support[0]-cStart]++] = support[1];
20070034214SMatthew G. Knepley             adj[tmp[support[1]-cStart]++] = support[0];
20170034214SMatthew G. Knepley           }
20270034214SMatthew G. Knepley         }
2039ffc88e4SToby Isaac       }
20470034214SMatthew G. Knepley #if defined(PETSC_USE_DEBUG)
20570034214SMatthew G. Knepley       for (c = 0; c < cEnd-cStart; ++c) if (tmp[c] != off[c+1]) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Offset %d != %d for cell %d", tmp[c], off[c], c+cStart);
20670034214SMatthew G. Knepley #endif
20770034214SMatthew G. Knepley       ierr = PetscFree(tmp);CHKERRQ(ierr);
20870034214SMatthew G. Knepley     }
20970034214SMatthew G. Knepley     if (numVertices) *numVertices = numCells;
21070034214SMatthew G. Knepley     if (offsets)   *offsets   = off;
21170034214SMatthew G. Knepley     if (adjacency) *adjacency = adj;
21270034214SMatthew G. Knepley     PetscFunctionReturn(0);
21370034214SMatthew G. Knepley   }
21470034214SMatthew G. Knepley   /* Setup face recognition */
21570034214SMatthew G. Knepley   if (faceDepth == 1) {
21670034214SMatthew G. Knepley     PetscInt cornersSeen[30] = {0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0,0}; /* Could use PetscBT */
21770034214SMatthew G. Knepley 
21870034214SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
21970034214SMatthew G. Knepley       PetscInt corners;
22070034214SMatthew G. Knepley 
22170034214SMatthew G. Knepley       ierr = DMPlexGetConeSize(dm, c, &corners);CHKERRQ(ierr);
22270034214SMatthew G. Knepley       if (!cornersSeen[corners]) {
22370034214SMatthew G. Knepley         PetscInt nFV;
22470034214SMatthew G. Knepley 
22570034214SMatthew G. Knepley         if (numFaceCases >= maxFaceCases) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Exceeded maximum number of face recognition cases");
22670034214SMatthew G. Knepley         cornersSeen[corners] = 1;
22770034214SMatthew G. Knepley 
22870034214SMatthew G. Knepley         ierr = DMPlexGetNumFaceVertices(dm, cellDim, corners, &nFV);CHKERRQ(ierr);
22970034214SMatthew G. Knepley 
23070034214SMatthew G. Knepley         numFaceVertices[numFaceCases++] = nFV;
23170034214SMatthew G. Knepley       }
23270034214SMatthew G. Knepley     }
23370034214SMatthew G. Knepley   }
23470034214SMatthew G. Knepley   ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr);
23570034214SMatthew G. Knepley   /* Count neighboring cells */
23670034214SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) {
23770034214SMatthew G. Knepley     PetscInt numNeighbors = PETSC_DETERMINE, n;
23870034214SMatthew G. Knepley 
2398b0b4c70SToby Isaac     ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &numNeighbors, &neighborCells);CHKERRQ(ierr);
24070034214SMatthew G. Knepley     /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */
24170034214SMatthew G. Knepley     for (n = 0; n < numNeighbors; ++n) {
24270034214SMatthew G. Knepley       PetscInt        cellPair[2];
24370034214SMatthew G. Knepley       PetscBool       found    = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE;
24470034214SMatthew G. Knepley       PetscInt        meetSize = 0;
24570034214SMatthew G. Knepley       const PetscInt *meet    = NULL;
24670034214SMatthew G. Knepley 
24770034214SMatthew G. Knepley       cellPair[0] = cell; cellPair[1] = neighborCells[n];
24870034214SMatthew G. Knepley       if (cellPair[0] == cellPair[1]) continue;
24970034214SMatthew G. Knepley       if (!found) {
25070034214SMatthew G. Knepley         ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
25170034214SMatthew G. Knepley         if (meetSize) {
25270034214SMatthew G. Knepley           PetscInt f;
25370034214SMatthew G. Knepley 
25470034214SMatthew G. Knepley           for (f = 0; f < numFaceCases; ++f) {
25570034214SMatthew G. Knepley             if (numFaceVertices[f] == meetSize) {
25670034214SMatthew G. Knepley               found = PETSC_TRUE;
25770034214SMatthew G. Knepley               break;
25870034214SMatthew G. Knepley             }
25970034214SMatthew G. Knepley           }
26070034214SMatthew G. Knepley         }
26170034214SMatthew G. Knepley         ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
26270034214SMatthew G. Knepley       }
26370034214SMatthew G. Knepley       if (found) ++off[cell-cStart+1];
26470034214SMatthew G. Knepley     }
26570034214SMatthew G. Knepley   }
26670034214SMatthew G. Knepley   /* Prefix sum */
26770034214SMatthew G. Knepley   for (cell = 1; cell <= numCells; ++cell) off[cell] += off[cell-1];
26870034214SMatthew G. Knepley 
26970034214SMatthew G. Knepley   if (adjacency) {
27070034214SMatthew G. Knepley     ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr);
27170034214SMatthew G. Knepley     /* Get neighboring cells */
27270034214SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
27370034214SMatthew G. Knepley       PetscInt numNeighbors = PETSC_DETERMINE, n;
27470034214SMatthew G. Knepley       PetscInt cellOffset   = 0;
27570034214SMatthew G. Knepley 
2768b0b4c70SToby Isaac       ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &numNeighbors, &neighborCells);CHKERRQ(ierr);
27770034214SMatthew G. Knepley       /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */
27870034214SMatthew G. Knepley       for (n = 0; n < numNeighbors; ++n) {
27970034214SMatthew G. Knepley         PetscInt        cellPair[2];
28070034214SMatthew G. Knepley         PetscBool       found    = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE;
28170034214SMatthew G. Knepley         PetscInt        meetSize = 0;
28270034214SMatthew G. Knepley         const PetscInt *meet    = NULL;
28370034214SMatthew G. Knepley 
28470034214SMatthew G. Knepley         cellPair[0] = cell; cellPair[1] = neighborCells[n];
28570034214SMatthew G. Knepley         if (cellPair[0] == cellPair[1]) continue;
28670034214SMatthew G. Knepley         if (!found) {
28770034214SMatthew G. Knepley           ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
28870034214SMatthew G. Knepley           if (meetSize) {
28970034214SMatthew G. Knepley             PetscInt f;
29070034214SMatthew G. Knepley 
29170034214SMatthew G. Knepley             for (f = 0; f < numFaceCases; ++f) {
29270034214SMatthew G. Knepley               if (numFaceVertices[f] == meetSize) {
29370034214SMatthew G. Knepley                 found = PETSC_TRUE;
29470034214SMatthew G. Knepley                 break;
29570034214SMatthew G. Knepley               }
29670034214SMatthew G. Knepley             }
29770034214SMatthew G. Knepley           }
29870034214SMatthew G. Knepley           ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
29970034214SMatthew G. Knepley         }
30070034214SMatthew G. Knepley         if (found) {
30170034214SMatthew G. Knepley           adj[off[cell-cStart]+cellOffset] = neighborCells[n];
30270034214SMatthew G. Knepley           ++cellOffset;
30370034214SMatthew G. Knepley         }
30470034214SMatthew G. Knepley       }
30570034214SMatthew G. Knepley     }
30670034214SMatthew G. Knepley   }
30770034214SMatthew G. Knepley   ierr = PetscFree(neighborCells);CHKERRQ(ierr);
30870034214SMatthew G. Knepley   if (numVertices) *numVertices = numCells;
30970034214SMatthew G. Knepley   if (offsets)   *offsets   = off;
31070034214SMatthew G. Knepley   if (adjacency) *adjacency = adj;
31170034214SMatthew G. Knepley   PetscFunctionReturn(0);
31270034214SMatthew G. Knepley }
31370034214SMatthew G. Knepley 
31477623264SMatthew G. Knepley /*@C
31577623264SMatthew G. Knepley   PetscPartitionerRegister - Adds a new PetscPartitioner implementation
31677623264SMatthew G. Knepley 
31777623264SMatthew G. Knepley   Not Collective
31877623264SMatthew G. Knepley 
31977623264SMatthew G. Knepley   Input Parameters:
32077623264SMatthew G. Knepley + name        - The name of a new user-defined creation routine
32177623264SMatthew G. Knepley - create_func - The creation routine itself
32277623264SMatthew G. Knepley 
32377623264SMatthew G. Knepley   Notes:
32477623264SMatthew G. Knepley   PetscPartitionerRegister() may be called multiple times to add several user-defined PetscPartitioners
32577623264SMatthew G. Knepley 
32677623264SMatthew G. Knepley   Sample usage:
32777623264SMatthew G. Knepley .vb
32877623264SMatthew G. Knepley     PetscPartitionerRegister("my_part", MyPetscPartitionerCreate);
32977623264SMatthew G. Knepley .ve
33077623264SMatthew G. Knepley 
33177623264SMatthew G. Knepley   Then, your PetscPartitioner type can be chosen with the procedural interface via
33277623264SMatthew G. Knepley .vb
33377623264SMatthew G. Knepley     PetscPartitionerCreate(MPI_Comm, PetscPartitioner *);
33477623264SMatthew G. Knepley     PetscPartitionerSetType(PetscPartitioner, "my_part");
33577623264SMatthew G. Knepley .ve
33677623264SMatthew G. Knepley    or at runtime via the option
33777623264SMatthew G. Knepley .vb
33877623264SMatthew G. Knepley     -petscpartitioner_type my_part
33977623264SMatthew G. Knepley .ve
34077623264SMatthew G. Knepley 
34177623264SMatthew G. Knepley   Level: advanced
34277623264SMatthew G. Knepley 
34377623264SMatthew G. Knepley .keywords: PetscPartitioner, register
34477623264SMatthew G. Knepley .seealso: PetscPartitionerRegisterAll(), PetscPartitionerRegisterDestroy()
34577623264SMatthew G. Knepley 
34677623264SMatthew G. Knepley @*/
34777623264SMatthew G. Knepley PetscErrorCode PetscPartitionerRegister(const char sname[], PetscErrorCode (*function)(PetscPartitioner))
34877623264SMatthew G. Knepley {
34977623264SMatthew G. Knepley   PetscErrorCode ierr;
35077623264SMatthew G. Knepley 
35177623264SMatthew G. Knepley   PetscFunctionBegin;
35277623264SMatthew G. Knepley   ierr = PetscFunctionListAdd(&PetscPartitionerList, sname, function);CHKERRQ(ierr);
35377623264SMatthew G. Knepley   PetscFunctionReturn(0);
35477623264SMatthew G. Knepley }
35577623264SMatthew G. Knepley 
35677623264SMatthew G. Knepley /*@C
35777623264SMatthew G. Knepley   PetscPartitionerSetType - Builds a particular PetscPartitioner
35877623264SMatthew G. Knepley 
35977623264SMatthew G. Knepley   Collective on PetscPartitioner
36077623264SMatthew G. Knepley 
36177623264SMatthew G. Knepley   Input Parameters:
36277623264SMatthew G. Knepley + part - The PetscPartitioner object
36377623264SMatthew G. Knepley - name - The kind of partitioner
36477623264SMatthew G. Knepley 
36577623264SMatthew G. Knepley   Options Database Key:
36677623264SMatthew G. Knepley . -petscpartitioner_type <type> - Sets the PetscPartitioner type; use -help for a list of available types
36777623264SMatthew G. Knepley 
36877623264SMatthew G. Knepley   Level: intermediate
36977623264SMatthew G. Knepley 
37077623264SMatthew G. Knepley .keywords: PetscPartitioner, set, type
37177623264SMatthew G. Knepley .seealso: PetscPartitionerGetType(), PetscPartitionerCreate()
37277623264SMatthew G. Knepley @*/
37377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetType(PetscPartitioner part, PetscPartitionerType name)
37477623264SMatthew G. Knepley {
37577623264SMatthew G. Knepley   PetscErrorCode (*r)(PetscPartitioner);
37677623264SMatthew G. Knepley   PetscBool      match;
37777623264SMatthew G. Knepley   PetscErrorCode ierr;
37877623264SMatthew G. Knepley 
37977623264SMatthew G. Knepley   PetscFunctionBegin;
38077623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
38177623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) part, name, &match);CHKERRQ(ierr);
38277623264SMatthew G. Knepley   if (match) PetscFunctionReturn(0);
38377623264SMatthew G. Knepley 
3840f51fdf8SToby Isaac   ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr);
38577623264SMatthew G. Knepley   ierr = PetscFunctionListFind(PetscPartitionerList, name, &r);CHKERRQ(ierr);
38677623264SMatthew G. Knepley   if (!r) SETERRQ1(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscPartitioner type: %s", name);
38777623264SMatthew G. Knepley 
38877623264SMatthew G. Knepley   if (part->ops->destroy) {
38977623264SMatthew G. Knepley     ierr              = (*part->ops->destroy)(part);CHKERRQ(ierr);
39077623264SMatthew G. Knepley     part->ops->destroy = NULL;
39177623264SMatthew G. Knepley   }
39277623264SMatthew G. Knepley   ierr = (*r)(part);CHKERRQ(ierr);
39377623264SMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) part, name);CHKERRQ(ierr);
39477623264SMatthew G. Knepley   PetscFunctionReturn(0);
39577623264SMatthew G. Knepley }
39677623264SMatthew G. Knepley 
39777623264SMatthew G. Knepley /*@C
39877623264SMatthew G. Knepley   PetscPartitionerGetType - Gets the PetscPartitioner type name (as a string) from the object.
39977623264SMatthew G. Knepley 
40077623264SMatthew G. Knepley   Not Collective
40177623264SMatthew G. Knepley 
40277623264SMatthew G. Knepley   Input Parameter:
40377623264SMatthew G. Knepley . part - The PetscPartitioner
40477623264SMatthew G. Knepley 
40577623264SMatthew G. Knepley   Output Parameter:
40677623264SMatthew G. Knepley . name - The PetscPartitioner type name
40777623264SMatthew G. Knepley 
40877623264SMatthew G. Knepley   Level: intermediate
40977623264SMatthew G. Knepley 
41077623264SMatthew G. Knepley .keywords: PetscPartitioner, get, type, name
41177623264SMatthew G. Knepley .seealso: PetscPartitionerSetType(), PetscPartitionerCreate()
41277623264SMatthew G. Knepley @*/
41377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerGetType(PetscPartitioner part, PetscPartitionerType *name)
41477623264SMatthew G. Knepley {
41577623264SMatthew G. Knepley   PetscErrorCode ierr;
41677623264SMatthew G. Knepley 
41777623264SMatthew G. Knepley   PetscFunctionBegin;
41877623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
41977623264SMatthew G. Knepley   PetscValidPointer(name, 2);
4200f51fdf8SToby Isaac   ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr);
42177623264SMatthew G. Knepley   *name = ((PetscObject) part)->type_name;
42277623264SMatthew G. Knepley   PetscFunctionReturn(0);
42377623264SMatthew G. Knepley }
42477623264SMatthew G. Knepley 
42577623264SMatthew G. Knepley /*@C
42677623264SMatthew G. Knepley   PetscPartitionerView - Views a PetscPartitioner
42777623264SMatthew G. Knepley 
42877623264SMatthew G. Knepley   Collective on PetscPartitioner
42977623264SMatthew G. Knepley 
43077623264SMatthew G. Knepley   Input Parameter:
43177623264SMatthew G. Knepley + part - the PetscPartitioner object to view
43277623264SMatthew G. Knepley - v    - the viewer
43377623264SMatthew G. Knepley 
43477623264SMatthew G. Knepley   Level: developer
43577623264SMatthew G. Knepley 
43677623264SMatthew G. Knepley .seealso: PetscPartitionerDestroy()
43777623264SMatthew G. Knepley @*/
43877623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView(PetscPartitioner part, PetscViewer v)
43977623264SMatthew G. Knepley {
44077623264SMatthew G. Knepley   PetscErrorCode ierr;
44177623264SMatthew G. Knepley 
44277623264SMatthew G. Knepley   PetscFunctionBegin;
44377623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
44477623264SMatthew G. Knepley   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) part), &v);CHKERRQ(ierr);}
44577623264SMatthew G. Knepley   if (part->ops->view) {ierr = (*part->ops->view)(part, v);CHKERRQ(ierr);}
44677623264SMatthew G. Knepley   PetscFunctionReturn(0);
44777623264SMatthew G. Knepley }
44877623264SMatthew G. Knepley 
44977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetTypeFromOptions_Internal(PetscPartitioner part)
45077623264SMatthew G. Knepley {
45177623264SMatthew G. Knepley   const char    *defaultType;
45277623264SMatthew G. Knepley   char           name[256];
45377623264SMatthew G. Knepley   PetscBool      flg;
45477623264SMatthew G. Knepley   PetscErrorCode ierr;
45577623264SMatthew G. Knepley 
45677623264SMatthew G. Knepley   PetscFunctionBegin;
45777623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
4587cf71cfdSLisandro Dalcin   if (!((PetscObject) part)->type_name)
4597cf71cfdSLisandro Dalcin #if defined(PETSC_HAVE_CHACO)
4607cf71cfdSLisandro Dalcin     defaultType = PETSCPARTITIONERCHACO;
4617cf71cfdSLisandro Dalcin #elif defined(PETSC_HAVE_PARMETIS)
4627cf71cfdSLisandro Dalcin     defaultType = PETSCPARTITIONERPARMETIS;
4637cf71cfdSLisandro Dalcin #else
4647cf71cfdSLisandro Dalcin     defaultType = PETSCPARTITIONERSIMPLE;
4657cf71cfdSLisandro Dalcin #endif
4667cf71cfdSLisandro Dalcin   else
4677cf71cfdSLisandro Dalcin     defaultType = ((PetscObject) part)->type_name;
4680f51fdf8SToby Isaac   ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr);
46977623264SMatthew G. Knepley 
47077623264SMatthew G. Knepley   ierr = PetscObjectOptionsBegin((PetscObject) part);CHKERRQ(ierr);
47177623264SMatthew G. Knepley   ierr = PetscOptionsFList("-petscpartitioner_type", "Graph partitioner", "PetscPartitionerSetType", PetscPartitionerList, defaultType, name, 256, &flg);CHKERRQ(ierr);
47277623264SMatthew G. Knepley   if (flg) {
47377623264SMatthew G. Knepley     ierr = PetscPartitionerSetType(part, name);CHKERRQ(ierr);
47477623264SMatthew G. Knepley   } else if (!((PetscObject) part)->type_name) {
47577623264SMatthew G. Knepley     ierr = PetscPartitionerSetType(part, defaultType);CHKERRQ(ierr);
47677623264SMatthew G. Knepley   }
47777623264SMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
47877623264SMatthew G. Knepley   PetscFunctionReturn(0);
47977623264SMatthew G. Knepley }
48077623264SMatthew G. Knepley 
48177623264SMatthew G. Knepley /*@
48277623264SMatthew G. Knepley   PetscPartitionerSetFromOptions - sets parameters in a PetscPartitioner from the options database
48377623264SMatthew G. Knepley 
48477623264SMatthew G. Knepley   Collective on PetscPartitioner
48577623264SMatthew G. Knepley 
48677623264SMatthew G. Knepley   Input Parameter:
48777623264SMatthew G. Knepley . part - the PetscPartitioner object to set options for
48877623264SMatthew G. Knepley 
48977623264SMatthew G. Knepley   Level: developer
49077623264SMatthew G. Knepley 
49177623264SMatthew G. Knepley .seealso: PetscPartitionerView()
49277623264SMatthew G. Knepley @*/
49377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetFromOptions(PetscPartitioner part)
49477623264SMatthew G. Knepley {
49577623264SMatthew G. Knepley   PetscErrorCode ierr;
49677623264SMatthew G. Knepley 
49777623264SMatthew G. Knepley   PetscFunctionBegin;
49877623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
49977623264SMatthew G. Knepley   ierr = PetscPartitionerSetTypeFromOptions_Internal(part);CHKERRQ(ierr);
50077623264SMatthew G. Knepley 
50177623264SMatthew G. Knepley   ierr = PetscObjectOptionsBegin((PetscObject) part);CHKERRQ(ierr);
502077101c0SMatthew G. Knepley   if (part->ops->setfromoptions) {ierr = (*part->ops->setfromoptions)(PetscOptionsObject,part);CHKERRQ(ierr);}
50377623264SMatthew G. Knepley   /* process any options handlers added with PetscObjectAddOptionsHandler() */
5040633abcbSJed Brown   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) part);CHKERRQ(ierr);
50577623264SMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
50677623264SMatthew G. Knepley   ierr = PetscPartitionerViewFromOptions(part, NULL, "-petscpartitioner_view");CHKERRQ(ierr);
50777623264SMatthew G. Knepley   PetscFunctionReturn(0);
50877623264SMatthew G. Knepley }
50977623264SMatthew G. Knepley 
51077623264SMatthew G. Knepley /*@C
51177623264SMatthew G. Knepley   PetscPartitionerSetUp - Construct data structures for the PetscPartitioner
51277623264SMatthew G. Knepley 
51377623264SMatthew G. Knepley   Collective on PetscPartitioner
51477623264SMatthew G. Knepley 
51577623264SMatthew G. Knepley   Input Parameter:
51677623264SMatthew G. Knepley . part - the PetscPartitioner object to setup
51777623264SMatthew G. Knepley 
51877623264SMatthew G. Knepley   Level: developer
51977623264SMatthew G. Knepley 
52077623264SMatthew G. Knepley .seealso: PetscPartitionerView(), PetscPartitionerDestroy()
52177623264SMatthew G. Knepley @*/
52277623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetUp(PetscPartitioner part)
52377623264SMatthew G. Knepley {
52477623264SMatthew G. Knepley   PetscErrorCode ierr;
52577623264SMatthew G. Knepley 
52677623264SMatthew G. Knepley   PetscFunctionBegin;
52777623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
52877623264SMatthew G. Knepley   if (part->ops->setup) {ierr = (*part->ops->setup)(part);CHKERRQ(ierr);}
52977623264SMatthew G. Knepley   PetscFunctionReturn(0);
53077623264SMatthew G. Knepley }
53177623264SMatthew G. Knepley 
53277623264SMatthew G. Knepley /*@
53377623264SMatthew G. Knepley   PetscPartitionerDestroy - Destroys a PetscPartitioner object
53477623264SMatthew G. Knepley 
53577623264SMatthew G. Knepley   Collective on PetscPartitioner
53677623264SMatthew G. Knepley 
53777623264SMatthew G. Knepley   Input Parameter:
53877623264SMatthew G. Knepley . part - the PetscPartitioner object to destroy
53977623264SMatthew G. Knepley 
54077623264SMatthew G. Knepley   Level: developer
54177623264SMatthew G. Knepley 
54277623264SMatthew G. Knepley .seealso: PetscPartitionerView()
54377623264SMatthew G. Knepley @*/
54477623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy(PetscPartitioner *part)
54577623264SMatthew G. Knepley {
54677623264SMatthew G. Knepley   PetscErrorCode ierr;
54777623264SMatthew G. Knepley 
54877623264SMatthew G. Knepley   PetscFunctionBegin;
54977623264SMatthew G. Knepley   if (!*part) PetscFunctionReturn(0);
55077623264SMatthew G. Knepley   PetscValidHeaderSpecific((*part), PETSCPARTITIONER_CLASSID, 1);
55177623264SMatthew G. Knepley 
55277623264SMatthew G. Knepley   if (--((PetscObject)(*part))->refct > 0) {*part = 0; PetscFunctionReturn(0);}
55377623264SMatthew G. Knepley   ((PetscObject) (*part))->refct = 0;
55477623264SMatthew G. Knepley 
55577623264SMatthew G. Knepley   if ((*part)->ops->destroy) {ierr = (*(*part)->ops->destroy)(*part);CHKERRQ(ierr);}
55677623264SMatthew G. Knepley   ierr = PetscHeaderDestroy(part);CHKERRQ(ierr);
55777623264SMatthew G. Knepley   PetscFunctionReturn(0);
55877623264SMatthew G. Knepley }
55977623264SMatthew G. Knepley 
56077623264SMatthew G. Knepley /*@
56177623264SMatthew G. Knepley   PetscPartitionerCreate - Creates an empty PetscPartitioner object. The type can then be set with PetscPartitionerSetType().
56277623264SMatthew G. Knepley 
56377623264SMatthew G. Knepley   Collective on MPI_Comm
56477623264SMatthew G. Knepley 
56577623264SMatthew G. Knepley   Input Parameter:
56677623264SMatthew G. Knepley . comm - The communicator for the PetscPartitioner object
56777623264SMatthew G. Knepley 
56877623264SMatthew G. Knepley   Output Parameter:
56977623264SMatthew G. Knepley . part - The PetscPartitioner object
57077623264SMatthew G. Knepley 
57177623264SMatthew G. Knepley   Level: beginner
57277623264SMatthew G. Knepley 
573dae52e14SToby Isaac .seealso: PetscPartitionerSetType(), PETSCPARTITIONERCHACO, PETSCPARTITIONERPARMETIS, PETSCPARTITIONERSHELL, PETSCPARTITIONERSIMPLE, PETSCPARTITIONERGATHER
57477623264SMatthew G. Knepley @*/
57577623264SMatthew G. Knepley PetscErrorCode PetscPartitionerCreate(MPI_Comm comm, PetscPartitioner *part)
57677623264SMatthew G. Knepley {
57777623264SMatthew G. Knepley   PetscPartitioner p;
57877623264SMatthew G. Knepley   PetscErrorCode   ierr;
57977623264SMatthew G. Knepley 
58077623264SMatthew G. Knepley   PetscFunctionBegin;
58177623264SMatthew G. Knepley   PetscValidPointer(part, 2);
58277623264SMatthew G. Knepley   *part = NULL;
58383cde681SMatthew G. Knepley   ierr = DMInitializePackage();CHKERRQ(ierr);
58477623264SMatthew G. Knepley 
58573107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(p, PETSCPARTITIONER_CLASSID, "PetscPartitioner", "Graph Partitioner", "PetscPartitioner", comm, PetscPartitionerDestroy, PetscPartitionerView);CHKERRQ(ierr);
58677623264SMatthew G. Knepley 
58777623264SMatthew G. Knepley   *part = p;
58877623264SMatthew G. Knepley   PetscFunctionReturn(0);
58977623264SMatthew G. Knepley }
59077623264SMatthew G. Knepley 
59177623264SMatthew G. Knepley /*@
59277623264SMatthew G. Knepley   PetscPartitionerPartition - Create a non-overlapping partition of the cells in the mesh
59377623264SMatthew G. Knepley 
59477623264SMatthew G. Knepley   Collective on DM
59577623264SMatthew G. Knepley 
59677623264SMatthew G. Knepley   Input Parameters:
59777623264SMatthew G. Knepley + part    - The PetscPartitioner
598f8987ae8SMichael Lange - dm      - The mesh DM
59977623264SMatthew G. Knepley 
60077623264SMatthew G. Knepley   Output Parameters:
60177623264SMatthew G. Knepley + partSection     - The PetscSection giving the division of points by partition
602f8987ae8SMichael Lange - partition       - The list of points by partition
60377623264SMatthew G. Knepley 
60477623264SMatthew G. Knepley   Note: Instead of cells, points at a given height can be partitioned by calling PetscPartitionerSetPointHeight()
60577623264SMatthew G. Knepley 
60677623264SMatthew G. Knepley   Level: developer
60777623264SMatthew G. Knepley 
60877623264SMatthew G. Knepley .seealso DMPlexDistribute(), PetscPartitionerSetPointHeight(), PetscPartitionerCreate()
6094b15ede2SMatthew G. Knepley @*/
610f8987ae8SMichael Lange PetscErrorCode PetscPartitionerPartition(PetscPartitioner part, DM dm, PetscSection partSection, IS *partition)
61177623264SMatthew G. Knepley {
61277623264SMatthew G. Knepley   PetscMPIInt    size;
61377623264SMatthew G. Knepley   PetscErrorCode ierr;
61477623264SMatthew G. Knepley 
61577623264SMatthew G. Knepley   PetscFunctionBegin;
61677623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
61777623264SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
61877623264SMatthew G. Knepley   PetscValidHeaderSpecific(partSection, PETSC_SECTION_CLASSID, 4);
61977623264SMatthew G. Knepley   PetscValidPointer(partition, 5);
62077623264SMatthew G. Knepley   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) part), &size);CHKERRQ(ierr);
62177623264SMatthew G. Knepley   if (size == 1) {
62277623264SMatthew G. Knepley     PetscInt *points;
62377623264SMatthew G. Knepley     PetscInt  cStart, cEnd, c;
62477623264SMatthew G. Knepley 
62577623264SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, part->height, &cStart, &cEnd);CHKERRQ(ierr);
62677623264SMatthew G. Knepley     ierr = PetscSectionSetChart(partSection, 0, size);CHKERRQ(ierr);
62777623264SMatthew G. Knepley     ierr = PetscSectionSetDof(partSection, 0, cEnd-cStart);CHKERRQ(ierr);
62877623264SMatthew G. Knepley     ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
62977623264SMatthew G. Knepley     ierr = PetscMalloc1(cEnd-cStart, &points);CHKERRQ(ierr);
63077623264SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) points[c] = c;
63177623264SMatthew G. Knepley     ierr = ISCreateGeneral(PetscObjectComm((PetscObject) part), cEnd-cStart, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
63277623264SMatthew G. Knepley     PetscFunctionReturn(0);
63377623264SMatthew G. Knepley   }
63477623264SMatthew G. Knepley   if (part->height == 0) {
63577623264SMatthew G. Knepley     PetscInt numVertices;
63677623264SMatthew G. Knepley     PetscInt *start     = NULL;
63777623264SMatthew G. Knepley     PetscInt *adjacency = NULL;
6383fa7752dSToby Isaac     IS       globalNumbering;
63977623264SMatthew G. Knepley 
6403fa7752dSToby Isaac     ierr = DMPlexCreatePartitionerGraph(dm, 0, &numVertices, &start, &adjacency, &globalNumbering);CHKERRQ(ierr);
64177623264SMatthew G. Knepley     if (!part->ops->partition) SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_WRONGSTATE, "PetscPartitioner has no type");
64277623264SMatthew G. Knepley     ierr = (*part->ops->partition)(part, dm, size, numVertices, start, adjacency, partSection, partition);CHKERRQ(ierr);
64377623264SMatthew G. Knepley     ierr = PetscFree(start);CHKERRQ(ierr);
64477623264SMatthew G. Knepley     ierr = PetscFree(adjacency);CHKERRQ(ierr);
6453fa7752dSToby Isaac     if (globalNumbering) { /* partition is wrt global unique numbering: change this to be wrt local numbering */
6463fa7752dSToby Isaac       const PetscInt *globalNum;
6473fa7752dSToby Isaac       const PetscInt *partIdx;
6483fa7752dSToby Isaac       PetscInt *map, cStart, cEnd;
6493fa7752dSToby Isaac       PetscInt *adjusted, i, localSize, offset;
6503fa7752dSToby Isaac       IS    newPartition;
6513fa7752dSToby Isaac 
6523fa7752dSToby Isaac       ierr = ISGetLocalSize(*partition,&localSize);CHKERRQ(ierr);
6533fa7752dSToby Isaac       ierr = PetscMalloc1(localSize,&adjusted);CHKERRQ(ierr);
6543fa7752dSToby Isaac       ierr = ISGetIndices(globalNumbering,&globalNum);CHKERRQ(ierr);
6553fa7752dSToby Isaac       ierr = ISGetIndices(*partition,&partIdx);CHKERRQ(ierr);
6563fa7752dSToby Isaac       ierr = PetscMalloc1(localSize,&map);CHKERRQ(ierr);
6573fa7752dSToby Isaac       ierr = DMPlexGetHeightStratum(dm, part->height, &cStart, &cEnd);CHKERRQ(ierr);
6583fa7752dSToby Isaac       for (i = cStart, offset = 0; i < cEnd; i++) {
6593fa7752dSToby Isaac         if (globalNum[i - cStart] >= 0) map[offset++] = i;
6603fa7752dSToby Isaac       }
6613fa7752dSToby Isaac       for (i = 0; i < localSize; i++) {
6623fa7752dSToby Isaac         adjusted[i] = map[partIdx[i]];
6633fa7752dSToby Isaac       }
6643fa7752dSToby Isaac       ierr = PetscFree(map);CHKERRQ(ierr);
6653fa7752dSToby Isaac       ierr = ISRestoreIndices(*partition,&partIdx);CHKERRQ(ierr);
6663fa7752dSToby Isaac       ierr = ISRestoreIndices(globalNumbering,&globalNum);CHKERRQ(ierr);
6673fa7752dSToby Isaac       ierr = ISCreateGeneral(PETSC_COMM_SELF,localSize,adjusted,PETSC_OWN_POINTER,&newPartition);CHKERRQ(ierr);
6683fa7752dSToby Isaac       ierr = ISDestroy(&globalNumbering);CHKERRQ(ierr);
6693fa7752dSToby Isaac       ierr = ISDestroy(partition);CHKERRQ(ierr);
6703fa7752dSToby Isaac       *partition = newPartition;
6713fa7752dSToby Isaac     }
67277623264SMatthew G. Knepley   } else SETERRQ1(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_OUTOFRANGE, "Invalid height %D for points to partition", part->height);
67377623264SMatthew G. Knepley   PetscFunctionReturn(0);
6743fa7752dSToby Isaac 
67577623264SMatthew G. Knepley }
67677623264SMatthew G. Knepley 
67777623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy_Shell(PetscPartitioner part)
67877623264SMatthew G. Knepley {
67977623264SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
68077623264SMatthew G. Knepley   PetscErrorCode          ierr;
68177623264SMatthew G. Knepley 
68277623264SMatthew G. Knepley   PetscFunctionBegin;
68377623264SMatthew G. Knepley   ierr = PetscSectionDestroy(&p->section);CHKERRQ(ierr);
68477623264SMatthew G. Knepley   ierr = ISDestroy(&p->partition);CHKERRQ(ierr);
68577623264SMatthew G. Knepley   ierr = PetscFree(p);CHKERRQ(ierr);
68677623264SMatthew G. Knepley   PetscFunctionReturn(0);
68777623264SMatthew G. Knepley }
68877623264SMatthew G. Knepley 
68977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Shell_Ascii(PetscPartitioner part, PetscViewer viewer)
69077623264SMatthew G. Knepley {
691077101c0SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
69277623264SMatthew G. Knepley   PetscViewerFormat       format;
69377623264SMatthew G. Knepley   PetscErrorCode          ierr;
69477623264SMatthew G. Knepley 
69577623264SMatthew G. Knepley   PetscFunctionBegin;
69677623264SMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
69777623264SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Shell Graph Partitioner:\n");CHKERRQ(ierr);
698077101c0SMatthew G. Knepley   if (p->random) {
699077101c0SMatthew G. Knepley     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
700077101c0SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "using random partition\n");CHKERRQ(ierr);
701077101c0SMatthew G. Knepley     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
702077101c0SMatthew G. Knepley   }
70377623264SMatthew G. Knepley   PetscFunctionReturn(0);
70477623264SMatthew G. Knepley }
70577623264SMatthew G. Knepley 
70677623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Shell(PetscPartitioner part, PetscViewer viewer)
70777623264SMatthew G. Knepley {
70877623264SMatthew G. Knepley   PetscBool      iascii;
70977623264SMatthew G. Knepley   PetscErrorCode ierr;
71077623264SMatthew G. Knepley 
71177623264SMatthew G. Knepley   PetscFunctionBegin;
71277623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
71377623264SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
71477623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
71577623264SMatthew G. Knepley   if (iascii) {ierr = PetscPartitionerView_Shell_Ascii(part, viewer);CHKERRQ(ierr);}
71677623264SMatthew G. Knepley   PetscFunctionReturn(0);
71777623264SMatthew G. Knepley }
71877623264SMatthew G. Knepley 
719077101c0SMatthew G. Knepley PetscErrorCode PetscPartitionerSetFromOptions_Shell(PetscOptionItems *PetscOptionsObject, PetscPartitioner part)
720077101c0SMatthew G. Knepley {
721077101c0SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
722077101c0SMatthew G. Knepley   PetscErrorCode          ierr;
723077101c0SMatthew G. Knepley 
724077101c0SMatthew G. Knepley   PetscFunctionBegin;
725077101c0SMatthew G. Knepley   ierr = PetscOptionsHead(PetscOptionsObject, "PetscPartitionerShell Options");CHKERRQ(ierr);
726077101c0SMatthew G. Knepley   ierr = PetscOptionsBool("-petscpartitioner_shell_random", "Use a random partition", "PetscPartitionerView", PETSC_FALSE, &p->random, NULL);CHKERRQ(ierr);
727077101c0SMatthew G. Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
728077101c0SMatthew G. Knepley   PetscFunctionReturn(0);
729077101c0SMatthew G. Knepley }
730077101c0SMatthew G. Knepley 
73177623264SMatthew G. Knepley PetscErrorCode PetscPartitionerPartition_Shell(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
73277623264SMatthew G. Knepley {
73377623264SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
73477623264SMatthew G. Knepley   PetscInt                np;
73577623264SMatthew G. Knepley   PetscErrorCode          ierr;
73677623264SMatthew G. Knepley 
73777623264SMatthew G. Knepley   PetscFunctionBegin;
738077101c0SMatthew G. Knepley   if (p->random) {
739077101c0SMatthew G. Knepley     PetscRandom r;
740077101c0SMatthew G. Knepley     PetscInt   *sizes, *points, v;
741077101c0SMatthew G. Knepley 
742077101c0SMatthew G. Knepley     ierr = PetscRandomCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
743351ac257SMatthew G. Knepley     ierr = PetscRandomSetInterval(r, 0.0, (PetscScalar) (nparts+1));CHKERRQ(ierr);
744077101c0SMatthew G. Knepley     ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr);
745077101c0SMatthew G. Knepley     ierr = PetscCalloc2(nparts, &sizes, numVertices, &points);CHKERRQ(ierr);
746077101c0SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
747077101c0SMatthew G. Knepley       PetscReal val;
748077101c0SMatthew G. Knepley       PetscInt  part;
749077101c0SMatthew G. Knepley 
750077101c0SMatthew G. Knepley       ierr = PetscRandomGetValueReal(r, &val);CHKERRQ(ierr);
751077101c0SMatthew G. Knepley       part = PetscFloorReal(val);
752077101c0SMatthew G. Knepley       points[v] = part;
753077101c0SMatthew G. Knepley       ++sizes[part];
754077101c0SMatthew G. Knepley     }
755077101c0SMatthew G. Knepley     ierr = PetscRandomDestroy(&r);CHKERRQ(ierr);
756077101c0SMatthew G. Knepley     ierr = PetscPartitionerShellSetPartition(part, nparts, sizes, points);CHKERRQ(ierr);
757077101c0SMatthew G. Knepley     ierr = PetscFree2(sizes, points);CHKERRQ(ierr);
758077101c0SMatthew G. Knepley   }
75977623264SMatthew G. Knepley   ierr = PetscSectionGetChart(p->section, NULL, &np);CHKERRQ(ierr);
76077623264SMatthew G. Knepley   if (nparts != np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of requested partitions %d != configured partitions %d", nparts, np);
76177623264SMatthew G. Knepley   ierr = ISGetLocalSize(p->partition, &np);CHKERRQ(ierr);
76277623264SMatthew G. Knepley   if (numVertices != np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of input vertices %d != configured vertices %d", numVertices, np);
7635680f57bSMatthew G. Knepley   ierr = PetscSectionCopy(p->section, partSection);CHKERRQ(ierr);
76477623264SMatthew G. Knepley   *partition = p->partition;
76577623264SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) p->partition);CHKERRQ(ierr);
76677623264SMatthew G. Knepley   PetscFunctionReturn(0);
76777623264SMatthew G. Knepley }
76877623264SMatthew G. Knepley 
76977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerInitialize_Shell(PetscPartitioner part)
77077623264SMatthew G. Knepley {
77177623264SMatthew G. Knepley   PetscFunctionBegin;
77277623264SMatthew G. Knepley   part->ops->view           = PetscPartitionerView_Shell;
773077101c0SMatthew G. Knepley   part->ops->setfromoptions = PetscPartitionerSetFromOptions_Shell;
77477623264SMatthew G. Knepley   part->ops->destroy        = PetscPartitionerDestroy_Shell;
77577623264SMatthew G. Knepley   part->ops->partition      = PetscPartitionerPartition_Shell;
77677623264SMatthew G. Knepley   PetscFunctionReturn(0);
77777623264SMatthew G. Knepley }
77877623264SMatthew G. Knepley 
77977623264SMatthew G. Knepley /*MC
78077623264SMatthew G. Knepley   PETSCPARTITIONERSHELL = "shell" - A PetscPartitioner object
78177623264SMatthew G. Knepley 
78277623264SMatthew G. Knepley   Level: intermediate
78377623264SMatthew G. Knepley 
78477623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
78577623264SMatthew G. Knepley M*/
78677623264SMatthew G. Knepley 
78777623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Shell(PetscPartitioner part)
78877623264SMatthew G. Knepley {
78977623264SMatthew G. Knepley   PetscPartitioner_Shell *p;
79077623264SMatthew G. Knepley   PetscErrorCode          ierr;
79177623264SMatthew G. Knepley 
79277623264SMatthew G. Knepley   PetscFunctionBegin;
79377623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
79477623264SMatthew G. Knepley   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
79577623264SMatthew G. Knepley   part->data = p;
79677623264SMatthew G. Knepley 
79777623264SMatthew G. Knepley   ierr = PetscPartitionerInitialize_Shell(part);CHKERRQ(ierr);
798077101c0SMatthew G. Knepley   p->random = PETSC_FALSE;
79977623264SMatthew G. Knepley   PetscFunctionReturn(0);
80077623264SMatthew G. Knepley }
80177623264SMatthew G. Knepley 
8025680f57bSMatthew G. Knepley /*@C
8035680f57bSMatthew G. Knepley   PetscPartitionerShellSetPartition - Set an artifical partition for a mesh
8045680f57bSMatthew G. Knepley 
805077101c0SMatthew G. Knepley   Collective on Part
8065680f57bSMatthew G. Knepley 
8075680f57bSMatthew G. Knepley   Input Parameters:
8085680f57bSMatthew G. Knepley + part     - The PetscPartitioner
809*9852e123SBarry Smith . size - The number of partitions
810*9852e123SBarry Smith . sizes    - array of size size (or NULL) providing the number of points in each partition
8119758bf69SToby Isaac - points   - array of size sum(sizes) (may be NULL iff sizes is NULL), a permutation of the points that groups those assigned to each partition in order (i.e., partition 0 first, partition 1 next, etc.)
8125680f57bSMatthew G. Knepley 
8135680f57bSMatthew G. Knepley   Level: developer
8145680f57bSMatthew G. Knepley 
815b7e49471SLawrence Mitchell   Notes:
816b7e49471SLawrence Mitchell 
817b7e49471SLawrence Mitchell     It is safe to free the sizes and points arrays after use in this routine.
818b7e49471SLawrence Mitchell 
8195680f57bSMatthew G. Knepley .seealso DMPlexDistribute(), PetscPartitionerCreate()
8205680f57bSMatthew G. Knepley @*/
821*9852e123SBarry Smith PetscErrorCode PetscPartitionerShellSetPartition(PetscPartitioner part, PetscInt size, const PetscInt sizes[], const PetscInt points[])
8225680f57bSMatthew G. Knepley {
8235680f57bSMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
8245680f57bSMatthew G. Knepley   PetscInt                proc, numPoints;
8255680f57bSMatthew G. Knepley   PetscErrorCode          ierr;
8265680f57bSMatthew G. Knepley 
8275680f57bSMatthew G. Knepley   PetscFunctionBegin;
8285680f57bSMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
8295680f57bSMatthew G. Knepley   if (sizes) {PetscValidPointer(sizes, 3);}
8305680f57bSMatthew G. Knepley   if (sizes) {PetscValidPointer(points, 4);}
8315680f57bSMatthew G. Knepley   ierr = PetscSectionDestroy(&p->section);CHKERRQ(ierr);
8325680f57bSMatthew G. Knepley   ierr = ISDestroy(&p->partition);CHKERRQ(ierr);
8335680f57bSMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) part), &p->section);CHKERRQ(ierr);
834*9852e123SBarry Smith   ierr = PetscSectionSetChart(p->section, 0, size);CHKERRQ(ierr);
8355680f57bSMatthew G. Knepley   if (sizes) {
836*9852e123SBarry Smith     for (proc = 0; proc < size; ++proc) {
8375680f57bSMatthew G. Knepley       ierr = PetscSectionSetDof(p->section, proc, sizes[proc]);CHKERRQ(ierr);
8385680f57bSMatthew G. Knepley     }
8395680f57bSMatthew G. Knepley   }
8405680f57bSMatthew G. Knepley   ierr = PetscSectionSetUp(p->section);CHKERRQ(ierr);
8415680f57bSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(p->section, &numPoints);CHKERRQ(ierr);
8425680f57bSMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) part), numPoints, points, PETSC_COPY_VALUES, &p->partition);CHKERRQ(ierr);
8435680f57bSMatthew G. Knepley   PetscFunctionReturn(0);
8445680f57bSMatthew G. Knepley }
8455680f57bSMatthew G. Knepley 
846077101c0SMatthew G. Knepley /*@
847077101c0SMatthew G. Knepley   PetscPartitionerShellSetRandom - Set the flag to use a random partition
848077101c0SMatthew G. Knepley 
849077101c0SMatthew G. Knepley   Collective on Part
850077101c0SMatthew G. Knepley 
851077101c0SMatthew G. Knepley   Input Parameters:
852077101c0SMatthew G. Knepley + part   - The PetscPartitioner
853077101c0SMatthew G. Knepley - random - The flag to use a random partition
854077101c0SMatthew G. Knepley 
855077101c0SMatthew G. Knepley   Level: intermediate
856077101c0SMatthew G. Knepley 
857077101c0SMatthew G. Knepley .seealso PetscPartitionerShellGetRandom(), PetscPartitionerCreate()
858077101c0SMatthew G. Knepley @*/
859077101c0SMatthew G. Knepley PetscErrorCode PetscPartitionerShellSetRandom(PetscPartitioner part, PetscBool random)
860077101c0SMatthew G. Knepley {
861077101c0SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
862077101c0SMatthew G. Knepley 
863077101c0SMatthew G. Knepley   PetscFunctionBegin;
864077101c0SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
865077101c0SMatthew G. Knepley   p->random = random;
866077101c0SMatthew G. Knepley   PetscFunctionReturn(0);
867077101c0SMatthew G. Knepley }
868077101c0SMatthew G. Knepley 
869077101c0SMatthew G. Knepley /*@
870077101c0SMatthew G. Knepley   PetscPartitionerShellGetRandom - get the flag to use a random partition
871077101c0SMatthew G. Knepley 
872077101c0SMatthew G. Knepley   Collective on Part
873077101c0SMatthew G. Knepley 
874077101c0SMatthew G. Knepley   Input Parameter:
875077101c0SMatthew G. Knepley . part   - The PetscPartitioner
876077101c0SMatthew G. Knepley 
877077101c0SMatthew G. Knepley   Output Parameter
878077101c0SMatthew G. Knepley . random - The flag to use a random partition
879077101c0SMatthew G. Knepley 
880077101c0SMatthew G. Knepley   Level: intermediate
881077101c0SMatthew G. Knepley 
882077101c0SMatthew G. Knepley .seealso PetscPartitionerShellSetRandom(), PetscPartitionerCreate()
883077101c0SMatthew G. Knepley @*/
884077101c0SMatthew G. Knepley PetscErrorCode PetscPartitionerShellGetRandom(PetscPartitioner part, PetscBool *random)
885077101c0SMatthew G. Knepley {
886077101c0SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
887077101c0SMatthew G. Knepley 
888077101c0SMatthew G. Knepley   PetscFunctionBegin;
889077101c0SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
890077101c0SMatthew G. Knepley   PetscValidPointer(random, 2);
891077101c0SMatthew G. Knepley   *random = p->random;
892077101c0SMatthew G. Knepley   PetscFunctionReturn(0);
893077101c0SMatthew G. Knepley }
894077101c0SMatthew G. Knepley 
895555a9cf8SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy_Simple(PetscPartitioner part)
896555a9cf8SMatthew G. Knepley {
897555a9cf8SMatthew G. Knepley   PetscPartitioner_Simple *p = (PetscPartitioner_Simple *) part->data;
898555a9cf8SMatthew G. Knepley   PetscErrorCode          ierr;
899555a9cf8SMatthew G. Knepley 
900555a9cf8SMatthew G. Knepley   PetscFunctionBegin;
901555a9cf8SMatthew G. Knepley   ierr = PetscFree(p);CHKERRQ(ierr);
902555a9cf8SMatthew G. Knepley   PetscFunctionReturn(0);
903555a9cf8SMatthew G. Knepley }
904555a9cf8SMatthew G. Knepley 
905555a9cf8SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Simple_Ascii(PetscPartitioner part, PetscViewer viewer)
906555a9cf8SMatthew G. Knepley {
907555a9cf8SMatthew G. Knepley   PetscViewerFormat format;
908555a9cf8SMatthew G. Knepley   PetscErrorCode    ierr;
909555a9cf8SMatthew G. Knepley 
910555a9cf8SMatthew G. Knepley   PetscFunctionBegin;
911555a9cf8SMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
912555a9cf8SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Simple Graph Partitioner:\n");CHKERRQ(ierr);
913555a9cf8SMatthew G. Knepley   PetscFunctionReturn(0);
914555a9cf8SMatthew G. Knepley }
915555a9cf8SMatthew G. Knepley 
916555a9cf8SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Simple(PetscPartitioner part, PetscViewer viewer)
917555a9cf8SMatthew G. Knepley {
918555a9cf8SMatthew G. Knepley   PetscBool      iascii;
919555a9cf8SMatthew G. Knepley   PetscErrorCode ierr;
920555a9cf8SMatthew G. Knepley 
921555a9cf8SMatthew G. Knepley   PetscFunctionBegin;
922555a9cf8SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
923555a9cf8SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
924555a9cf8SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
925555a9cf8SMatthew G. Knepley   if (iascii) {ierr = PetscPartitionerView_Simple_Ascii(part, viewer);CHKERRQ(ierr);}
926555a9cf8SMatthew G. Knepley   PetscFunctionReturn(0);
927555a9cf8SMatthew G. Knepley }
928555a9cf8SMatthew G. Knepley 
929555a9cf8SMatthew G. Knepley PetscErrorCode PetscPartitionerPartition_Simple(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
930555a9cf8SMatthew G. Knepley {
931cead94edSToby Isaac   MPI_Comm       comm;
932555a9cf8SMatthew G. Knepley   PetscInt       np;
933cead94edSToby Isaac   PetscMPIInt    size;
934555a9cf8SMatthew G. Knepley   PetscErrorCode ierr;
935555a9cf8SMatthew G. Knepley 
936555a9cf8SMatthew G. Knepley   PetscFunctionBegin;
937cead94edSToby Isaac   comm = PetscObjectComm((PetscObject)dm);
938cead94edSToby Isaac   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
939555a9cf8SMatthew G. Knepley   ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
940555a9cf8SMatthew G. Knepley   ierr = ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition);CHKERRQ(ierr);
941cead94edSToby Isaac   if (size == 1) {
942cead94edSToby Isaac     for (np = 0; np < nparts; ++np) {ierr = PetscSectionSetDof(partSection, np, numVertices/nparts + ((numVertices % nparts) > np));CHKERRQ(ierr);}
943cead94edSToby Isaac   }
944cead94edSToby Isaac   else {
945cead94edSToby Isaac     PetscMPIInt rank;
946cead94edSToby Isaac     PetscInt nvGlobal, *offsets, myFirst, myLast;
947cead94edSToby Isaac 
948a679a563SToby Isaac     ierr = PetscMalloc1(size+1,&offsets);CHKERRQ(ierr);
949cead94edSToby Isaac     offsets[0] = 0;
950cead94edSToby Isaac     ierr = MPI_Allgather(&numVertices,1,MPIU_INT,&offsets[1],1,MPIU_INT,comm);CHKERRQ(ierr);
951cead94edSToby Isaac     for (np = 2; np <= size; np++) {
952cead94edSToby Isaac       offsets[np] += offsets[np-1];
953cead94edSToby Isaac     }
954cead94edSToby Isaac     nvGlobal = offsets[size];
955cead94edSToby Isaac     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
956cead94edSToby Isaac     myFirst = offsets[rank];
957cead94edSToby Isaac     myLast  = offsets[rank + 1] - 1;
958cead94edSToby Isaac     ierr = PetscFree(offsets);CHKERRQ(ierr);
959cead94edSToby Isaac     if (numVertices) {
960cead94edSToby Isaac       PetscInt firstPart = 0, firstLargePart = 0;
961cead94edSToby Isaac       PetscInt lastPart = 0, lastLargePart = 0;
962cead94edSToby Isaac       PetscInt rem = nvGlobal % nparts;
963cead94edSToby Isaac       PetscInt pSmall = nvGlobal/nparts;
964cead94edSToby Isaac       PetscInt pBig = nvGlobal/nparts + 1;
965cead94edSToby Isaac 
966cead94edSToby Isaac 
967cead94edSToby Isaac       if (rem) {
968cead94edSToby Isaac         firstLargePart = myFirst / pBig;
969cead94edSToby Isaac         lastLargePart  = myLast  / pBig;
970cead94edSToby Isaac 
971cead94edSToby Isaac         if (firstLargePart < rem) {
972cead94edSToby Isaac           firstPart = firstLargePart;
973cead94edSToby Isaac         }
974cead94edSToby Isaac         else {
975cead94edSToby Isaac           firstPart = rem + (myFirst - (rem * pBig)) / pSmall;
976cead94edSToby Isaac         }
977cead94edSToby Isaac         if (lastLargePart < rem) {
978cead94edSToby Isaac           lastPart = lastLargePart;
979cead94edSToby Isaac         }
980cead94edSToby Isaac         else {
981cead94edSToby Isaac           lastPart = rem + (myLast - (rem * pBig)) / pSmall;
982cead94edSToby Isaac         }
983cead94edSToby Isaac       }
984cead94edSToby Isaac       else {
985cead94edSToby Isaac         firstPart = myFirst / (nvGlobal/nparts);
986cead94edSToby Isaac         lastPart  = myLast  / (nvGlobal/nparts);
987cead94edSToby Isaac       }
988cead94edSToby Isaac 
989cead94edSToby Isaac       for (np = firstPart; np <= lastPart; np++) {
990cead94edSToby Isaac         PetscInt PartStart =  np    * (nvGlobal/nparts) + PetscMin(nvGlobal % nparts,np);
991cead94edSToby Isaac         PetscInt PartEnd   = (np+1) * (nvGlobal/nparts) + PetscMin(nvGlobal % nparts,np+1);
992cead94edSToby Isaac 
993cead94edSToby Isaac         PartStart = PetscMax(PartStart,myFirst);
994cead94edSToby Isaac         PartEnd   = PetscMin(PartEnd,myLast+1);
995cead94edSToby Isaac         ierr = PetscSectionSetDof(partSection,np,PartEnd-PartStart);CHKERRQ(ierr);
996cead94edSToby Isaac       }
997cead94edSToby Isaac     }
998cead94edSToby Isaac   }
999cead94edSToby Isaac   ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
1000555a9cf8SMatthew G. Knepley   PetscFunctionReturn(0);
1001555a9cf8SMatthew G. Knepley }
1002555a9cf8SMatthew G. Knepley 
1003555a9cf8SMatthew G. Knepley PetscErrorCode PetscPartitionerInitialize_Simple(PetscPartitioner part)
1004555a9cf8SMatthew G. Knepley {
1005555a9cf8SMatthew G. Knepley   PetscFunctionBegin;
1006555a9cf8SMatthew G. Knepley   part->ops->view      = PetscPartitionerView_Simple;
1007555a9cf8SMatthew G. Knepley   part->ops->destroy   = PetscPartitionerDestroy_Simple;
1008555a9cf8SMatthew G. Knepley   part->ops->partition = PetscPartitionerPartition_Simple;
1009555a9cf8SMatthew G. Knepley   PetscFunctionReturn(0);
1010555a9cf8SMatthew G. Knepley }
1011555a9cf8SMatthew G. Knepley 
1012555a9cf8SMatthew G. Knepley /*MC
1013555a9cf8SMatthew G. Knepley   PETSCPARTITIONERSIMPLE = "simple" - A PetscPartitioner object
1014555a9cf8SMatthew G. Knepley 
1015555a9cf8SMatthew G. Knepley   Level: intermediate
1016555a9cf8SMatthew G. Knepley 
1017555a9cf8SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
1018555a9cf8SMatthew G. Knepley M*/
1019555a9cf8SMatthew G. Knepley 
1020555a9cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Simple(PetscPartitioner part)
1021555a9cf8SMatthew G. Knepley {
1022555a9cf8SMatthew G. Knepley   PetscPartitioner_Simple *p;
1023555a9cf8SMatthew G. Knepley   PetscErrorCode           ierr;
1024555a9cf8SMatthew G. Knepley 
1025555a9cf8SMatthew G. Knepley   PetscFunctionBegin;
1026555a9cf8SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
1027555a9cf8SMatthew G. Knepley   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
1028555a9cf8SMatthew G. Knepley   part->data = p;
1029555a9cf8SMatthew G. Knepley 
1030555a9cf8SMatthew G. Knepley   ierr = PetscPartitionerInitialize_Simple(part);CHKERRQ(ierr);
1031555a9cf8SMatthew G. Knepley   PetscFunctionReturn(0);
1032555a9cf8SMatthew G. Knepley }
1033555a9cf8SMatthew G. Knepley 
1034dae52e14SToby Isaac PetscErrorCode PetscPartitionerDestroy_Gather(PetscPartitioner part)
1035dae52e14SToby Isaac {
1036dae52e14SToby Isaac   PetscPartitioner_Gather *p = (PetscPartitioner_Gather *) part->data;
1037dae52e14SToby Isaac   PetscErrorCode          ierr;
1038dae52e14SToby Isaac 
1039dae52e14SToby Isaac   PetscFunctionBegin;
1040dae52e14SToby Isaac   ierr = PetscFree(p);CHKERRQ(ierr);
1041dae52e14SToby Isaac   PetscFunctionReturn(0);
1042dae52e14SToby Isaac }
1043dae52e14SToby Isaac 
1044dae52e14SToby Isaac PetscErrorCode PetscPartitionerView_Gather_Ascii(PetscPartitioner part, PetscViewer viewer)
1045dae52e14SToby Isaac {
1046dae52e14SToby Isaac   PetscViewerFormat format;
1047dae52e14SToby Isaac   PetscErrorCode    ierr;
1048dae52e14SToby Isaac 
1049dae52e14SToby Isaac   PetscFunctionBegin;
1050dae52e14SToby Isaac   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
1051dae52e14SToby Isaac   ierr = PetscViewerASCIIPrintf(viewer, "Gather Graph Partitioner:\n");CHKERRQ(ierr);
1052dae52e14SToby Isaac   PetscFunctionReturn(0);
1053dae52e14SToby Isaac }
1054dae52e14SToby Isaac 
1055dae52e14SToby Isaac PetscErrorCode PetscPartitionerView_Gather(PetscPartitioner part, PetscViewer viewer)
1056dae52e14SToby Isaac {
1057dae52e14SToby Isaac   PetscBool      iascii;
1058dae52e14SToby Isaac   PetscErrorCode ierr;
1059dae52e14SToby Isaac 
1060dae52e14SToby Isaac   PetscFunctionBegin;
1061dae52e14SToby Isaac   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
1062dae52e14SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1063dae52e14SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
1064dae52e14SToby Isaac   if (iascii) {ierr = PetscPartitionerView_Gather_Ascii(part, viewer);CHKERRQ(ierr);}
1065dae52e14SToby Isaac   PetscFunctionReturn(0);
1066dae52e14SToby Isaac }
1067dae52e14SToby Isaac 
1068dae52e14SToby Isaac PetscErrorCode PetscPartitionerPartition_Gather(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
1069dae52e14SToby Isaac {
1070dae52e14SToby Isaac   PetscInt       np;
1071dae52e14SToby Isaac   PetscErrorCode ierr;
1072dae52e14SToby Isaac 
1073dae52e14SToby Isaac   PetscFunctionBegin;
1074dae52e14SToby Isaac   ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
1075dae52e14SToby Isaac   ierr = ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition);CHKERRQ(ierr);
1076dae52e14SToby Isaac   ierr = PetscSectionSetDof(partSection,0,numVertices);CHKERRQ(ierr);
1077dae52e14SToby Isaac   for (np = 1; np < nparts; ++np) {ierr = PetscSectionSetDof(partSection, np, 0);CHKERRQ(ierr);}
1078dae52e14SToby Isaac   ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
1079dae52e14SToby Isaac   PetscFunctionReturn(0);
1080dae52e14SToby Isaac }
1081dae52e14SToby Isaac 
1082dae52e14SToby Isaac PetscErrorCode PetscPartitionerInitialize_Gather(PetscPartitioner part)
1083dae52e14SToby Isaac {
1084dae52e14SToby Isaac   PetscFunctionBegin;
1085dae52e14SToby Isaac   part->ops->view      = PetscPartitionerView_Gather;
1086dae52e14SToby Isaac   part->ops->destroy   = PetscPartitionerDestroy_Gather;
1087dae52e14SToby Isaac   part->ops->partition = PetscPartitionerPartition_Gather;
1088dae52e14SToby Isaac   PetscFunctionReturn(0);
1089dae52e14SToby Isaac }
1090dae52e14SToby Isaac 
1091dae52e14SToby Isaac /*MC
1092dae52e14SToby Isaac   PETSCPARTITIONERGATHER = "gather" - A PetscPartitioner object
1093dae52e14SToby Isaac 
1094dae52e14SToby Isaac   Level: intermediate
1095dae52e14SToby Isaac 
1096dae52e14SToby Isaac .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
1097dae52e14SToby Isaac M*/
1098dae52e14SToby Isaac 
1099dae52e14SToby Isaac PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Gather(PetscPartitioner part)
1100dae52e14SToby Isaac {
1101dae52e14SToby Isaac   PetscPartitioner_Gather *p;
1102dae52e14SToby Isaac   PetscErrorCode           ierr;
1103dae52e14SToby Isaac 
1104dae52e14SToby Isaac   PetscFunctionBegin;
1105dae52e14SToby Isaac   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
1106dae52e14SToby Isaac   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
1107dae52e14SToby Isaac   part->data = p;
1108dae52e14SToby Isaac 
1109dae52e14SToby Isaac   ierr = PetscPartitionerInitialize_Gather(part);CHKERRQ(ierr);
1110dae52e14SToby Isaac   PetscFunctionReturn(0);
1111dae52e14SToby Isaac }
1112dae52e14SToby Isaac 
1113dae52e14SToby Isaac 
111477623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy_Chaco(PetscPartitioner part)
111577623264SMatthew G. Knepley {
111677623264SMatthew G. Knepley   PetscPartitioner_Chaco *p = (PetscPartitioner_Chaco *) part->data;
111777623264SMatthew G. Knepley   PetscErrorCode          ierr;
111877623264SMatthew G. Knepley 
111977623264SMatthew G. Knepley   PetscFunctionBegin;
112077623264SMatthew G. Knepley   ierr = PetscFree(p);CHKERRQ(ierr);
112177623264SMatthew G. Knepley   PetscFunctionReturn(0);
112277623264SMatthew G. Knepley }
112377623264SMatthew G. Knepley 
112477623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Chaco_Ascii(PetscPartitioner part, PetscViewer viewer)
112577623264SMatthew G. Knepley {
112677623264SMatthew G. Knepley   PetscViewerFormat format;
112777623264SMatthew G. Knepley   PetscErrorCode    ierr;
112877623264SMatthew G. Knepley 
112977623264SMatthew G. Knepley   PetscFunctionBegin;
113077623264SMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
113177623264SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Chaco Graph Partitioner:\n");CHKERRQ(ierr);
113277623264SMatthew G. Knepley   PetscFunctionReturn(0);
113377623264SMatthew G. Knepley }
113477623264SMatthew G. Knepley 
113577623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Chaco(PetscPartitioner part, PetscViewer viewer)
113677623264SMatthew G. Knepley {
113777623264SMatthew G. Knepley   PetscBool      iascii;
113877623264SMatthew G. Knepley   PetscErrorCode ierr;
113977623264SMatthew G. Knepley 
114077623264SMatthew G. Knepley   PetscFunctionBegin;
114177623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
114277623264SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
114377623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
114477623264SMatthew G. Knepley   if (iascii) {ierr = PetscPartitionerView_Chaco_Ascii(part, viewer);CHKERRQ(ierr);}
114577623264SMatthew G. Knepley   PetscFunctionReturn(0);
114677623264SMatthew G. Knepley }
114777623264SMatthew G. Knepley 
114870034214SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO)
114970034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
115070034214SMatthew G. Knepley #include <unistd.h>
115170034214SMatthew G. Knepley #endif
115211d1e910SBarry Smith #if defined(PETSC_HAVE_CHACO_INT_ASSIGNMENT)
115311d1e910SBarry Smith #include <chaco.h>
115411d1e910SBarry Smith #else
115511d1e910SBarry Smith /* Older versions of Chaco do not have an include file */
115670034214SMatthew G. Knepley PETSC_EXTERN int interface(int nvtxs, int *start, int *adjacency, int *vwgts,
115770034214SMatthew G. Knepley                        float *ewgts, float *x, float *y, float *z, char *outassignname,
115870034214SMatthew G. Knepley                        char *outfilename, short *assignment, int architecture, int ndims_tot,
115970034214SMatthew G. Knepley                        int mesh_dims[3], double *goal, int global_method, int local_method,
116070034214SMatthew G. Knepley                        int rqi_flag, int vmax, int ndims, double eigtol, long seed);
116111d1e910SBarry Smith #endif
116270034214SMatthew G. Knepley extern int FREE_GRAPH;
116377623264SMatthew G. Knepley #endif
116470034214SMatthew G. Knepley 
116577623264SMatthew G. Knepley PetscErrorCode PetscPartitionerPartition_Chaco(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
116670034214SMatthew G. Knepley {
116777623264SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO)
116870034214SMatthew G. Knepley   enum {DEFAULT_METHOD = 1, INERTIAL_METHOD = 3};
116970034214SMatthew G. Knepley   MPI_Comm       comm;
117070034214SMatthew G. Knepley   int            nvtxs          = numVertices; /* number of vertices in full graph */
117170034214SMatthew G. Knepley   int           *vwgts          = NULL;   /* weights for all vertices */
117270034214SMatthew G. Knepley   float         *ewgts          = NULL;   /* weights for all edges */
117370034214SMatthew G. Knepley   float         *x              = NULL, *y = NULL, *z = NULL; /* coordinates for inertial method */
117470034214SMatthew G. Knepley   char          *outassignname  = NULL;   /*  name of assignment output file */
117570034214SMatthew G. Knepley   char          *outfilename    = NULL;   /* output file name */
117670034214SMatthew G. Knepley   int            architecture   = 1;      /* 0 => hypercube, d => d-dimensional mesh */
117770034214SMatthew G. Knepley   int            ndims_tot      = 0;      /* total number of cube dimensions to divide */
117870034214SMatthew G. Knepley   int            mesh_dims[3];            /* dimensions of mesh of processors */
117970034214SMatthew G. Knepley   double        *goal          = NULL;    /* desired set sizes for each set */
118070034214SMatthew G. Knepley   int            global_method = 1;       /* global partitioning algorithm */
118170034214SMatthew G. Knepley   int            local_method  = 1;       /* local partitioning algorithm */
118270034214SMatthew G. Knepley   int            rqi_flag      = 0;       /* should I use RQI/Symmlq eigensolver? */
118370034214SMatthew G. Knepley   int            vmax          = 200;     /* how many vertices to coarsen down to? */
118470034214SMatthew G. Knepley   int            ndims         = 1;       /* number of eigenvectors (2^d sets) */
118570034214SMatthew G. Knepley   double         eigtol        = 0.001;   /* tolerance on eigenvectors */
118670034214SMatthew G. Knepley   long           seed          = 123636512; /* for random graph mutations */
118711d1e910SBarry Smith #if defined(PETSC_HAVE_CHACO_INT_ASSIGNMENT)
118811d1e910SBarry Smith   int           *assignment;              /* Output partition */
118911d1e910SBarry Smith #else
119070034214SMatthew G. Knepley   short int     *assignment;              /* Output partition */
119111d1e910SBarry Smith #endif
119270034214SMatthew G. Knepley   int            fd_stdout, fd_pipe[2];
119370034214SMatthew G. Knepley   PetscInt      *points;
119470034214SMatthew G. Knepley   int            i, v, p;
119570034214SMatthew G. Knepley   PetscErrorCode ierr;
119670034214SMatthew G. Knepley 
119770034214SMatthew G. Knepley   PetscFunctionBegin;
119870034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
119970034214SMatthew G. Knepley   if (!numVertices) {
120077623264SMatthew G. Knepley     ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
120177623264SMatthew G. Knepley     ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
120270034214SMatthew G. Knepley     ierr = ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
120370034214SMatthew G. Knepley     PetscFunctionReturn(0);
120470034214SMatthew G. Knepley   }
120570034214SMatthew G. Knepley   FREE_GRAPH = 0;                         /* Do not let Chaco free my memory */
120670034214SMatthew G. Knepley   for (i = 0; i < start[numVertices]; ++i) ++adjacency[i];
120770034214SMatthew G. Knepley 
120870034214SMatthew G. Knepley   if (global_method == INERTIAL_METHOD) {
120970034214SMatthew G. Knepley     /* manager.createCellCoordinates(nvtxs, &x, &y, &z); */
121070034214SMatthew G. Knepley     SETERRQ(comm, PETSC_ERR_SUP, "Inertial partitioning not yet supported");
121170034214SMatthew G. Knepley   }
121277623264SMatthew G. Knepley   mesh_dims[0] = nparts;
121370034214SMatthew G. Knepley   mesh_dims[1] = 1;
121470034214SMatthew G. Knepley   mesh_dims[2] = 1;
121570034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &assignment);CHKERRQ(ierr);
121670034214SMatthew G. Knepley   /* Chaco outputs to stdout. We redirect this to a buffer. */
121770034214SMatthew G. Knepley   /* TODO: check error codes for UNIX calls */
121870034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
121970034214SMatthew G. Knepley   {
122070034214SMatthew G. Knepley     int piperet;
122170034214SMatthew G. Knepley     piperet = pipe(fd_pipe);
122270034214SMatthew G. Knepley     if (piperet) SETERRQ(comm,PETSC_ERR_SYS,"Could not create pipe");
122370034214SMatthew G. Knepley     fd_stdout = dup(1);
122470034214SMatthew G. Knepley     close(1);
122570034214SMatthew G. Knepley     dup2(fd_pipe[1], 1);
122670034214SMatthew G. Knepley   }
122770034214SMatthew G. Knepley #endif
122870034214SMatthew G. Knepley   ierr = interface(nvtxs, (int*) start, (int*) adjacency, vwgts, ewgts, x, y, z, outassignname, outfilename,
122970034214SMatthew G. Knepley                    assignment, architecture, ndims_tot, mesh_dims, goal, global_method, local_method, rqi_flag,
123070034214SMatthew G. Knepley                    vmax, ndims, eigtol, seed);
123170034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
123270034214SMatthew G. Knepley   {
123370034214SMatthew G. Knepley     char msgLog[10000];
123470034214SMatthew G. Knepley     int  count;
123570034214SMatthew G. Knepley 
123670034214SMatthew G. Knepley     fflush(stdout);
123770034214SMatthew G. Knepley     count = read(fd_pipe[0], msgLog, (10000-1)*sizeof(char));
123870034214SMatthew G. Knepley     if (count < 0) count = 0;
123970034214SMatthew G. Knepley     msgLog[count] = 0;
124070034214SMatthew G. Knepley     close(1);
124170034214SMatthew G. Knepley     dup2(fd_stdout, 1);
124270034214SMatthew G. Knepley     close(fd_stdout);
124370034214SMatthew G. Knepley     close(fd_pipe[0]);
124470034214SMatthew G. Knepley     close(fd_pipe[1]);
124570034214SMatthew G. Knepley     if (ierr) SETERRQ1(comm, PETSC_ERR_LIB, "Error in Chaco library: %s", msgLog);
124670034214SMatthew G. Knepley   }
124770034214SMatthew G. Knepley #endif
124870034214SMatthew G. Knepley   /* Convert to PetscSection+IS */
124977623264SMatthew G. Knepley   ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
125070034214SMatthew G. Knepley   for (v = 0; v < nvtxs; ++v) {
125177623264SMatthew G. Knepley     ierr = PetscSectionAddDof(partSection, assignment[v], 1);CHKERRQ(ierr);
125270034214SMatthew G. Knepley   }
125377623264SMatthew G. Knepley   ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
125470034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr);
125577623264SMatthew G. Knepley   for (p = 0, i = 0; p < nparts; ++p) {
125670034214SMatthew G. Knepley     for (v = 0; v < nvtxs; ++v) {
125770034214SMatthew G. Knepley       if (assignment[v] == p) points[i++] = v;
125870034214SMatthew G. Knepley     }
125970034214SMatthew G. Knepley   }
126070034214SMatthew G. Knepley   if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs);
126170034214SMatthew G. Knepley   ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
126270034214SMatthew G. Knepley   if (global_method == INERTIAL_METHOD) {
126370034214SMatthew G. Knepley     /* manager.destroyCellCoordinates(nvtxs, &x, &y, &z); */
126470034214SMatthew G. Knepley   }
126570034214SMatthew G. Knepley   ierr = PetscFree(assignment);CHKERRQ(ierr);
126670034214SMatthew G. Knepley   for (i = 0; i < start[numVertices]; ++i) --adjacency[i];
126770034214SMatthew G. Knepley   PetscFunctionReturn(0);
126877623264SMatthew G. Knepley #else
126977623264SMatthew G. Knepley   SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-chaco.");
127070034214SMatthew G. Knepley #endif
127177623264SMatthew G. Knepley }
127277623264SMatthew G. Knepley 
127377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerInitialize_Chaco(PetscPartitioner part)
127477623264SMatthew G. Knepley {
127577623264SMatthew G. Knepley   PetscFunctionBegin;
127677623264SMatthew G. Knepley   part->ops->view      = PetscPartitionerView_Chaco;
127777623264SMatthew G. Knepley   part->ops->destroy   = PetscPartitionerDestroy_Chaco;
127877623264SMatthew G. Knepley   part->ops->partition = PetscPartitionerPartition_Chaco;
127977623264SMatthew G. Knepley   PetscFunctionReturn(0);
128077623264SMatthew G. Knepley }
128177623264SMatthew G. Knepley 
128277623264SMatthew G. Knepley /*MC
128377623264SMatthew G. Knepley   PETSCPARTITIONERCHACO = "chaco" - A PetscPartitioner object using the Chaco library
128477623264SMatthew G. Knepley 
128577623264SMatthew G. Knepley   Level: intermediate
128677623264SMatthew G. Knepley 
128777623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
128877623264SMatthew G. Knepley M*/
128977623264SMatthew G. Knepley 
129077623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Chaco(PetscPartitioner part)
129177623264SMatthew G. Knepley {
129277623264SMatthew G. Knepley   PetscPartitioner_Chaco *p;
129377623264SMatthew G. Knepley   PetscErrorCode          ierr;
129477623264SMatthew G. Knepley 
129577623264SMatthew G. Knepley   PetscFunctionBegin;
129677623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
129777623264SMatthew G. Knepley   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
129877623264SMatthew G. Knepley   part->data = p;
129977623264SMatthew G. Knepley 
130077623264SMatthew G. Knepley   ierr = PetscPartitionerInitialize_Chaco(part);CHKERRQ(ierr);
130177623264SMatthew G. Knepley   ierr = PetscCitationsRegister(ChacoPartitionerCitation, &ChacoPartitionercite);CHKERRQ(ierr);
130277623264SMatthew G. Knepley   PetscFunctionReturn(0);
130377623264SMatthew G. Knepley }
130477623264SMatthew G. Knepley 
130577623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy_ParMetis(PetscPartitioner part)
130677623264SMatthew G. Knepley {
130777623264SMatthew G. Knepley   PetscPartitioner_ParMetis *p = (PetscPartitioner_ParMetis *) part->data;
130877623264SMatthew G. Knepley   PetscErrorCode             ierr;
130977623264SMatthew G. Knepley 
131077623264SMatthew G. Knepley   PetscFunctionBegin;
131177623264SMatthew G. Knepley   ierr = PetscFree(p);CHKERRQ(ierr);
131277623264SMatthew G. Knepley   PetscFunctionReturn(0);
131377623264SMatthew G. Knepley }
131477623264SMatthew G. Knepley 
131577623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_ParMetis_Ascii(PetscPartitioner part, PetscViewer viewer)
131677623264SMatthew G. Knepley {
131777623264SMatthew G. Knepley   PetscViewerFormat format;
131877623264SMatthew G. Knepley   PetscErrorCode    ierr;
131977623264SMatthew G. Knepley 
132077623264SMatthew G. Knepley   PetscFunctionBegin;
132177623264SMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
132277623264SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "ParMetis Graph Partitioner:\n");CHKERRQ(ierr);
132377623264SMatthew G. Knepley   PetscFunctionReturn(0);
132477623264SMatthew G. Knepley }
132577623264SMatthew G. Knepley 
132677623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_ParMetis(PetscPartitioner part, PetscViewer viewer)
132777623264SMatthew G. Knepley {
132877623264SMatthew G. Knepley   PetscBool      iascii;
132977623264SMatthew G. Knepley   PetscErrorCode ierr;
133077623264SMatthew G. Knepley 
133177623264SMatthew G. Knepley   PetscFunctionBegin;
133277623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
133377623264SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
133477623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
133577623264SMatthew G. Knepley   if (iascii) {ierr = PetscPartitionerView_ParMetis_Ascii(part, viewer);CHKERRQ(ierr);}
133677623264SMatthew G. Knepley   PetscFunctionReturn(0);
133777623264SMatthew G. Knepley }
133870034214SMatthew G. Knepley 
133970034214SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS)
134070034214SMatthew G. Knepley #include <parmetis.h>
134177623264SMatthew G. Knepley #endif
134270034214SMatthew G. Knepley 
134377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerPartition_ParMetis(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
134470034214SMatthew G. Knepley {
134577623264SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS)
134670034214SMatthew G. Knepley   MPI_Comm       comm;
1347deea36a5SMatthew G. Knepley   PetscSection   section;
134870034214SMatthew G. Knepley   PetscInt       nvtxs      = numVertices; /* The number of vertices in full graph */
134970034214SMatthew G. Knepley   PetscInt      *vtxdist;                  /* Distribution of vertices across processes */
135070034214SMatthew G. Knepley   PetscInt      *xadj       = start;       /* Start of edge list for each vertex */
135170034214SMatthew G. Knepley   PetscInt      *adjncy     = adjacency;   /* Edge lists for all vertices */
135270034214SMatthew G. Knepley   PetscInt      *vwgt       = NULL;        /* Vertex weights */
135370034214SMatthew G. Knepley   PetscInt      *adjwgt     = NULL;        /* Edge weights */
135470034214SMatthew G. Knepley   PetscInt       wgtflag    = 0;           /* Indicates which weights are present */
135570034214SMatthew G. Knepley   PetscInt       numflag    = 0;           /* Indicates initial offset (0 or 1) */
135670034214SMatthew G. Knepley   PetscInt       ncon       = 1;           /* The number of weights per vertex */
135770034214SMatthew G. Knepley   PetscReal     *tpwgts;                   /* The fraction of vertex weights assigned to each partition */
135870034214SMatthew G. Knepley   PetscReal     *ubvec;                    /* The balance intolerance for vertex weights */
135970034214SMatthew G. Knepley   PetscInt       options[5];               /* Options */
136070034214SMatthew G. Knepley   /* Outputs */
136170034214SMatthew G. Knepley   PetscInt       edgeCut;                  /* The number of edges cut by the partition */
136270034214SMatthew G. Knepley   PetscInt      *assignment, *points;
136377623264SMatthew G. Knepley   PetscMPIInt    rank, p, v, i;
136470034214SMatthew G. Knepley   PetscErrorCode ierr;
136570034214SMatthew G. Knepley 
136670034214SMatthew G. Knepley   PetscFunctionBegin;
136777623264SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) part, &comm);CHKERRQ(ierr);
136870034214SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
136970034214SMatthew G. Knepley   options[0] = 0; /* Use all defaults */
137070034214SMatthew G. Knepley   /* Calculate vertex distribution */
1371cd0de0f2SShri   ierr = PetscMalloc5(nparts+1,&vtxdist,nparts*ncon,&tpwgts,ncon,&ubvec,nvtxs,&assignment,nvtxs,&vwgt);CHKERRQ(ierr);
137270034214SMatthew G. Knepley   vtxdist[0] = 0;
137370034214SMatthew G. Knepley   ierr = MPI_Allgather(&nvtxs, 1, MPIU_INT, &vtxdist[1], 1, MPIU_INT, comm);CHKERRQ(ierr);
137470034214SMatthew G. Knepley   for (p = 2; p <= nparts; ++p) {
137570034214SMatthew G. Knepley     vtxdist[p] += vtxdist[p-1];
137670034214SMatthew G. Knepley   }
137770034214SMatthew G. Knepley   /* Calculate weights */
137870034214SMatthew G. Knepley   for (p = 0; p < nparts; ++p) {
137970034214SMatthew G. Knepley     tpwgts[p] = 1.0/nparts;
138070034214SMatthew G. Knepley   }
138170034214SMatthew G. Knepley   ubvec[0] = 1.05;
1382deea36a5SMatthew G. Knepley   /* Weight cells by dofs on cell by default */
1383deea36a5SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1384deea36a5SMatthew G. Knepley   if (section) {
1385deea36a5SMatthew G. Knepley     PetscInt cStart, cEnd, dof;
138670034214SMatthew G. Knepley 
1387deea36a5SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1388deea36a5SMatthew G. Knepley     for (v = cStart; v < cEnd; ++v) {
1389deea36a5SMatthew G. Knepley       ierr = PetscSectionGetDof(section, v, &dof);CHKERRQ(ierr);
1390deea36a5SMatthew G. Knepley       vwgt[v-cStart] = PetscMax(dof, 1);
1391deea36a5SMatthew G. Knepley     }
1392deea36a5SMatthew G. Knepley   } else {
1393deea36a5SMatthew G. Knepley     for (v = 0; v < nvtxs; ++v) vwgt[v] = 1;
1394cd0de0f2SShri   }
1395cd0de0f2SShri 
139670034214SMatthew G. Knepley   if (nparts == 1) {
13979fc93327SToby Isaac     ierr = PetscMemzero(assignment, nvtxs * sizeof(PetscInt));CHKERRQ(ierr);
139870034214SMatthew G. Knepley   } else {
139970034214SMatthew G. Knepley     if (vtxdist[1] == vtxdist[nparts]) {
140070034214SMatthew G. Knepley       if (!rank) {
140170034214SMatthew G. Knepley         PetscStackPush("METIS_PartGraphKway");
140270034214SMatthew G. Knepley         ierr = METIS_PartGraphKway(&nvtxs, &ncon, xadj, adjncy, vwgt, NULL, adjwgt, &nparts, tpwgts, ubvec, NULL, &edgeCut, assignment);
140370034214SMatthew G. Knepley         PetscStackPop;
140470034214SMatthew G. Knepley         if (ierr != METIS_OK) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in METIS_PartGraphKway()");
140570034214SMatthew G. Knepley       }
140670034214SMatthew G. Knepley     } else {
140770034214SMatthew G. Knepley       PetscStackPush("ParMETIS_V3_PartKway");
140870034214SMatthew G. Knepley       ierr = ParMETIS_V3_PartKway(vtxdist, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgeCut, assignment, &comm);
140970034214SMatthew G. Knepley       PetscStackPop;
141070034214SMatthew G. Knepley       if (ierr != METIS_OK) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in ParMETIS_V3_PartKway()");
141170034214SMatthew G. Knepley     }
141270034214SMatthew G. Knepley   }
141370034214SMatthew G. Knepley   /* Convert to PetscSection+IS */
141477623264SMatthew G. Knepley   ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
141577623264SMatthew G. Knepley   for (v = 0; v < nvtxs; ++v) {ierr = PetscSectionAddDof(partSection, assignment[v], 1);CHKERRQ(ierr);}
141677623264SMatthew G. Knepley   ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
141770034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr);
141877623264SMatthew G. Knepley   for (p = 0, i = 0; p < nparts; ++p) {
141970034214SMatthew G. Knepley     for (v = 0; v < nvtxs; ++v) {
142070034214SMatthew G. Knepley       if (assignment[v] == p) points[i++] = v;
142170034214SMatthew G. Knepley     }
142270034214SMatthew G. Knepley   }
142370034214SMatthew G. Knepley   if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs);
142470034214SMatthew G. Knepley   ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
1425cd0de0f2SShri   ierr = PetscFree5(vtxdist,tpwgts,ubvec,assignment,vwgt);CHKERRQ(ierr);
14269b80ac48SMatthew G. Knepley   PetscFunctionReturn(0);
142770034214SMatthew G. Knepley #else
142877623264SMatthew G. Knepley   SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-parmetis.");
142970034214SMatthew G. Knepley #endif
143070034214SMatthew G. Knepley }
143170034214SMatthew G. Knepley 
143277623264SMatthew G. Knepley PetscErrorCode PetscPartitionerInitialize_ParMetis(PetscPartitioner part)
143377623264SMatthew G. Knepley {
143477623264SMatthew G. Knepley   PetscFunctionBegin;
143577623264SMatthew G. Knepley   part->ops->view      = PetscPartitionerView_ParMetis;
143677623264SMatthew G. Knepley   part->ops->destroy   = PetscPartitionerDestroy_ParMetis;
143777623264SMatthew G. Knepley   part->ops->partition = PetscPartitionerPartition_ParMetis;
143877623264SMatthew G. Knepley   PetscFunctionReturn(0);
143977623264SMatthew G. Knepley }
144077623264SMatthew G. Knepley 
144177623264SMatthew G. Knepley /*MC
144277623264SMatthew G. Knepley   PETSCPARTITIONERPARMETIS = "parmetis" - A PetscPartitioner object using the ParMetis library
144377623264SMatthew G. Knepley 
144477623264SMatthew G. Knepley   Level: intermediate
144577623264SMatthew G. Knepley 
144677623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
144777623264SMatthew G. Knepley M*/
144877623264SMatthew G. Knepley 
144977623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_ParMetis(PetscPartitioner part)
145077623264SMatthew G. Knepley {
145177623264SMatthew G. Knepley   PetscPartitioner_ParMetis *p;
145277623264SMatthew G. Knepley   PetscErrorCode          ierr;
145377623264SMatthew G. Knepley 
145477623264SMatthew G. Knepley   PetscFunctionBegin;
145577623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
145677623264SMatthew G. Knepley   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
145777623264SMatthew G. Knepley   part->data = p;
145877623264SMatthew G. Knepley 
145977623264SMatthew G. Knepley   ierr = PetscPartitionerInitialize_ParMetis(part);CHKERRQ(ierr);
146077623264SMatthew G. Knepley   ierr = PetscCitationsRegister(ParMetisPartitionerCitation, &ParMetisPartitionercite);CHKERRQ(ierr);
146170034214SMatthew G. Knepley   PetscFunctionReturn(0);
146270034214SMatthew G. Knepley }
146370034214SMatthew G. Knepley 
14645680f57bSMatthew G. Knepley /*@
14655680f57bSMatthew G. Knepley   DMPlexGetPartitioner - Get the mesh partitioner
14665680f57bSMatthew G. Knepley 
14675680f57bSMatthew G. Knepley   Not collective
14685680f57bSMatthew G. Knepley 
14695680f57bSMatthew G. Knepley   Input Parameter:
14705680f57bSMatthew G. Knepley . dm - The DM
14715680f57bSMatthew G. Knepley 
14725680f57bSMatthew G. Knepley   Output Parameter:
14735680f57bSMatthew G. Knepley . part - The PetscPartitioner
14745680f57bSMatthew G. Knepley 
14755680f57bSMatthew G. Knepley   Level: developer
14765680f57bSMatthew G. Knepley 
147798599a47SLawrence Mitchell   Note: This gets a borrowed reference, so the user should not destroy this PetscPartitioner.
147898599a47SLawrence Mitchell 
147998599a47SLawrence Mitchell .seealso DMPlexDistribute(), DMPlexSetPartitioner(), PetscPartitionerCreate()
14805680f57bSMatthew G. Knepley @*/
14815680f57bSMatthew G. Knepley PetscErrorCode DMPlexGetPartitioner(DM dm, PetscPartitioner *part)
14825680f57bSMatthew G. Knepley {
14835680f57bSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
14845680f57bSMatthew G. Knepley 
14855680f57bSMatthew G. Knepley   PetscFunctionBegin;
14865680f57bSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14875680f57bSMatthew G. Knepley   PetscValidPointer(part, 2);
14885680f57bSMatthew G. Knepley   *part = mesh->partitioner;
14895680f57bSMatthew G. Knepley   PetscFunctionReturn(0);
14905680f57bSMatthew G. Knepley }
14915680f57bSMatthew G. Knepley 
149271bb2955SLawrence Mitchell /*@
149371bb2955SLawrence Mitchell   DMPlexSetPartitioner - Set the mesh partitioner
149471bb2955SLawrence Mitchell 
149571bb2955SLawrence Mitchell   logically collective on dm and part
149671bb2955SLawrence Mitchell 
149771bb2955SLawrence Mitchell   Input Parameters:
149871bb2955SLawrence Mitchell + dm - The DM
149971bb2955SLawrence Mitchell - part - The partitioner
150071bb2955SLawrence Mitchell 
150171bb2955SLawrence Mitchell   Level: developer
150271bb2955SLawrence Mitchell 
150371bb2955SLawrence Mitchell   Note: Any existing PetscPartitioner will be destroyed.
150471bb2955SLawrence Mitchell 
150571bb2955SLawrence Mitchell .seealso DMPlexDistribute(), DMPlexGetPartitioner(), PetscPartitionerCreate()
150671bb2955SLawrence Mitchell @*/
150771bb2955SLawrence Mitchell PetscErrorCode DMPlexSetPartitioner(DM dm, PetscPartitioner part)
150871bb2955SLawrence Mitchell {
150971bb2955SLawrence Mitchell   DM_Plex       *mesh = (DM_Plex *) dm->data;
151071bb2955SLawrence Mitchell   PetscErrorCode ierr;
151171bb2955SLawrence Mitchell 
151271bb2955SLawrence Mitchell   PetscFunctionBegin;
151371bb2955SLawrence Mitchell   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
151471bb2955SLawrence Mitchell   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 2);
151571bb2955SLawrence Mitchell   ierr = PetscObjectReference((PetscObject)part);CHKERRQ(ierr);
151671bb2955SLawrence Mitchell   ierr = PetscPartitionerDestroy(&mesh->partitioner);CHKERRQ(ierr);
151771bb2955SLawrence Mitchell   mesh->partitioner = part;
151871bb2955SLawrence Mitchell   PetscFunctionReturn(0);
151971bb2955SLawrence Mitchell }
152071bb2955SLawrence Mitchell 
15216a5a2ffdSToby Isaac static PetscErrorCode DMPlexPartitionLabelClosure_Tree(DM dm, DMLabel label, PetscInt rank, PetscInt point, PetscBool up, PetscBool down)
1522270bba0cSToby Isaac {
1523270bba0cSToby Isaac   PetscErrorCode ierr;
1524270bba0cSToby Isaac 
1525270bba0cSToby Isaac   PetscFunctionBegin;
15266a5a2ffdSToby Isaac   if (up) {
15276a5a2ffdSToby Isaac     PetscInt parent;
15286a5a2ffdSToby Isaac 
1529270bba0cSToby Isaac     ierr = DMPlexGetTreeParent(dm,point,&parent,NULL);CHKERRQ(ierr);
15306a5a2ffdSToby Isaac     if (parent != point) {
15316a5a2ffdSToby Isaac       PetscInt closureSize, *closure = NULL, i;
15326a5a2ffdSToby Isaac 
1533270bba0cSToby Isaac       ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1534270bba0cSToby Isaac       for (i = 0; i < closureSize; i++) {
1535270bba0cSToby Isaac         PetscInt cpoint = closure[2*i];
1536270bba0cSToby Isaac 
1537270bba0cSToby Isaac         ierr = DMLabelSetValue(label,cpoint,rank);CHKERRQ(ierr);
15386a5a2ffdSToby Isaac         ierr = DMPlexPartitionLabelClosure_Tree(dm,label,rank,cpoint,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1539270bba0cSToby Isaac       }
1540270bba0cSToby Isaac       ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
15416a5a2ffdSToby Isaac     }
15426a5a2ffdSToby Isaac   }
15436a5a2ffdSToby Isaac   if (down) {
15446a5a2ffdSToby Isaac     PetscInt numChildren;
15456a5a2ffdSToby Isaac     const PetscInt *children;
15466a5a2ffdSToby Isaac 
15476a5a2ffdSToby Isaac     ierr = DMPlexGetTreeChildren(dm,point,&numChildren,&children);CHKERRQ(ierr);
15486a5a2ffdSToby Isaac     if (numChildren) {
15496a5a2ffdSToby Isaac       PetscInt i;
15506a5a2ffdSToby Isaac 
15516a5a2ffdSToby Isaac       for (i = 0; i < numChildren; i++) {
15526a5a2ffdSToby Isaac         PetscInt cpoint = children[i];
15536a5a2ffdSToby Isaac 
15546a5a2ffdSToby Isaac         ierr = DMLabelSetValue(label,cpoint,rank);CHKERRQ(ierr);
15556a5a2ffdSToby Isaac         ierr = DMPlexPartitionLabelClosure_Tree(dm,label,rank,cpoint,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr);
15566a5a2ffdSToby Isaac       }
15576a5a2ffdSToby Isaac     }
15586a5a2ffdSToby Isaac   }
1559270bba0cSToby Isaac   PetscFunctionReturn(0);
1560270bba0cSToby Isaac }
1561270bba0cSToby Isaac 
15625abbe4feSMichael Lange /*@
15635abbe4feSMichael Lange   DMPlexPartitionLabelClosure - Add the closure of all points to the partition label
15645abbe4feSMichael Lange 
15655abbe4feSMichael Lange   Input Parameters:
15665abbe4feSMichael Lange + dm     - The DM
15675abbe4feSMichael Lange - label  - DMLabel assinging ranks to remote roots
15685abbe4feSMichael Lange 
15695abbe4feSMichael Lange   Level: developer
15705abbe4feSMichael Lange 
15715abbe4feSMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap
15725abbe4feSMichael Lange @*/
15735abbe4feSMichael Lange PetscErrorCode DMPlexPartitionLabelClosure(DM dm, DMLabel label)
15749ffc88e4SToby Isaac {
15755abbe4feSMichael Lange   IS              rankIS,   pointIS;
15765abbe4feSMichael Lange   const PetscInt *ranks,   *points;
15775abbe4feSMichael Lange   PetscInt        numRanks, numPoints, r, p, c, closureSize;
15785abbe4feSMichael Lange   PetscInt       *closure = NULL;
15799ffc88e4SToby Isaac   PetscErrorCode  ierr;
15809ffc88e4SToby Isaac 
15819ffc88e4SToby Isaac   PetscFunctionBegin;
15825abbe4feSMichael Lange   ierr = DMLabelGetValueIS(label, &rankIS);CHKERRQ(ierr);
15835abbe4feSMichael Lange   ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr);
15845abbe4feSMichael Lange   ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr);
15855abbe4feSMichael Lange   for (r = 0; r < numRanks; ++r) {
15865abbe4feSMichael Lange     const PetscInt rank = ranks[r];
15879ffc88e4SToby Isaac 
15885abbe4feSMichael Lange     ierr = DMLabelGetStratumIS(label, rank, &pointIS);CHKERRQ(ierr);
15895abbe4feSMichael Lange     ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
15905abbe4feSMichael Lange     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
15915abbe4feSMichael Lange     for (p = 0; p < numPoints; ++p) {
15925abbe4feSMichael Lange       ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1593270bba0cSToby Isaac       for (c = 0; c < closureSize*2; c += 2) {
1594270bba0cSToby Isaac         ierr = DMLabelSetValue(label, closure[c], rank);CHKERRQ(ierr);
15956a5a2ffdSToby Isaac         ierr = DMPlexPartitionLabelClosure_Tree(dm,label,rank,closure[c],PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr);
1596270bba0cSToby Isaac       }
15979ffc88e4SToby Isaac     }
15985abbe4feSMichael Lange     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
15995abbe4feSMichael Lange     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
16009ffc88e4SToby Isaac   }
16017de78196SMichael Lange   if (closure) {ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);}
16025abbe4feSMichael Lange   ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr);
16035abbe4feSMichael Lange   ierr = ISDestroy(&rankIS);CHKERRQ(ierr);
16049ffc88e4SToby Isaac   PetscFunctionReturn(0);
16059ffc88e4SToby Isaac }
16069ffc88e4SToby Isaac 
160724d039d7SMichael Lange /*@
160824d039d7SMichael Lange   DMPlexPartitionLabelAdjacency - Add one level of adjacent points to the partition label
160924d039d7SMichael Lange 
161024d039d7SMichael Lange   Input Parameters:
161124d039d7SMichael Lange + dm     - The DM
161224d039d7SMichael Lange - label  - DMLabel assinging ranks to remote roots
161324d039d7SMichael Lange 
161424d039d7SMichael Lange   Level: developer
161524d039d7SMichael Lange 
161624d039d7SMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap
161724d039d7SMichael Lange @*/
161824d039d7SMichael Lange PetscErrorCode DMPlexPartitionLabelAdjacency(DM dm, DMLabel label)
161970034214SMatthew G. Knepley {
162024d039d7SMichael Lange   IS              rankIS,   pointIS;
162124d039d7SMichael Lange   const PetscInt *ranks,   *points;
162224d039d7SMichael Lange   PetscInt        numRanks, numPoints, r, p, a, adjSize;
162324d039d7SMichael Lange   PetscInt       *adj = NULL;
162470034214SMatthew G. Knepley   PetscErrorCode  ierr;
162570034214SMatthew G. Knepley 
162670034214SMatthew G. Knepley   PetscFunctionBegin;
162724d039d7SMichael Lange   ierr = DMLabelGetValueIS(label, &rankIS);CHKERRQ(ierr);
162824d039d7SMichael Lange   ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr);
162924d039d7SMichael Lange   ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr);
163024d039d7SMichael Lange   for (r = 0; r < numRanks; ++r) {
163124d039d7SMichael Lange     const PetscInt rank = ranks[r];
163270034214SMatthew G. Knepley 
163324d039d7SMichael Lange     ierr = DMLabelGetStratumIS(label, rank, &pointIS);CHKERRQ(ierr);
163424d039d7SMichael Lange     ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
163524d039d7SMichael Lange     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
163670034214SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
163724d039d7SMichael Lange       adjSize = PETSC_DETERMINE;
163824d039d7SMichael Lange       ierr = DMPlexGetAdjacency(dm, points[p], &adjSize, &adj);CHKERRQ(ierr);
163924d039d7SMichael Lange       for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(label, adj[a], rank);CHKERRQ(ierr);}
164070034214SMatthew G. Knepley     }
164124d039d7SMichael Lange     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
164224d039d7SMichael Lange     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
164370034214SMatthew G. Knepley   }
164424d039d7SMichael Lange   ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr);
164524d039d7SMichael Lange   ierr = ISDestroy(&rankIS);CHKERRQ(ierr);
164624d039d7SMichael Lange   ierr = PetscFree(adj);CHKERRQ(ierr);
164724d039d7SMichael Lange   PetscFunctionReturn(0);
164870034214SMatthew G. Knepley }
164970034214SMatthew G. Knepley 
1650be200f8dSMichael Lange /*@
1651be200f8dSMichael Lange   DMPlexPartitionLabelPropagate - Propagate points in a partition label over the point SF
1652be200f8dSMichael Lange 
1653be200f8dSMichael Lange   Input Parameters:
1654be200f8dSMichael Lange + dm     - The DM
1655be200f8dSMichael Lange - label  - DMLabel assinging ranks to remote roots
1656be200f8dSMichael Lange 
1657be200f8dSMichael Lange   Level: developer
1658be200f8dSMichael Lange 
1659be200f8dSMichael Lange   Note: This is required when generating multi-level overlaps to capture
1660be200f8dSMichael Lange   overlap points from non-neighbouring partitions.
1661be200f8dSMichael Lange 
1662be200f8dSMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap
1663be200f8dSMichael Lange @*/
1664be200f8dSMichael Lange PetscErrorCode DMPlexPartitionLabelPropagate(DM dm, DMLabel label)
1665be200f8dSMichael Lange {
1666be200f8dSMichael Lange   MPI_Comm        comm;
1667be200f8dSMichael Lange   PetscMPIInt     rank;
1668be200f8dSMichael Lange   PetscSF         sfPoint;
16695d04f6ebSMichael Lange   DMLabel         lblRoots, lblLeaves;
1670be200f8dSMichael Lange   IS              rankIS, pointIS;
1671be200f8dSMichael Lange   const PetscInt *ranks;
1672be200f8dSMichael Lange   PetscInt        numRanks, r;
1673be200f8dSMichael Lange   PetscErrorCode  ierr;
1674be200f8dSMichael Lange 
1675be200f8dSMichael Lange   PetscFunctionBegin;
1676be200f8dSMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
1677be200f8dSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1678be200f8dSMichael Lange   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
16795d04f6ebSMichael Lange   /* Pull point contributions from remote leaves into local roots */
16805d04f6ebSMichael Lange   ierr = DMLabelGather(label, sfPoint, &lblLeaves);CHKERRQ(ierr);
16815d04f6ebSMichael Lange   ierr = DMLabelGetValueIS(lblLeaves, &rankIS);CHKERRQ(ierr);
16825d04f6ebSMichael Lange   ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr);
16835d04f6ebSMichael Lange   ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr);
16845d04f6ebSMichael Lange   for (r = 0; r < numRanks; ++r) {
16855d04f6ebSMichael Lange     const PetscInt remoteRank = ranks[r];
16865d04f6ebSMichael Lange     if (remoteRank == rank) continue;
16875d04f6ebSMichael Lange     ierr = DMLabelGetStratumIS(lblLeaves, remoteRank, &pointIS);CHKERRQ(ierr);
16885d04f6ebSMichael Lange     ierr = DMLabelInsertIS(label, pointIS, remoteRank);CHKERRQ(ierr);
16895d04f6ebSMichael Lange     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
16905d04f6ebSMichael Lange   }
16915d04f6ebSMichael Lange   ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr);
16925d04f6ebSMichael Lange   ierr = ISDestroy(&rankIS);CHKERRQ(ierr);
16935d04f6ebSMichael Lange   ierr = DMLabelDestroy(&lblLeaves);CHKERRQ(ierr);
1694be200f8dSMichael Lange   /* Push point contributions from roots into remote leaves */
1695be200f8dSMichael Lange   ierr = DMLabelDistribute(label, sfPoint, &lblRoots);CHKERRQ(ierr);
1696be200f8dSMichael Lange   ierr = DMLabelGetValueIS(lblRoots, &rankIS);CHKERRQ(ierr);
1697be200f8dSMichael Lange   ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr);
1698be200f8dSMichael Lange   ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr);
1699be200f8dSMichael Lange   for (r = 0; r < numRanks; ++r) {
1700be200f8dSMichael Lange     const PetscInt remoteRank = ranks[r];
1701be200f8dSMichael Lange     if (remoteRank == rank) continue;
1702be200f8dSMichael Lange     ierr = DMLabelGetStratumIS(lblRoots, remoteRank, &pointIS);CHKERRQ(ierr);
1703be200f8dSMichael Lange     ierr = DMLabelInsertIS(label, pointIS, remoteRank);CHKERRQ(ierr);
1704be200f8dSMichael Lange     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1705be200f8dSMichael Lange   }
1706be200f8dSMichael Lange   ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr);
1707be200f8dSMichael Lange   ierr = ISDestroy(&rankIS);CHKERRQ(ierr);
1708be200f8dSMichael Lange   ierr = DMLabelDestroy(&lblRoots);CHKERRQ(ierr);
1709be200f8dSMichael Lange   PetscFunctionReturn(0);
1710be200f8dSMichael Lange }
1711be200f8dSMichael Lange 
17121fd9873aSMichael Lange /*@
17131fd9873aSMichael Lange   DMPlexPartitionLabelInvert - Create a partition label of remote roots from a local root label
17141fd9873aSMichael Lange 
17151fd9873aSMichael Lange   Input Parameters:
17161fd9873aSMichael Lange + dm        - The DM
17171fd9873aSMichael Lange . rootLabel - DMLabel assinging ranks to local roots
17181fd9873aSMichael Lange . processSF - A star forest mapping into the local index on each remote rank
17191fd9873aSMichael Lange 
17201fd9873aSMichael Lange   Output Parameter:
17211fd9873aSMichael Lange - leafLabel - DMLabel assinging ranks to remote roots
17221fd9873aSMichael Lange 
17231fd9873aSMichael Lange   Note: The rootLabel defines a send pattern by mapping local points to remote target ranks. The
17241fd9873aSMichael Lange   resulting leafLabel is a receiver mapping of remote roots to their parent rank.
17251fd9873aSMichael Lange 
17261fd9873aSMichael Lange   Level: developer
17271fd9873aSMichael Lange 
17281fd9873aSMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap
17291fd9873aSMichael Lange @*/
17301fd9873aSMichael Lange PetscErrorCode DMPlexPartitionLabelInvert(DM dm, DMLabel rootLabel, PetscSF processSF, DMLabel leafLabel)
17311fd9873aSMichael Lange {
17321fd9873aSMichael Lange   MPI_Comm           comm;
1733*9852e123SBarry Smith   PetscMPIInt        rank, size;
1734*9852e123SBarry Smith   PetscInt           p, n, numNeighbors, ssize, l, nleaves;
17351fd9873aSMichael Lange   PetscSF            sfPoint;
17361fd9873aSMichael Lange   PetscSFNode       *rootPoints, *leafPoints;
17371fd9873aSMichael Lange   PetscSection       rootSection, leafSection;
17381fd9873aSMichael Lange   const PetscSFNode *remote;
17391fd9873aSMichael Lange   const PetscInt    *local, *neighbors;
17401fd9873aSMichael Lange   IS                 valueIS;
17411fd9873aSMichael Lange   PetscErrorCode     ierr;
17421fd9873aSMichael Lange 
17431fd9873aSMichael Lange   PetscFunctionBegin;
17441fd9873aSMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
17451fd9873aSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1746*9852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
17471fd9873aSMichael Lange   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
17481fd9873aSMichael Lange 
17491fd9873aSMichael Lange   /* Convert to (point, rank) and use actual owners */
17501fd9873aSMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
1751*9852e123SBarry Smith   ierr = PetscSectionSetChart(rootSection, 0, size);CHKERRQ(ierr);
17521fd9873aSMichael Lange   ierr = DMLabelGetValueIS(rootLabel, &valueIS);CHKERRQ(ierr);
17531fd9873aSMichael Lange   ierr = ISGetLocalSize(valueIS, &numNeighbors);CHKERRQ(ierr);
17541fd9873aSMichael Lange   ierr = ISGetIndices(valueIS, &neighbors);CHKERRQ(ierr);
17551fd9873aSMichael Lange   for (n = 0; n < numNeighbors; ++n) {
17561fd9873aSMichael Lange     PetscInt numPoints;
17571fd9873aSMichael Lange 
17581fd9873aSMichael Lange     ierr = DMLabelGetStratumSize(rootLabel, neighbors[n], &numPoints);CHKERRQ(ierr);
17591fd9873aSMichael Lange     ierr = PetscSectionAddDof(rootSection, neighbors[n], numPoints);CHKERRQ(ierr);
17601fd9873aSMichael Lange   }
17611fd9873aSMichael Lange   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
1762*9852e123SBarry Smith   ierr = PetscSectionGetStorageSize(rootSection, &ssize);CHKERRQ(ierr);
1763*9852e123SBarry Smith   ierr = PetscMalloc1(ssize, &rootPoints);CHKERRQ(ierr);
17641fd9873aSMichael Lange   ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr);
17651fd9873aSMichael Lange   for (n = 0; n < numNeighbors; ++n) {
17661fd9873aSMichael Lange     IS              pointIS;
17671fd9873aSMichael Lange     const PetscInt *points;
17681fd9873aSMichael Lange     PetscInt        off, numPoints, p;
17691fd9873aSMichael Lange 
17701fd9873aSMichael Lange     ierr = PetscSectionGetOffset(rootSection, neighbors[n], &off);CHKERRQ(ierr);
17711fd9873aSMichael Lange     ierr = DMLabelGetStratumIS(rootLabel, neighbors[n], &pointIS);CHKERRQ(ierr);
17721fd9873aSMichael Lange     ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
17731fd9873aSMichael Lange     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
17741fd9873aSMichael Lange     for (p = 0; p < numPoints; ++p) {
1775f8987ae8SMichael Lange       if (local) {ierr = PetscFindInt(points[p], nleaves, local, &l);CHKERRQ(ierr);}
1776f8987ae8SMichael Lange       else       {l = -1;}
17771fd9873aSMichael Lange       if (l >= 0) {rootPoints[off+p] = remote[l];}
17781fd9873aSMichael Lange       else        {rootPoints[off+p].index = points[p]; rootPoints[off+p].rank = rank;}
17791fd9873aSMichael Lange     }
17801fd9873aSMichael Lange     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
17811fd9873aSMichael Lange     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
17821fd9873aSMichael Lange   }
17831fd9873aSMichael Lange   ierr = ISRestoreIndices(valueIS, &neighbors);CHKERRQ(ierr);
17841fd9873aSMichael Lange   ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
17851fd9873aSMichael Lange   /* Communicate overlap */
17861fd9873aSMichael Lange   ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr);
17871fd9873aSMichael Lange   ierr = DMPlexDistributeData(dm, processSF, rootSection, MPIU_2INT, rootPoints, leafSection, (void**) &leafPoints);CHKERRQ(ierr);
17881fd9873aSMichael Lange   /* Filter remote contributions (ovLeafPoints) into the overlapSF */
1789*9852e123SBarry Smith   ierr = PetscSectionGetStorageSize(leafSection, &ssize);CHKERRQ(ierr);
1790*9852e123SBarry Smith   for (p = 0; p < ssize; p++) {
17911fd9873aSMichael Lange     ierr = DMLabelSetValue(leafLabel, leafPoints[p].index, leafPoints[p].rank);CHKERRQ(ierr);
17921fd9873aSMichael Lange   }
17931fd9873aSMichael Lange   ierr = PetscFree(rootPoints);CHKERRQ(ierr);
17941fd9873aSMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
17951fd9873aSMichael Lange   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
17961fd9873aSMichael Lange   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
17971fd9873aSMichael Lange   PetscFunctionReturn(0);
17981fd9873aSMichael Lange }
17991fd9873aSMichael Lange 
1800aa3148a8SMichael Lange /*@
1801aa3148a8SMichael Lange   DMPlexPartitionLabelCreateSF - Create a star forest from a label that assigns ranks to points
1802aa3148a8SMichael Lange 
1803aa3148a8SMichael Lange   Input Parameters:
1804aa3148a8SMichael Lange + dm    - The DM
1805aa3148a8SMichael Lange . label - DMLabel assinging ranks to remote roots
1806aa3148a8SMichael Lange 
1807aa3148a8SMichael Lange   Output Parameter:
1808aa3148a8SMichael Lange - sf    - The star forest communication context encapsulating the defined mapping
1809aa3148a8SMichael Lange 
1810aa3148a8SMichael Lange   Note: The incoming label is a receiver mapping of remote points to their parent rank.
1811aa3148a8SMichael Lange 
1812aa3148a8SMichael Lange   Level: developer
1813aa3148a8SMichael Lange 
1814aa3148a8SMichael Lange .seealso: DMPlexDistribute(), DMPlexCreateOverlap
1815aa3148a8SMichael Lange @*/
1816aa3148a8SMichael Lange PetscErrorCode DMPlexPartitionLabelCreateSF(DM dm, DMLabel label, PetscSF *sf)
1817aa3148a8SMichael Lange {
1818*9852e123SBarry Smith   PetscMPIInt     rank, size;
181943f7d02bSMichael Lange   PetscInt        n, numRemote, p, numPoints, pStart, pEnd, idx = 0;
1820aa3148a8SMichael Lange   PetscSFNode    *remotePoints;
182143f7d02bSMichael Lange   IS              remoteRootIS;
182243f7d02bSMichael Lange   const PetscInt *remoteRoots;
1823aa3148a8SMichael Lange   PetscErrorCode ierr;
1824aa3148a8SMichael Lange 
1825aa3148a8SMichael Lange   PetscFunctionBegin;
182643f7d02bSMichael Lange   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
1827*9852e123SBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr);
1828aa3148a8SMichael Lange 
1829*9852e123SBarry Smith   for (numRemote = 0, n = 0; n < size; ++n) {
1830aa3148a8SMichael Lange     ierr = DMLabelGetStratumSize(label, n, &numPoints);CHKERRQ(ierr);
1831aa3148a8SMichael Lange     numRemote += numPoints;
1832aa3148a8SMichael Lange   }
1833aa3148a8SMichael Lange   ierr = PetscMalloc1(numRemote, &remotePoints);CHKERRQ(ierr);
183443f7d02bSMichael Lange   /* Put owned points first */
183543f7d02bSMichael Lange   ierr = DMLabelGetStratumSize(label, rank, &numPoints);CHKERRQ(ierr);
183643f7d02bSMichael Lange   if (numPoints > 0) {
183743f7d02bSMichael Lange     ierr = DMLabelGetStratumIS(label, rank, &remoteRootIS);CHKERRQ(ierr);
183843f7d02bSMichael Lange     ierr = ISGetIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
183943f7d02bSMichael Lange     for (p = 0; p < numPoints; p++) {
184043f7d02bSMichael Lange       remotePoints[idx].index = remoteRoots[p];
184143f7d02bSMichael Lange       remotePoints[idx].rank = rank;
184243f7d02bSMichael Lange       idx++;
184343f7d02bSMichael Lange     }
184443f7d02bSMichael Lange     ierr = ISRestoreIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
184543f7d02bSMichael Lange     ierr = ISDestroy(&remoteRootIS);CHKERRQ(ierr);
184643f7d02bSMichael Lange   }
184743f7d02bSMichael Lange   /* Now add remote points */
1848*9852e123SBarry Smith   for (n = 0; n < size; ++n) {
1849aa3148a8SMichael Lange     ierr = DMLabelGetStratumSize(label, n, &numPoints);CHKERRQ(ierr);
185043f7d02bSMichael Lange     if (numPoints <= 0 || n == rank) continue;
1851aa3148a8SMichael Lange     ierr = DMLabelGetStratumIS(label, n, &remoteRootIS);CHKERRQ(ierr);
1852aa3148a8SMichael Lange     ierr = ISGetIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
1853aa3148a8SMichael Lange     for (p = 0; p < numPoints; p++) {
1854aa3148a8SMichael Lange       remotePoints[idx].index = remoteRoots[p];
1855aa3148a8SMichael Lange       remotePoints[idx].rank = n;
1856aa3148a8SMichael Lange       idx++;
1857aa3148a8SMichael Lange     }
1858aa3148a8SMichael Lange     ierr = ISRestoreIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
1859aa3148a8SMichael Lange     ierr = ISDestroy(&remoteRootIS);CHKERRQ(ierr);
1860aa3148a8SMichael Lange   }
1861aa3148a8SMichael Lange   ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), sf);CHKERRQ(ierr);
1862aa3148a8SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
1863aa3148a8SMichael Lange   ierr = PetscSFSetGraph(*sf, pEnd-pStart, numRemote, NULL, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
186470034214SMatthew G. Knepley   PetscFunctionReturn(0);
186570034214SMatthew G. Knepley }
1866