xref: /petsc/src/dm/impls/plex/plexpartition.c (revision 8cfe4c1f06d0c32e59a0d6d8dd36f10663364657)
170034214SMatthew G. Knepley #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 
2970034214SMatthew G. Knepley #undef __FUNCT__
30532c4e7dSMichael Lange #define __FUNCT__ "DMPlexCreatePartitionerGraph"
31532c4e7dSMichael Lange /*@C
32532c4e7dSMichael Lange   DMPlexCreatePartitionerGraph - Create a CSR graph of point connections for the partitioner
33532c4e7dSMichael Lange 
34532c4e7dSMichael Lange   Input Parameters:
35532c4e7dSMichael Lange + dm      - The mesh DM dm
36532c4e7dSMichael Lange - height  - Height of the strata from which to construct the graph
37532c4e7dSMichael Lange 
38532c4e7dSMichael Lange   Output Parameter:
39532c4e7dSMichael Lange + numVertices - Number of vertices in the graph
40532c4e7dSMichael Lange - offsets     - Point offsets in the graph
41532c4e7dSMichael Lange - adjacency   - Point connectivity in the graph
42532c4e7dSMichael Lange 
43532c4e7dSMichael Lange   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
44532c4e7dSMichael Lange   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
45532c4e7dSMichael Lange   representation on the mesh.
46532c4e7dSMichael Lange 
47532c4e7dSMichael Lange   Level: developer
48532c4e7dSMichael Lange 
49532c4e7dSMichael Lange .seealso: PetscPartitionerGetType(), PetscPartitionerCreate(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
50532c4e7dSMichael Lange @*/
51532c4e7dSMichael Lange PetscErrorCode DMPlexCreatePartitionerGraph(DM dm, PetscInt height, PetscInt *numVertices, PetscInt **offsets, PetscInt **adjacency)
52532c4e7dSMichael Lange {
53*8cfe4c1fSMichael Lange   PetscInt       p, pStart, pEnd, a, adjSize, size, nroots;
54532c4e7dSMichael Lange   PetscInt      *adj = NULL, *vOffsets = NULL;
55*8cfe4c1fSMichael Lange   IS             cellNumbering;
56*8cfe4c1fSMichael Lange   const PetscInt *cellNum;
57532c4e7dSMichael Lange   PetscBool      useCone, useClosure;
58532c4e7dSMichael Lange   PetscSection   section;
59532c4e7dSMichael Lange   PetscSegBuffer adjBuffer;
60*8cfe4c1fSMichael Lange   PetscSF        sfPoint;
61532c4e7dSMichael Lange   PetscMPIInt    rank;
62532c4e7dSMichael Lange   PetscErrorCode ierr;
63532c4e7dSMichael Lange 
64532c4e7dSMichael Lange   PetscFunctionBegin;
65532c4e7dSMichael Lange   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
66532c4e7dSMichael Lange   ierr = DMPlexGetHeightStratum(dm, height, &pStart, &pEnd);CHKERRQ(ierr);
67*8cfe4c1fSMichael Lange   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
68*8cfe4c1fSMichael Lange   ierr = PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
69532c4e7dSMichael Lange   /* Build adjacency graph via a section/segbuffer */
70532c4e7dSMichael Lange   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &section);CHKERRQ(ierr);
71532c4e7dSMichael Lange   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
72532c4e7dSMichael Lange   ierr = PetscSegBufferCreate(sizeof(PetscInt),1000,&adjBuffer);CHKERRQ(ierr);
73532c4e7dSMichael Lange   /* Always use FVM adjacency to create partitioner graph */
74532c4e7dSMichael Lange   ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
75532c4e7dSMichael Lange   ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr);
76532c4e7dSMichael Lange   ierr = DMPlexSetAdjacencyUseCone(dm, PETSC_TRUE);CHKERRQ(ierr);
77532c4e7dSMichael Lange   ierr = DMPlexSetAdjacencyUseClosure(dm, PETSC_FALSE);CHKERRQ(ierr);
78*8cfe4c1fSMichael Lange   if (nroots > 0) {
79*8cfe4c1fSMichael Lange     ierr = DMPlexGetCellNumbering(dm, &cellNumbering);CHKERRQ(ierr);
80*8cfe4c1fSMichael Lange     ierr = ISGetIndices(cellNumbering, &cellNum);CHKERRQ(ierr);
81*8cfe4c1fSMichael Lange   }
82*8cfe4c1fSMichael Lange   for (*numVertices = 0, p = pStart; p < pEnd; p++) {
83*8cfe4c1fSMichael Lange     /* Skip non-owned cells in parallel (ParMetis expects no overlap) */
84*8cfe4c1fSMichael Lange     if (nroots > 0) {if (cellNum[p] < 0) continue;}
85532c4e7dSMichael Lange     adjSize = PETSC_DETERMINE;
86532c4e7dSMichael Lange     ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
87532c4e7dSMichael Lange     for (a = 0; a < adjSize; ++a) {
88532c4e7dSMichael Lange       const PetscInt point = adj[a];
89532c4e7dSMichael Lange       if (point != p && pStart <= point && point < pEnd) {
90532c4e7dSMichael Lange         PetscInt *PETSC_RESTRICT pBuf;
91532c4e7dSMichael Lange         ierr = PetscSectionAddDof(section, p, 1);CHKERRQ(ierr);
92532c4e7dSMichael Lange         ierr = PetscSegBufferGetInts(adjBuffer, 1, &pBuf);CHKERRQ(ierr);
93532c4e7dSMichael Lange         *pBuf = point;
94532c4e7dSMichael Lange       }
95532c4e7dSMichael Lange     }
96*8cfe4c1fSMichael Lange     (*numVertices)++;
97532c4e7dSMichael Lange   }
98532c4e7dSMichael Lange   ierr = DMPlexSetAdjacencyUseCone(dm, useCone);CHKERRQ(ierr);
99532c4e7dSMichael Lange   ierr = DMPlexSetAdjacencyUseClosure(dm, useClosure);CHKERRQ(ierr);
100532c4e7dSMichael Lange   /* Derive CSR graph from section/segbuffer */
101532c4e7dSMichael Lange   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
102532c4e7dSMichael Lange   ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
103532c4e7dSMichael Lange   ierr = PetscMalloc2(*numVertices+1, &vOffsets, size, adjacency);CHKERRQ(ierr);
104*8cfe4c1fSMichael Lange   for (p = 0; p < *numVertices; p++) {ierr = PetscSectionGetOffset(section, p, &(vOffsets[p]));CHKERRQ(ierr);}
105532c4e7dSMichael Lange   vOffsets[*numVertices] = size;
106532c4e7dSMichael Lange   if (offsets) *offsets = vOffsets;
107*8cfe4c1fSMichael Lange   if (adjacency) {
108*8cfe4c1fSMichael Lange     ierr = PetscSegBufferExtractTo(adjBuffer, *adjacency);CHKERRQ(ierr);
109*8cfe4c1fSMichael Lange     if (nroots > 0) {
110*8cfe4c1fSMichael Lange       ISLocalToGlobalMapping ltogCells;
111*8cfe4c1fSMichael Lange       PetscInt n, size, *cells_arr;
112*8cfe4c1fSMichael Lange       /* In parallel, apply a global cell numbering to the graph */
113*8cfe4c1fSMichael Lange       ierr = ISRestoreIndices(cellNumbering, &cellNum);CHKERRQ(ierr);
114*8cfe4c1fSMichael Lange       ierr = ISLocalToGlobalMappingCreateIS(cellNumbering, &ltogCells);CHKERRQ(ierr);
115*8cfe4c1fSMichael Lange       ierr = ISLocalToGlobalMappingGetSize(ltogCells, &size);CHKERRQ(ierr);
116*8cfe4c1fSMichael Lange       ierr = ISLocalToGlobalMappingGetIndices(ltogCells, (const PetscInt**)&cells_arr);CHKERRQ(ierr);
117*8cfe4c1fSMichael Lange       /* Convert to positive global cell numbers */
118*8cfe4c1fSMichael Lange       for (n=0; n<size; n++) {if (cells_arr[n] < 0) cells_arr[n] = -(cells_arr[n]+1);}
119*8cfe4c1fSMichael Lange       ierr = ISLocalToGlobalMappingRestoreIndices(ltogCells, (const PetscInt**)&cells_arr);CHKERRQ(ierr);
120*8cfe4c1fSMichael Lange       ierr = ISLocalToGlobalMappingApplyBlock(ltogCells, vOffsets[*numVertices], *adjacency, *adjacency);CHKERRQ(ierr);
121*8cfe4c1fSMichael Lange       ierr = ISLocalToGlobalMappingDestroy(&ltogCells);CHKERRQ(ierr);
122*8cfe4c1fSMichael Lange     }
123*8cfe4c1fSMichael Lange   }
124532c4e7dSMichael Lange   /* Clean up */
125532c4e7dSMichael Lange   ierr = PetscSegBufferDestroy(&adjBuffer);CHKERRQ(ierr);
126532c4e7dSMichael Lange   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
127532c4e7dSMichael Lange   ierr = PetscFree(adj);CHKERRQ(ierr);
128532c4e7dSMichael Lange   PetscFunctionReturn(0);
129532c4e7dSMichael Lange }
130532c4e7dSMichael Lange 
131532c4e7dSMichael Lange #undef __FUNCT__
13270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateNeighborCSR"
13370034214SMatthew G. Knepley PetscErrorCode DMPlexCreateNeighborCSR(DM dm, PetscInt cellHeight, PetscInt *numVertices, PetscInt **offsets, PetscInt **adjacency)
13470034214SMatthew G. Knepley {
13570034214SMatthew G. Knepley   const PetscInt maxFaceCases = 30;
13670034214SMatthew G. Knepley   PetscInt       numFaceCases = 0;
13770034214SMatthew G. Knepley   PetscInt       numFaceVertices[30]; /* maxFaceCases, C89 sucks sucks sucks */
13870034214SMatthew G. Knepley   PetscInt      *off, *adj;
13970034214SMatthew G. Knepley   PetscInt      *neighborCells = NULL;
14070034214SMatthew G. Knepley   PetscInt       dim, cellDim, depth = 0, faceDepth, cStart, cEnd, c, numCells, cell;
14170034214SMatthew G. Knepley   PetscErrorCode ierr;
14270034214SMatthew G. Knepley 
14370034214SMatthew G. Knepley   PetscFunctionBegin;
14470034214SMatthew G. Knepley   /* For parallel partitioning, I think you have to communicate supports */
145c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
14670034214SMatthew G. Knepley   cellDim = dim - cellHeight;
14770034214SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
14870034214SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
14970034214SMatthew G. Knepley   if (cEnd - cStart == 0) {
15070034214SMatthew G. Knepley     if (numVertices) *numVertices = 0;
15170034214SMatthew G. Knepley     if (offsets)   *offsets   = NULL;
15270034214SMatthew G. Knepley     if (adjacency) *adjacency = NULL;
15370034214SMatthew G. Knepley     PetscFunctionReturn(0);
15470034214SMatthew G. Knepley   }
15570034214SMatthew G. Knepley   numCells  = cEnd - cStart;
15670034214SMatthew G. Knepley   faceDepth = depth - cellHeight;
15770034214SMatthew G. Knepley   if (dim == depth) {
15870034214SMatthew G. Knepley     PetscInt f, fStart, fEnd;
15970034214SMatthew G. Knepley 
16070034214SMatthew G. Knepley     ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr);
16170034214SMatthew G. Knepley     /* Count neighboring cells */
16270034214SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, cellHeight+1, &fStart, &fEnd);CHKERRQ(ierr);
16370034214SMatthew G. Knepley     for (f = fStart; f < fEnd; ++f) {
16470034214SMatthew G. Knepley       const PetscInt *support;
16570034214SMatthew G. Knepley       PetscInt        supportSize;
16670034214SMatthew G. Knepley       ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr);
16770034214SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr);
16870034214SMatthew G. Knepley       if (supportSize == 2) {
1699ffc88e4SToby Isaac         PetscInt numChildren;
1709ffc88e4SToby Isaac 
1719ffc88e4SToby Isaac         ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr);
1729ffc88e4SToby Isaac         if (!numChildren) {
17370034214SMatthew G. Knepley           ++off[support[0]-cStart+1];
17470034214SMatthew G. Knepley           ++off[support[1]-cStart+1];
17570034214SMatthew G. Knepley         }
17670034214SMatthew G. Knepley       }
1779ffc88e4SToby Isaac     }
17870034214SMatthew G. Knepley     /* Prefix sum */
17970034214SMatthew G. Knepley     for (c = 1; c <= numCells; ++c) off[c] += off[c-1];
18070034214SMatthew G. Knepley     if (adjacency) {
18170034214SMatthew G. Knepley       PetscInt *tmp;
18270034214SMatthew G. Knepley 
18370034214SMatthew G. Knepley       ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr);
18470034214SMatthew G. Knepley       ierr = PetscMalloc1((numCells+1), &tmp);CHKERRQ(ierr);
18570034214SMatthew G. Knepley       ierr = PetscMemcpy(tmp, off, (numCells+1) * sizeof(PetscInt));CHKERRQ(ierr);
18670034214SMatthew G. Knepley       /* Get neighboring cells */
18770034214SMatthew G. Knepley       for (f = fStart; f < fEnd; ++f) {
18870034214SMatthew G. Knepley         const PetscInt *support;
18970034214SMatthew G. Knepley         PetscInt        supportSize;
19070034214SMatthew G. Knepley         ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr);
19170034214SMatthew G. Knepley         ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr);
19270034214SMatthew G. Knepley         if (supportSize == 2) {
1939ffc88e4SToby Isaac           PetscInt numChildren;
1949ffc88e4SToby Isaac 
1959ffc88e4SToby Isaac           ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr);
1969ffc88e4SToby Isaac           if (!numChildren) {
19770034214SMatthew G. Knepley             adj[tmp[support[0]-cStart]++] = support[1];
19870034214SMatthew G. Knepley             adj[tmp[support[1]-cStart]++] = support[0];
19970034214SMatthew G. Knepley           }
20070034214SMatthew G. Knepley         }
2019ffc88e4SToby Isaac       }
20270034214SMatthew G. Knepley #if defined(PETSC_USE_DEBUG)
20370034214SMatthew 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);
20470034214SMatthew G. Knepley #endif
20570034214SMatthew G. Knepley       ierr = PetscFree(tmp);CHKERRQ(ierr);
20670034214SMatthew G. Knepley     }
20770034214SMatthew G. Knepley     if (numVertices) *numVertices = numCells;
20870034214SMatthew G. Knepley     if (offsets)   *offsets   = off;
20970034214SMatthew G. Knepley     if (adjacency) *adjacency = adj;
21070034214SMatthew G. Knepley     PetscFunctionReturn(0);
21170034214SMatthew G. Knepley   }
21270034214SMatthew G. Knepley   /* Setup face recognition */
21370034214SMatthew G. Knepley   if (faceDepth == 1) {
21470034214SMatthew 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 */
21570034214SMatthew G. Knepley 
21670034214SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
21770034214SMatthew G. Knepley       PetscInt corners;
21870034214SMatthew G. Knepley 
21970034214SMatthew G. Knepley       ierr = DMPlexGetConeSize(dm, c, &corners);CHKERRQ(ierr);
22070034214SMatthew G. Knepley       if (!cornersSeen[corners]) {
22170034214SMatthew G. Knepley         PetscInt nFV;
22270034214SMatthew G. Knepley 
22370034214SMatthew G. Knepley         if (numFaceCases >= maxFaceCases) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Exceeded maximum number of face recognition cases");
22470034214SMatthew G. Knepley         cornersSeen[corners] = 1;
22570034214SMatthew G. Knepley 
22670034214SMatthew G. Knepley         ierr = DMPlexGetNumFaceVertices(dm, cellDim, corners, &nFV);CHKERRQ(ierr);
22770034214SMatthew G. Knepley 
22870034214SMatthew G. Knepley         numFaceVertices[numFaceCases++] = nFV;
22970034214SMatthew G. Knepley       }
23070034214SMatthew G. Knepley     }
23170034214SMatthew G. Knepley   }
23270034214SMatthew G. Knepley   ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr);
23370034214SMatthew G. Knepley   /* Count neighboring cells */
23470034214SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) {
23570034214SMatthew G. Knepley     PetscInt numNeighbors = PETSC_DETERMINE, n;
23670034214SMatthew G. Knepley 
2378b0b4c70SToby Isaac     ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &numNeighbors, &neighborCells);CHKERRQ(ierr);
23870034214SMatthew G. Knepley     /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */
23970034214SMatthew G. Knepley     for (n = 0; n < numNeighbors; ++n) {
24070034214SMatthew G. Knepley       PetscInt        cellPair[2];
24170034214SMatthew G. Knepley       PetscBool       found    = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE;
24270034214SMatthew G. Knepley       PetscInt        meetSize = 0;
24370034214SMatthew G. Knepley       const PetscInt *meet    = NULL;
24470034214SMatthew G. Knepley 
24570034214SMatthew G. Knepley       cellPair[0] = cell; cellPair[1] = neighborCells[n];
24670034214SMatthew G. Knepley       if (cellPair[0] == cellPair[1]) continue;
24770034214SMatthew G. Knepley       if (!found) {
24870034214SMatthew G. Knepley         ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
24970034214SMatthew G. Knepley         if (meetSize) {
25070034214SMatthew G. Knepley           PetscInt f;
25170034214SMatthew G. Knepley 
25270034214SMatthew G. Knepley           for (f = 0; f < numFaceCases; ++f) {
25370034214SMatthew G. Knepley             if (numFaceVertices[f] == meetSize) {
25470034214SMatthew G. Knepley               found = PETSC_TRUE;
25570034214SMatthew G. Knepley               break;
25670034214SMatthew G. Knepley             }
25770034214SMatthew G. Knepley           }
25870034214SMatthew G. Knepley         }
25970034214SMatthew G. Knepley         ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
26070034214SMatthew G. Knepley       }
26170034214SMatthew G. Knepley       if (found) ++off[cell-cStart+1];
26270034214SMatthew G. Knepley     }
26370034214SMatthew G. Knepley   }
26470034214SMatthew G. Knepley   /* Prefix sum */
26570034214SMatthew G. Knepley   for (cell = 1; cell <= numCells; ++cell) off[cell] += off[cell-1];
26670034214SMatthew G. Knepley 
26770034214SMatthew G. Knepley   if (adjacency) {
26870034214SMatthew G. Knepley     ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr);
26970034214SMatthew G. Knepley     /* Get neighboring cells */
27070034214SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
27170034214SMatthew G. Knepley       PetscInt numNeighbors = PETSC_DETERMINE, n;
27270034214SMatthew G. Knepley       PetscInt cellOffset   = 0;
27370034214SMatthew G. Knepley 
2748b0b4c70SToby Isaac       ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &numNeighbors, &neighborCells);CHKERRQ(ierr);
27570034214SMatthew G. Knepley       /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */
27670034214SMatthew G. Knepley       for (n = 0; n < numNeighbors; ++n) {
27770034214SMatthew G. Knepley         PetscInt        cellPair[2];
27870034214SMatthew G. Knepley         PetscBool       found    = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE;
27970034214SMatthew G. Knepley         PetscInt        meetSize = 0;
28070034214SMatthew G. Knepley         const PetscInt *meet    = NULL;
28170034214SMatthew G. Knepley 
28270034214SMatthew G. Knepley         cellPair[0] = cell; cellPair[1] = neighborCells[n];
28370034214SMatthew G. Knepley         if (cellPair[0] == cellPair[1]) continue;
28470034214SMatthew G. Knepley         if (!found) {
28570034214SMatthew G. Knepley           ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
28670034214SMatthew G. Knepley           if (meetSize) {
28770034214SMatthew G. Knepley             PetscInt f;
28870034214SMatthew G. Knepley 
28970034214SMatthew G. Knepley             for (f = 0; f < numFaceCases; ++f) {
29070034214SMatthew G. Knepley               if (numFaceVertices[f] == meetSize) {
29170034214SMatthew G. Knepley                 found = PETSC_TRUE;
29270034214SMatthew G. Knepley                 break;
29370034214SMatthew G. Knepley               }
29470034214SMatthew G. Knepley             }
29570034214SMatthew G. Knepley           }
29670034214SMatthew G. Knepley           ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
29770034214SMatthew G. Knepley         }
29870034214SMatthew G. Knepley         if (found) {
29970034214SMatthew G. Knepley           adj[off[cell-cStart]+cellOffset] = neighborCells[n];
30070034214SMatthew G. Knepley           ++cellOffset;
30170034214SMatthew G. Knepley         }
30270034214SMatthew G. Knepley       }
30370034214SMatthew G. Knepley     }
30470034214SMatthew G. Knepley   }
30570034214SMatthew G. Knepley   ierr = PetscFree(neighborCells);CHKERRQ(ierr);
30670034214SMatthew G. Knepley   if (numVertices) *numVertices = numCells;
30770034214SMatthew G. Knepley   if (offsets)   *offsets   = off;
30870034214SMatthew G. Knepley   if (adjacency) *adjacency = adj;
30970034214SMatthew G. Knepley   PetscFunctionReturn(0);
31070034214SMatthew G. Knepley }
31170034214SMatthew G. Knepley 
31277623264SMatthew G. Knepley #undef __FUNCT__
31377623264SMatthew G. Knepley #define __FUNCT__ "DMPlexEnlargePartition"
31477623264SMatthew G. Knepley /* Expand the partition by BFS on the adjacency graph */
31577623264SMatthew G. Knepley PetscErrorCode DMPlexEnlargePartition(DM dm, const PetscInt start[], const PetscInt adjacency[], PetscSection origPartSection, IS origPartition, PetscSection partSection, IS *partition)
31677623264SMatthew G. Knepley {
31777623264SMatthew G. Knepley   PetscHashI      h;
31877623264SMatthew G. Knepley   const PetscInt *points;
31977623264SMatthew G. Knepley   PetscInt      **tmpPoints, *newPoints, totPoints = 0;
32077623264SMatthew G. Knepley   PetscInt        pStart, pEnd, part, q;
32177623264SMatthew G. Knepley   PetscBool       useCone;
32277623264SMatthew G. Knepley   PetscErrorCode  ierr;
32377623264SMatthew G. Knepley 
32477623264SMatthew G. Knepley   PetscFunctionBegin;
32577623264SMatthew G. Knepley   PetscHashICreate(h);
32677623264SMatthew G. Knepley   ierr = PetscSectionGetChart(origPartSection, &pStart, &pEnd);CHKERRQ(ierr);
32777623264SMatthew G. Knepley   ierr = PetscSectionSetChart(partSection, pStart, pEnd);CHKERRQ(ierr);
32877623264SMatthew G. Knepley   ierr = ISGetIndices(origPartition, &points);CHKERRQ(ierr);
32977623264SMatthew G. Knepley   ierr = PetscMalloc1((pEnd - pStart), &tmpPoints);CHKERRQ(ierr);
33077623264SMatthew G. Knepley   ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
33177623264SMatthew G. Knepley   ierr = DMPlexSetAdjacencyUseCone(dm, PETSC_TRUE);CHKERRQ(ierr);
33277623264SMatthew G. Knepley   for (part = pStart; part < pEnd; ++part) {
33377623264SMatthew G. Knepley     PetscInt *adj = NULL;
33477623264SMatthew G. Knepley     PetscInt  numPoints, nP, numNewPoints, off, p, n = 0;
33577623264SMatthew G. Knepley 
33677623264SMatthew G. Knepley     PetscHashIClear(h);
33777623264SMatthew G. Knepley     ierr = PetscSectionGetDof(origPartSection, part, &numPoints);CHKERRQ(ierr);
33877623264SMatthew G. Knepley     ierr = PetscSectionGetOffset(origPartSection, part, &off);CHKERRQ(ierr);
33977623264SMatthew G. Knepley     /* Add all existing points to h */
34077623264SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
34177623264SMatthew G. Knepley       const PetscInt point = points[off+p];
34277623264SMatthew G. Knepley       PetscHashIAdd(h, point, 1);
34377623264SMatthew G. Knepley     }
34477623264SMatthew G. Knepley     PetscHashISize(h, nP);
34577623264SMatthew G. Knepley     if (nP != numPoints) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid partition has %d points, but only %d were unique", numPoints, nP);
34677623264SMatthew G. Knepley     /* Add all points in next BFS level */
34777623264SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
34877623264SMatthew G. Knepley       const PetscInt point   = points[off+p];
34977623264SMatthew G. Knepley       PetscInt       adjSize = PETSC_DETERMINE, a;
35077623264SMatthew G. Knepley 
35177623264SMatthew G. Knepley       ierr = DMPlexGetAdjacency(dm, point, &adjSize, &adj);CHKERRQ(ierr);
35277623264SMatthew G. Knepley       for (a = 0; a < adjSize; ++a) PetscHashIAdd(h, adj[a], 1);
35377623264SMatthew G. Knepley     }
35477623264SMatthew G. Knepley     PetscHashISize(h, numNewPoints);
35577623264SMatthew G. Knepley     ierr = PetscSectionSetDof(partSection, part, numNewPoints);CHKERRQ(ierr);
35677623264SMatthew G. Knepley     ierr = PetscMalloc1(numNewPoints, &tmpPoints[part]);CHKERRQ(ierr);
35777623264SMatthew G. Knepley     ierr = PetscHashIGetKeys(h, &n, tmpPoints[part]);CHKERRQ(ierr);
35877623264SMatthew G. Knepley     ierr = PetscFree(adj);CHKERRQ(ierr);
35977623264SMatthew G. Knepley     totPoints += numNewPoints;
36077623264SMatthew G. Knepley   }
36177623264SMatthew G. Knepley   ierr = DMPlexSetAdjacencyUseCone(dm, useCone);CHKERRQ(ierr);
36277623264SMatthew G. Knepley   ierr = ISRestoreIndices(origPartition, &points);CHKERRQ(ierr);
36377623264SMatthew G. Knepley   PetscHashIDestroy(h);
36477623264SMatthew G. Knepley   ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
36577623264SMatthew G. Knepley   ierr = PetscMalloc1(totPoints, &newPoints);CHKERRQ(ierr);
36677623264SMatthew G. Knepley   for (part = pStart, q = 0; part < pEnd; ++part) {
36777623264SMatthew G. Knepley     PetscInt numPoints, p;
36877623264SMatthew G. Knepley 
36977623264SMatthew G. Knepley     ierr = PetscSectionGetDof(partSection, part, &numPoints);CHKERRQ(ierr);
37077623264SMatthew G. Knepley     for (p = 0; p < numPoints; ++p, ++q) newPoints[q] = tmpPoints[part][p];
37177623264SMatthew G. Knepley     ierr = PetscFree(tmpPoints[part]);CHKERRQ(ierr);
37277623264SMatthew G. Knepley   }
37377623264SMatthew G. Knepley   ierr = PetscFree(tmpPoints);CHKERRQ(ierr);
37477623264SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), totPoints, newPoints, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
37577623264SMatthew G. Knepley   PetscFunctionReturn(0);
37677623264SMatthew G. Knepley }
37777623264SMatthew G. Knepley 
37877623264SMatthew G. Knepley #undef __FUNCT__
37977623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerRegister"
38077623264SMatthew G. Knepley /*@C
38177623264SMatthew G. Knepley   PetscPartitionerRegister - Adds a new PetscPartitioner implementation
38277623264SMatthew G. Knepley 
38377623264SMatthew G. Knepley   Not Collective
38477623264SMatthew G. Knepley 
38577623264SMatthew G. Knepley   Input Parameters:
38677623264SMatthew G. Knepley + name        - The name of a new user-defined creation routine
38777623264SMatthew G. Knepley - create_func - The creation routine itself
38877623264SMatthew G. Knepley 
38977623264SMatthew G. Knepley   Notes:
39077623264SMatthew G. Knepley   PetscPartitionerRegister() may be called multiple times to add several user-defined PetscPartitioners
39177623264SMatthew G. Knepley 
39277623264SMatthew G. Knepley   Sample usage:
39377623264SMatthew G. Knepley .vb
39477623264SMatthew G. Knepley     PetscPartitionerRegister("my_part", MyPetscPartitionerCreate);
39577623264SMatthew G. Knepley .ve
39677623264SMatthew G. Knepley 
39777623264SMatthew G. Knepley   Then, your PetscPartitioner type can be chosen with the procedural interface via
39877623264SMatthew G. Knepley .vb
39977623264SMatthew G. Knepley     PetscPartitionerCreate(MPI_Comm, PetscPartitioner *);
40077623264SMatthew G. Knepley     PetscPartitionerSetType(PetscPartitioner, "my_part");
40177623264SMatthew G. Knepley .ve
40277623264SMatthew G. Knepley    or at runtime via the option
40377623264SMatthew G. Knepley .vb
40477623264SMatthew G. Knepley     -petscpartitioner_type my_part
40577623264SMatthew G. Knepley .ve
40677623264SMatthew G. Knepley 
40777623264SMatthew G. Knepley   Level: advanced
40877623264SMatthew G. Knepley 
40977623264SMatthew G. Knepley .keywords: PetscPartitioner, register
41077623264SMatthew G. Knepley .seealso: PetscPartitionerRegisterAll(), PetscPartitionerRegisterDestroy()
41177623264SMatthew G. Knepley 
41277623264SMatthew G. Knepley @*/
41377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerRegister(const char sname[], PetscErrorCode (*function)(PetscPartitioner))
41477623264SMatthew G. Knepley {
41577623264SMatthew G. Knepley   PetscErrorCode ierr;
41677623264SMatthew G. Knepley 
41777623264SMatthew G. Knepley   PetscFunctionBegin;
41877623264SMatthew G. Knepley   ierr = PetscFunctionListAdd(&PetscPartitionerList, sname, function);CHKERRQ(ierr);
41977623264SMatthew G. Knepley   PetscFunctionReturn(0);
42077623264SMatthew G. Knepley }
42177623264SMatthew G. Knepley 
42277623264SMatthew G. Knepley #undef __FUNCT__
42377623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerSetType"
42477623264SMatthew G. Knepley /*@C
42577623264SMatthew G. Knepley   PetscPartitionerSetType - Builds a particular PetscPartitioner
42677623264SMatthew G. Knepley 
42777623264SMatthew G. Knepley   Collective on PetscPartitioner
42877623264SMatthew G. Knepley 
42977623264SMatthew G. Knepley   Input Parameters:
43077623264SMatthew G. Knepley + part - The PetscPartitioner object
43177623264SMatthew G. Knepley - name - The kind of partitioner
43277623264SMatthew G. Knepley 
43377623264SMatthew G. Knepley   Options Database Key:
43477623264SMatthew G. Knepley . -petscpartitioner_type <type> - Sets the PetscPartitioner type; use -help for a list of available types
43577623264SMatthew G. Knepley 
43677623264SMatthew G. Knepley   Level: intermediate
43777623264SMatthew G. Knepley 
43877623264SMatthew G. Knepley .keywords: PetscPartitioner, set, type
43977623264SMatthew G. Knepley .seealso: PetscPartitionerGetType(), PetscPartitionerCreate()
44077623264SMatthew G. Knepley @*/
44177623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetType(PetscPartitioner part, PetscPartitionerType name)
44277623264SMatthew G. Knepley {
44377623264SMatthew G. Knepley   PetscErrorCode (*r)(PetscPartitioner);
44477623264SMatthew G. Knepley   PetscBool      match;
44577623264SMatthew G. Knepley   PetscErrorCode ierr;
44677623264SMatthew G. Knepley 
44777623264SMatthew G. Knepley   PetscFunctionBegin;
44877623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
44977623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) part, name, &match);CHKERRQ(ierr);
45077623264SMatthew G. Knepley   if (match) PetscFunctionReturn(0);
45177623264SMatthew G. Knepley 
45277623264SMatthew G. Knepley   if (!PetscPartitionerRegisterAllCalled) {ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr);}
45377623264SMatthew G. Knepley   ierr = PetscFunctionListFind(PetscPartitionerList, name, &r);CHKERRQ(ierr);
45477623264SMatthew G. Knepley   if (!r) SETERRQ1(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscPartitioner type: %s", name);
45577623264SMatthew G. Knepley 
45677623264SMatthew G. Knepley   if (part->ops->destroy) {
45777623264SMatthew G. Knepley     ierr              = (*part->ops->destroy)(part);CHKERRQ(ierr);
45877623264SMatthew G. Knepley     part->ops->destroy = NULL;
45977623264SMatthew G. Knepley   }
46077623264SMatthew G. Knepley   ierr = (*r)(part);CHKERRQ(ierr);
46177623264SMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) part, name);CHKERRQ(ierr);
46277623264SMatthew G. Knepley   PetscFunctionReturn(0);
46377623264SMatthew G. Knepley }
46477623264SMatthew G. Knepley 
46577623264SMatthew G. Knepley #undef __FUNCT__
46677623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerGetType"
46777623264SMatthew G. Knepley /*@C
46877623264SMatthew G. Knepley   PetscPartitionerGetType - Gets the PetscPartitioner type name (as a string) from the object.
46977623264SMatthew G. Knepley 
47077623264SMatthew G. Knepley   Not Collective
47177623264SMatthew G. Knepley 
47277623264SMatthew G. Knepley   Input Parameter:
47377623264SMatthew G. Knepley . part - The PetscPartitioner
47477623264SMatthew G. Knepley 
47577623264SMatthew G. Knepley   Output Parameter:
47677623264SMatthew G. Knepley . name - The PetscPartitioner type name
47777623264SMatthew G. Knepley 
47877623264SMatthew G. Knepley   Level: intermediate
47977623264SMatthew G. Knepley 
48077623264SMatthew G. Knepley .keywords: PetscPartitioner, get, type, name
48177623264SMatthew G. Knepley .seealso: PetscPartitionerSetType(), PetscPartitionerCreate()
48277623264SMatthew G. Knepley @*/
48377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerGetType(PetscPartitioner part, PetscPartitionerType *name)
48477623264SMatthew G. Knepley {
48577623264SMatthew G. Knepley   PetscErrorCode ierr;
48677623264SMatthew G. Knepley 
48777623264SMatthew G. Knepley   PetscFunctionBegin;
48877623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
48977623264SMatthew G. Knepley   PetscValidPointer(name, 2);
49077623264SMatthew G. Knepley   if (!PetscPartitionerRegisterAllCalled) {ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr);}
49177623264SMatthew G. Knepley   *name = ((PetscObject) part)->type_name;
49277623264SMatthew G. Knepley   PetscFunctionReturn(0);
49377623264SMatthew G. Knepley }
49477623264SMatthew G. Knepley 
49577623264SMatthew G. Knepley #undef __FUNCT__
49677623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView"
49777623264SMatthew G. Knepley /*@C
49877623264SMatthew G. Knepley   PetscPartitionerView - Views a PetscPartitioner
49977623264SMatthew G. Knepley 
50077623264SMatthew G. Knepley   Collective on PetscPartitioner
50177623264SMatthew G. Knepley 
50277623264SMatthew G. Knepley   Input Parameter:
50377623264SMatthew G. Knepley + part - the PetscPartitioner object to view
50477623264SMatthew G. Knepley - v    - the viewer
50577623264SMatthew G. Knepley 
50677623264SMatthew G. Knepley   Level: developer
50777623264SMatthew G. Knepley 
50877623264SMatthew G. Knepley .seealso: PetscPartitionerDestroy()
50977623264SMatthew G. Knepley @*/
51077623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView(PetscPartitioner part, PetscViewer v)
51177623264SMatthew G. Knepley {
51277623264SMatthew G. Knepley   PetscErrorCode ierr;
51377623264SMatthew G. Knepley 
51477623264SMatthew G. Knepley   PetscFunctionBegin;
51577623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
51677623264SMatthew G. Knepley   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) part), &v);CHKERRQ(ierr);}
51777623264SMatthew G. Knepley   if (part->ops->view) {ierr = (*part->ops->view)(part, v);CHKERRQ(ierr);}
51877623264SMatthew G. Knepley   PetscFunctionReturn(0);
51977623264SMatthew G. Knepley }
52077623264SMatthew G. Knepley 
52177623264SMatthew G. Knepley #undef __FUNCT__
52277623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerViewFromOptions"
52377623264SMatthew G. Knepley /*
52477623264SMatthew G. Knepley   PetscPartitionerViewFromOptions - Processes command line options to determine if/how a PetscPartitioner is to be viewed.
52577623264SMatthew G. Knepley 
52677623264SMatthew G. Knepley   Collective on PetscPartitioner
52777623264SMatthew G. Knepley 
52877623264SMatthew G. Knepley   Input Parameters:
52977623264SMatthew G. Knepley + part   - the PetscPartitioner
53077623264SMatthew G. Knepley . prefix - prefix to use for viewing, or NULL
53177623264SMatthew G. Knepley - optionname - option to activate viewing
53277623264SMatthew G. Knepley 
53377623264SMatthew G. Knepley   Level: intermediate
53477623264SMatthew G. Knepley 
53577623264SMatthew G. Knepley .keywords: PetscPartitioner, view, options, database
53677623264SMatthew G. Knepley .seealso: VecViewFromOptions(), MatViewFromOptions()
53777623264SMatthew G. Knepley */
53877623264SMatthew G. Knepley PetscErrorCode PetscPartitionerViewFromOptions(PetscPartitioner part, const char prefix[], const char optionname[])
53977623264SMatthew G. Knepley {
54077623264SMatthew G. Knepley   PetscViewer       viewer;
54177623264SMatthew G. Knepley   PetscViewerFormat format;
54277623264SMatthew G. Knepley   PetscBool         flg;
54377623264SMatthew G. Knepley   PetscErrorCode    ierr;
54477623264SMatthew G. Knepley 
54577623264SMatthew G. Knepley   PetscFunctionBegin;
54677623264SMatthew G. Knepley   if (prefix) {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) part), prefix,                       optionname, &viewer, &format, &flg);CHKERRQ(ierr);}
54777623264SMatthew G. Knepley   else        {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) part), ((PetscObject) part)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);}
54877623264SMatthew G. Knepley   if (flg) {
54977623264SMatthew G. Knepley     ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr);
55077623264SMatthew G. Knepley     ierr = PetscPartitionerView(part, viewer);CHKERRQ(ierr);
55177623264SMatthew G. Knepley     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
55277623264SMatthew G. Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
55377623264SMatthew G. Knepley   }
55477623264SMatthew G. Knepley   PetscFunctionReturn(0);
55577623264SMatthew G. Knepley }
55677623264SMatthew G. Knepley #undef __FUNCT__
55777623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerSetTypeFromOptions_Internal"
55877623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetTypeFromOptions_Internal(PetscPartitioner part)
55977623264SMatthew G. Knepley {
56077623264SMatthew G. Knepley   const char    *defaultType;
56177623264SMatthew G. Knepley   char           name[256];
56277623264SMatthew G. Knepley   PetscBool      flg;
56377623264SMatthew G. Knepley   PetscErrorCode ierr;
56477623264SMatthew G. Knepley 
56577623264SMatthew G. Knepley   PetscFunctionBegin;
56677623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
56777623264SMatthew G. Knepley   if (!((PetscObject) part)->type_name) defaultType = PETSCPARTITIONERCHACO;
56877623264SMatthew G. Knepley   else                                  defaultType = ((PetscObject) part)->type_name;
56977623264SMatthew G. Knepley   if (!PetscPartitionerRegisterAllCalled) {ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr);}
57077623264SMatthew G. Knepley 
57177623264SMatthew G. Knepley   ierr = PetscObjectOptionsBegin((PetscObject) part);CHKERRQ(ierr);
57277623264SMatthew G. Knepley   ierr = PetscOptionsFList("-petscpartitioner_type", "Graph partitioner", "PetscPartitionerSetType", PetscPartitionerList, defaultType, name, 256, &flg);CHKERRQ(ierr);
57377623264SMatthew G. Knepley   if (flg) {
57477623264SMatthew G. Knepley     ierr = PetscPartitionerSetType(part, name);CHKERRQ(ierr);
57577623264SMatthew G. Knepley   } else if (!((PetscObject) part)->type_name) {
57677623264SMatthew G. Knepley     ierr = PetscPartitionerSetType(part, defaultType);CHKERRQ(ierr);
57777623264SMatthew G. Knepley   }
57877623264SMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
57977623264SMatthew G. Knepley   PetscFunctionReturn(0);
58077623264SMatthew G. Knepley }
58177623264SMatthew G. Knepley 
58277623264SMatthew G. Knepley #undef __FUNCT__
58377623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerSetFromOptions"
58477623264SMatthew G. Knepley /*@
58577623264SMatthew G. Knepley   PetscPartitionerSetFromOptions - sets parameters in a PetscPartitioner from the options database
58677623264SMatthew G. Knepley 
58777623264SMatthew G. Knepley   Collective on PetscPartitioner
58877623264SMatthew G. Knepley 
58977623264SMatthew G. Knepley   Input Parameter:
59077623264SMatthew G. Knepley . part - the PetscPartitioner object to set options for
59177623264SMatthew G. Knepley 
59277623264SMatthew G. Knepley   Level: developer
59377623264SMatthew G. Knepley 
59477623264SMatthew G. Knepley .seealso: PetscPartitionerView()
59577623264SMatthew G. Knepley @*/
59677623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetFromOptions(PetscPartitioner part)
59777623264SMatthew G. Knepley {
59877623264SMatthew G. Knepley   PetscErrorCode ierr;
59977623264SMatthew G. Knepley 
60077623264SMatthew G. Knepley   PetscFunctionBegin;
60177623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
60277623264SMatthew G. Knepley   ierr = PetscPartitionerSetTypeFromOptions_Internal(part);CHKERRQ(ierr);
60377623264SMatthew G. Knepley 
60477623264SMatthew G. Knepley   ierr = PetscObjectOptionsBegin((PetscObject) part);CHKERRQ(ierr);
60577623264SMatthew G. Knepley   if (part->ops->setfromoptions) {ierr = (*part->ops->setfromoptions)(part);CHKERRQ(ierr);}
60677623264SMatthew G. Knepley   /* process any options handlers added with PetscObjectAddOptionsHandler() */
60777623264SMatthew G. Knepley   ierr = PetscObjectProcessOptionsHandlers((PetscObject) part);CHKERRQ(ierr);
60877623264SMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
60977623264SMatthew G. Knepley   ierr = PetscPartitionerViewFromOptions(part, NULL, "-petscpartitioner_view");CHKERRQ(ierr);
61077623264SMatthew G. Knepley   PetscFunctionReturn(0);
61177623264SMatthew G. Knepley }
61277623264SMatthew G. Knepley 
61377623264SMatthew G. Knepley #undef __FUNCT__
61477623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerSetUp"
61577623264SMatthew G. Knepley /*@C
61677623264SMatthew G. Knepley   PetscPartitionerSetUp - Construct data structures for the PetscPartitioner
61777623264SMatthew G. Knepley 
61877623264SMatthew G. Knepley   Collective on PetscPartitioner
61977623264SMatthew G. Knepley 
62077623264SMatthew G. Knepley   Input Parameter:
62177623264SMatthew G. Knepley . part - the PetscPartitioner object to setup
62277623264SMatthew G. Knepley 
62377623264SMatthew G. Knepley   Level: developer
62477623264SMatthew G. Knepley 
62577623264SMatthew G. Knepley .seealso: PetscPartitionerView(), PetscPartitionerDestroy()
62677623264SMatthew G. Knepley @*/
62777623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetUp(PetscPartitioner part)
62877623264SMatthew G. Knepley {
62977623264SMatthew G. Knepley   PetscErrorCode ierr;
63077623264SMatthew G. Knepley 
63177623264SMatthew G. Knepley   PetscFunctionBegin;
63277623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
63377623264SMatthew G. Knepley   if (part->ops->setup) {ierr = (*part->ops->setup)(part);CHKERRQ(ierr);}
63477623264SMatthew G. Knepley   PetscFunctionReturn(0);
63577623264SMatthew G. Knepley }
63677623264SMatthew G. Knepley 
63777623264SMatthew G. Knepley #undef __FUNCT__
63877623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerDestroy"
63977623264SMatthew G. Knepley /*@
64077623264SMatthew G. Knepley   PetscPartitionerDestroy - Destroys a PetscPartitioner object
64177623264SMatthew G. Knepley 
64277623264SMatthew G. Knepley   Collective on PetscPartitioner
64377623264SMatthew G. Knepley 
64477623264SMatthew G. Knepley   Input Parameter:
64577623264SMatthew G. Knepley . part - the PetscPartitioner object to destroy
64677623264SMatthew G. Knepley 
64777623264SMatthew G. Knepley   Level: developer
64877623264SMatthew G. Knepley 
64977623264SMatthew G. Knepley .seealso: PetscPartitionerView()
65077623264SMatthew G. Knepley @*/
65177623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy(PetscPartitioner *part)
65277623264SMatthew G. Knepley {
65377623264SMatthew G. Knepley   PetscErrorCode ierr;
65477623264SMatthew G. Knepley 
65577623264SMatthew G. Knepley   PetscFunctionBegin;
65677623264SMatthew G. Knepley   if (!*part) PetscFunctionReturn(0);
65777623264SMatthew G. Knepley   PetscValidHeaderSpecific((*part), PETSCPARTITIONER_CLASSID, 1);
65877623264SMatthew G. Knepley 
65977623264SMatthew G. Knepley   if (--((PetscObject)(*part))->refct > 0) {*part = 0; PetscFunctionReturn(0);}
66077623264SMatthew G. Knepley   ((PetscObject) (*part))->refct = 0;
66177623264SMatthew G. Knepley 
66277623264SMatthew G. Knepley   if ((*part)->ops->destroy) {ierr = (*(*part)->ops->destroy)(*part);CHKERRQ(ierr);}
66377623264SMatthew G. Knepley   ierr = PetscHeaderDestroy(part);CHKERRQ(ierr);
66477623264SMatthew G. Knepley   PetscFunctionReturn(0);
66577623264SMatthew G. Knepley }
66677623264SMatthew G. Knepley 
66777623264SMatthew G. Knepley #undef __FUNCT__
66877623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerCreate"
66977623264SMatthew G. Knepley /*@
67077623264SMatthew G. Knepley   PetscPartitionerCreate - Creates an empty PetscPartitioner object. The type can then be set with PetscPartitionerSetType().
67177623264SMatthew G. Knepley 
67277623264SMatthew G. Knepley   Collective on MPI_Comm
67377623264SMatthew G. Knepley 
67477623264SMatthew G. Knepley   Input Parameter:
67577623264SMatthew G. Knepley . comm - The communicator for the PetscPartitioner object
67677623264SMatthew G. Knepley 
67777623264SMatthew G. Knepley   Output Parameter:
67877623264SMatthew G. Knepley . part - The PetscPartitioner object
67977623264SMatthew G. Knepley 
68077623264SMatthew G. Knepley   Level: beginner
68177623264SMatthew G. Knepley 
68277623264SMatthew G. Knepley .seealso: PetscPartitionerSetType(), PETSCPARTITIONERCHACO, PETSCPARTITIONERPARMETIS, PETSCPARTITIONERSHELL
68377623264SMatthew G. Knepley @*/
68477623264SMatthew G. Knepley PetscErrorCode PetscPartitionerCreate(MPI_Comm comm, PetscPartitioner *part)
68577623264SMatthew G. Knepley {
68677623264SMatthew G. Knepley   PetscPartitioner p;
68777623264SMatthew G. Knepley   PetscErrorCode   ierr;
68877623264SMatthew G. Knepley 
68977623264SMatthew G. Knepley   PetscFunctionBegin;
69077623264SMatthew G. Knepley   PetscValidPointer(part, 2);
69177623264SMatthew G. Knepley   *part = NULL;
69277623264SMatthew G. Knepley   ierr = PetscFVInitializePackage();CHKERRQ(ierr);
69377623264SMatthew G. Knepley 
69477623264SMatthew G. Knepley   ierr = PetscHeaderCreate(p, _p_PetscPartitioner, struct _PetscPartitionerOps, PETSCPARTITIONER_CLASSID, "PetscPartitioner", "Graph Partitioner", "PetscPartitioner", comm, PetscPartitionerDestroy, PetscPartitionerView);CHKERRQ(ierr);
69577623264SMatthew G. Knepley   ierr = PetscMemzero(p->ops, sizeof(struct _PetscPartitionerOps));CHKERRQ(ierr);
69677623264SMatthew G. Knepley 
69777623264SMatthew G. Knepley   *part = p;
69877623264SMatthew G. Knepley   PetscFunctionReturn(0);
69977623264SMatthew G. Knepley }
70077623264SMatthew G. Knepley 
70177623264SMatthew G. Knepley #undef __FUNCT__
70277623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerPartition"
70377623264SMatthew G. Knepley /*@
70477623264SMatthew G. Knepley   PetscPartitionerPartition - Create a non-overlapping partition of the cells in the mesh
70577623264SMatthew G. Knepley 
70677623264SMatthew G. Knepley   Collective on DM
70777623264SMatthew G. Knepley 
70877623264SMatthew G. Knepley   Input Parameters:
70977623264SMatthew G. Knepley + part    - The PetscPartitioner
710f8987ae8SMichael Lange - dm      - The mesh DM
71177623264SMatthew G. Knepley 
71277623264SMatthew G. Knepley   Output Parameters:
71377623264SMatthew G. Knepley + partSection     - The PetscSection giving the division of points by partition
714f8987ae8SMichael Lange - partition       - The list of points by partition
71577623264SMatthew G. Knepley 
71677623264SMatthew G. Knepley   Note: Instead of cells, points at a given height can be partitioned by calling PetscPartitionerSetPointHeight()
71777623264SMatthew G. Knepley 
71877623264SMatthew G. Knepley   Level: developer
71977623264SMatthew G. Knepley 
72077623264SMatthew G. Knepley .seealso DMPlexDistribute(), PetscPartitionerSetPointHeight(), PetscPartitionerCreate()
7214b15ede2SMatthew G. Knepley @*/
722f8987ae8SMichael Lange PetscErrorCode PetscPartitionerPartition(PetscPartitioner part, DM dm, PetscSection partSection, IS *partition)
72377623264SMatthew G. Knepley {
72477623264SMatthew G. Knepley   PetscMPIInt    size;
72577623264SMatthew G. Knepley   PetscErrorCode ierr;
72677623264SMatthew G. Knepley 
72777623264SMatthew G. Knepley   PetscFunctionBegin;
72877623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
72977623264SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
73077623264SMatthew G. Knepley   PetscValidHeaderSpecific(partSection, PETSC_SECTION_CLASSID, 4);
73177623264SMatthew G. Knepley   PetscValidPointer(partition, 5);
73277623264SMatthew G. Knepley   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) part), &size);CHKERRQ(ierr);
73377623264SMatthew G. Knepley   if (size == 1) {
73477623264SMatthew G. Knepley     PetscInt *points;
73577623264SMatthew G. Knepley     PetscInt  cStart, cEnd, c;
73677623264SMatthew G. Knepley 
73777623264SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, part->height, &cStart, &cEnd);CHKERRQ(ierr);
73877623264SMatthew G. Knepley     ierr = PetscSectionSetChart(partSection, 0, size);CHKERRQ(ierr);
73977623264SMatthew G. Knepley     ierr = PetscSectionSetDof(partSection, 0, cEnd-cStart);CHKERRQ(ierr);
74077623264SMatthew G. Knepley     ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
74177623264SMatthew G. Knepley     ierr = PetscMalloc1(cEnd-cStart, &points);CHKERRQ(ierr);
74277623264SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) points[c] = c;
74377623264SMatthew G. Knepley     ierr = ISCreateGeneral(PetscObjectComm((PetscObject) part), cEnd-cStart, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
74477623264SMatthew G. Knepley     PetscFunctionReturn(0);
74577623264SMatthew G. Knepley   }
74677623264SMatthew G. Knepley   if (part->height == 0) {
74777623264SMatthew G. Knepley     PetscInt  numVertices;
74877623264SMatthew G. Knepley     PetscInt *start     = NULL;
74977623264SMatthew G. Knepley     PetscInt *adjacency = NULL;
75077623264SMatthew G. Knepley 
751532c4e7dSMichael Lange     ierr = DMPlexCreatePartitionerGraph(dm, 0, &numVertices, &start, &adjacency);CHKERRQ(ierr);
75277623264SMatthew G. Knepley     if (!part->ops->partition) SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_WRONGSTATE, "PetscPartitioner has no type");
75377623264SMatthew G. Knepley     ierr = (*part->ops->partition)(part, dm, size, numVertices, start, adjacency, partSection, partition);CHKERRQ(ierr);
75477623264SMatthew G. Knepley     ierr = PetscFree(start);CHKERRQ(ierr);
75577623264SMatthew G. Knepley     ierr = PetscFree(adjacency);CHKERRQ(ierr);
75677623264SMatthew G. Knepley   } else SETERRQ1(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_OUTOFRANGE, "Invalid height %D for points to partition", part->height);
75777623264SMatthew G. Knepley   PetscFunctionReturn(0);
75877623264SMatthew G. Knepley }
75977623264SMatthew G. Knepley 
76077623264SMatthew G. Knepley #undef __FUNCT__
76177623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerDestroy_Shell"
76277623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy_Shell(PetscPartitioner part)
76377623264SMatthew G. Knepley {
76477623264SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
76577623264SMatthew G. Knepley   PetscErrorCode          ierr;
76677623264SMatthew G. Knepley 
76777623264SMatthew G. Knepley   PetscFunctionBegin;
76877623264SMatthew G. Knepley   ierr = PetscSectionDestroy(&p->section);CHKERRQ(ierr);
76977623264SMatthew G. Knepley   ierr = ISDestroy(&p->partition);CHKERRQ(ierr);
77077623264SMatthew G. Knepley   ierr = PetscFree(p);CHKERRQ(ierr);
77177623264SMatthew G. Knepley   PetscFunctionReturn(0);
77277623264SMatthew G. Knepley }
77377623264SMatthew G. Knepley 
77477623264SMatthew G. Knepley #undef __FUNCT__
77577623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_Shell_Ascii"
77677623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Shell_Ascii(PetscPartitioner part, PetscViewer viewer)
77777623264SMatthew G. Knepley {
77877623264SMatthew G. Knepley   PetscViewerFormat format;
77977623264SMatthew G. Knepley   PetscErrorCode    ierr;
78077623264SMatthew G. Knepley 
78177623264SMatthew G. Knepley   PetscFunctionBegin;
78277623264SMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
78377623264SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Shell Graph Partitioner:\n");CHKERRQ(ierr);
78477623264SMatthew G. Knepley   PetscFunctionReturn(0);
78577623264SMatthew G. Knepley }
78677623264SMatthew G. Knepley 
78777623264SMatthew G. Knepley #undef __FUNCT__
78877623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_Shell"
78977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Shell(PetscPartitioner part, PetscViewer viewer)
79077623264SMatthew G. Knepley {
79177623264SMatthew G. Knepley   PetscBool      iascii;
79277623264SMatthew G. Knepley   PetscErrorCode ierr;
79377623264SMatthew G. Knepley 
79477623264SMatthew G. Knepley   PetscFunctionBegin;
79577623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
79677623264SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
79777623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
79877623264SMatthew G. Knepley   if (iascii) {ierr = PetscPartitionerView_Shell_Ascii(part, viewer);CHKERRQ(ierr);}
79977623264SMatthew G. Knepley   PetscFunctionReturn(0);
80077623264SMatthew G. Knepley }
80177623264SMatthew G. Knepley 
80277623264SMatthew G. Knepley #undef __FUNCT__
80377623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerPartition_Shell"
80477623264SMatthew G. Knepley PetscErrorCode PetscPartitionerPartition_Shell(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
80577623264SMatthew G. Knepley {
80677623264SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
80777623264SMatthew G. Knepley   PetscInt                np;
80877623264SMatthew G. Knepley   PetscErrorCode          ierr;
80977623264SMatthew G. Knepley 
81077623264SMatthew G. Knepley   PetscFunctionBegin;
81177623264SMatthew G. Knepley   ierr = PetscSectionGetChart(p->section, NULL, &np);CHKERRQ(ierr);
81277623264SMatthew G. Knepley   if (nparts != np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of requested partitions %d != configured partitions %d", nparts, np);
81377623264SMatthew G. Knepley   ierr = ISGetLocalSize(p->partition, &np);CHKERRQ(ierr);
81477623264SMatthew G. Knepley   if (numVertices != np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of input vertices %d != configured vertices %d", numVertices, np);
8155680f57bSMatthew G. Knepley   ierr = PetscSectionCopy(p->section, partSection);CHKERRQ(ierr);
81677623264SMatthew G. Knepley   *partition = p->partition;
81777623264SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) p->partition);CHKERRQ(ierr);
81877623264SMatthew G. Knepley   PetscFunctionReturn(0);
81977623264SMatthew G. Knepley }
82077623264SMatthew G. Knepley 
82177623264SMatthew G. Knepley #undef __FUNCT__
82277623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerInitialize_Shell"
82377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerInitialize_Shell(PetscPartitioner part)
82477623264SMatthew G. Knepley {
82577623264SMatthew G. Knepley   PetscFunctionBegin;
82677623264SMatthew G. Knepley   part->ops->view      = PetscPartitionerView_Shell;
82777623264SMatthew G. Knepley   part->ops->destroy   = PetscPartitionerDestroy_Shell;
82877623264SMatthew G. Knepley   part->ops->partition = PetscPartitionerPartition_Shell;
82977623264SMatthew G. Knepley   PetscFunctionReturn(0);
83077623264SMatthew G. Knepley }
83177623264SMatthew G. Knepley 
83277623264SMatthew G. Knepley /*MC
83377623264SMatthew G. Knepley   PETSCPARTITIONERSHELL = "shell" - A PetscPartitioner object
83477623264SMatthew G. Knepley 
83577623264SMatthew G. Knepley   Level: intermediate
83677623264SMatthew G. Knepley 
83777623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
83877623264SMatthew G. Knepley M*/
83977623264SMatthew G. Knepley 
84077623264SMatthew G. Knepley #undef __FUNCT__
84177623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerCreate_Shell"
84277623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Shell(PetscPartitioner part)
84377623264SMatthew G. Knepley {
84477623264SMatthew G. Knepley   PetscPartitioner_Shell *p;
84577623264SMatthew G. Knepley   PetscErrorCode          ierr;
84677623264SMatthew G. Knepley 
84777623264SMatthew G. Knepley   PetscFunctionBegin;
84877623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
84977623264SMatthew G. Knepley   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
85077623264SMatthew G. Knepley   part->data = p;
85177623264SMatthew G. Knepley 
85277623264SMatthew G. Knepley   ierr = PetscPartitionerInitialize_Shell(part);CHKERRQ(ierr);
85377623264SMatthew G. Knepley   PetscFunctionReturn(0);
85477623264SMatthew G. Knepley }
85577623264SMatthew G. Knepley 
85677623264SMatthew G. Knepley #undef __FUNCT__
8575680f57bSMatthew G. Knepley #define __FUNCT__ "PetscPartitionerShellSetPartition"
8585680f57bSMatthew G. Knepley /*@C
8595680f57bSMatthew G. Knepley   PetscPartitionerShellSetPartition - Set an artifical partition for a mesh
8605680f57bSMatthew G. Knepley 
8615680f57bSMatthew G. Knepley   Collective on DM
8625680f57bSMatthew G. Knepley 
8635680f57bSMatthew G. Knepley   Input Parameters:
8645680f57bSMatthew G. Knepley + part    - The PetscPartitioner
8655680f57bSMatthew G. Knepley . dm      - The mesh DM
8665680f57bSMatthew G. Knepley - enlarge - Expand each partition with neighbors
8675680f57bSMatthew G. Knepley 
8685680f57bSMatthew G. Knepley   Level: developer
8695680f57bSMatthew G. Knepley 
8705680f57bSMatthew G. Knepley .seealso DMPlexDistribute(), PetscPartitionerCreate()
8715680f57bSMatthew G. Knepley @*/
8725680f57bSMatthew G. Knepley PetscErrorCode PetscPartitionerShellSetPartition(PetscPartitioner part, PetscInt numProcs, const PetscInt sizes[], const PetscInt points[])
8735680f57bSMatthew G. Knepley {
8745680f57bSMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
8755680f57bSMatthew G. Knepley   PetscInt                proc, numPoints;
8765680f57bSMatthew G. Knepley   PetscErrorCode          ierr;
8775680f57bSMatthew G. Knepley 
8785680f57bSMatthew G. Knepley   PetscFunctionBegin;
8795680f57bSMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
8805680f57bSMatthew G. Knepley   if (sizes) {PetscValidPointer(sizes, 3);}
8815680f57bSMatthew G. Knepley   if (sizes) {PetscValidPointer(points, 4);}
8825680f57bSMatthew G. Knepley   ierr = PetscSectionDestroy(&p->section);CHKERRQ(ierr);
8835680f57bSMatthew G. Knepley   ierr = ISDestroy(&p->partition);CHKERRQ(ierr);
8845680f57bSMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) part), &p->section);CHKERRQ(ierr);
8855680f57bSMatthew G. Knepley   ierr = PetscSectionSetChart(p->section, 0, numProcs);CHKERRQ(ierr);
8865680f57bSMatthew G. Knepley   if (sizes) {
8875680f57bSMatthew G. Knepley     for (proc = 0; proc < numProcs; ++proc) {
8885680f57bSMatthew G. Knepley       ierr = PetscSectionSetDof(p->section, proc, sizes[proc]);CHKERRQ(ierr);
8895680f57bSMatthew G. Knepley     }
8905680f57bSMatthew G. Knepley   }
8915680f57bSMatthew G. Knepley   ierr = PetscSectionSetUp(p->section);CHKERRQ(ierr);
8925680f57bSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(p->section, &numPoints);CHKERRQ(ierr);
8935680f57bSMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) part), numPoints, points, PETSC_COPY_VALUES, &p->partition);CHKERRQ(ierr);
8945680f57bSMatthew G. Knepley   PetscFunctionReturn(0);
8955680f57bSMatthew G. Knepley }
8965680f57bSMatthew G. Knepley 
8975680f57bSMatthew G. Knepley #undef __FUNCT__
89877623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerDestroy_Chaco"
89977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy_Chaco(PetscPartitioner part)
90077623264SMatthew G. Knepley {
90177623264SMatthew G. Knepley   PetscPartitioner_Chaco *p = (PetscPartitioner_Chaco *) part->data;
90277623264SMatthew G. Knepley   PetscErrorCode          ierr;
90377623264SMatthew G. Knepley 
90477623264SMatthew G. Knepley   PetscFunctionBegin;
90577623264SMatthew G. Knepley   ierr = PetscFree(p);CHKERRQ(ierr);
90677623264SMatthew G. Knepley   PetscFunctionReturn(0);
90777623264SMatthew G. Knepley }
90877623264SMatthew G. Knepley 
90977623264SMatthew G. Knepley #undef __FUNCT__
91077623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_Chaco_Ascii"
91177623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Chaco_Ascii(PetscPartitioner part, PetscViewer viewer)
91277623264SMatthew G. Knepley {
91377623264SMatthew G. Knepley   PetscViewerFormat format;
91477623264SMatthew G. Knepley   PetscErrorCode    ierr;
91577623264SMatthew G. Knepley 
91677623264SMatthew G. Knepley   PetscFunctionBegin;
91777623264SMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
91877623264SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Chaco Graph Partitioner:\n");CHKERRQ(ierr);
91977623264SMatthew G. Knepley   PetscFunctionReturn(0);
92077623264SMatthew G. Knepley }
92177623264SMatthew G. Knepley 
92277623264SMatthew G. Knepley #undef __FUNCT__
92377623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_Chaco"
92477623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Chaco(PetscPartitioner part, PetscViewer viewer)
92577623264SMatthew G. Knepley {
92677623264SMatthew G. Knepley   PetscBool      iascii;
92777623264SMatthew G. Knepley   PetscErrorCode ierr;
92877623264SMatthew G. Knepley 
92977623264SMatthew G. Knepley   PetscFunctionBegin;
93077623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
93177623264SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
93277623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
93377623264SMatthew G. Knepley   if (iascii) {ierr = PetscPartitionerView_Chaco_Ascii(part, viewer);CHKERRQ(ierr);}
93477623264SMatthew G. Knepley   PetscFunctionReturn(0);
93577623264SMatthew G. Knepley }
93677623264SMatthew G. Knepley 
93770034214SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO)
93870034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
93970034214SMatthew G. Knepley #include <unistd.h>
94070034214SMatthew G. Knepley #endif
94170034214SMatthew G. Knepley /* Chaco does not have an include file */
94270034214SMatthew G. Knepley PETSC_EXTERN int interface(int nvtxs, int *start, int *adjacency, int *vwgts,
94370034214SMatthew G. Knepley                        float *ewgts, float *x, float *y, float *z, char *outassignname,
94470034214SMatthew G. Knepley                        char *outfilename, short *assignment, int architecture, int ndims_tot,
94570034214SMatthew G. Knepley                        int mesh_dims[3], double *goal, int global_method, int local_method,
94670034214SMatthew G. Knepley                        int rqi_flag, int vmax, int ndims, double eigtol, long seed);
94770034214SMatthew G. Knepley 
94870034214SMatthew G. Knepley extern int FREE_GRAPH;
94977623264SMatthew G. Knepley #endif
95070034214SMatthew G. Knepley 
95170034214SMatthew G. Knepley #undef __FUNCT__
95277623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerPartition_Chaco"
95377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerPartition_Chaco(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
95470034214SMatthew G. Knepley {
95577623264SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO)
95670034214SMatthew G. Knepley   enum {DEFAULT_METHOD = 1, INERTIAL_METHOD = 3};
95770034214SMatthew G. Knepley   MPI_Comm       comm;
95870034214SMatthew G. Knepley   int            nvtxs          = numVertices; /* number of vertices in full graph */
95970034214SMatthew G. Knepley   int           *vwgts          = NULL;   /* weights for all vertices */
96070034214SMatthew G. Knepley   float         *ewgts          = NULL;   /* weights for all edges */
96170034214SMatthew G. Knepley   float         *x              = NULL, *y = NULL, *z = NULL; /* coordinates for inertial method */
96270034214SMatthew G. Knepley   char          *outassignname  = NULL;   /*  name of assignment output file */
96370034214SMatthew G. Knepley   char          *outfilename    = NULL;   /* output file name */
96470034214SMatthew G. Knepley   int            architecture   = 1;      /* 0 => hypercube, d => d-dimensional mesh */
96570034214SMatthew G. Knepley   int            ndims_tot      = 0;      /* total number of cube dimensions to divide */
96670034214SMatthew G. Knepley   int            mesh_dims[3];            /* dimensions of mesh of processors */
96770034214SMatthew G. Knepley   double        *goal          = NULL;    /* desired set sizes for each set */
96870034214SMatthew G. Knepley   int            global_method = 1;       /* global partitioning algorithm */
96970034214SMatthew G. Knepley   int            local_method  = 1;       /* local partitioning algorithm */
97070034214SMatthew G. Knepley   int            rqi_flag      = 0;       /* should I use RQI/Symmlq eigensolver? */
97170034214SMatthew G. Knepley   int            vmax          = 200;     /* how many vertices to coarsen down to? */
97270034214SMatthew G. Knepley   int            ndims         = 1;       /* number of eigenvectors (2^d sets) */
97370034214SMatthew G. Knepley   double         eigtol        = 0.001;   /* tolerance on eigenvectors */
97470034214SMatthew G. Knepley   long           seed          = 123636512; /* for random graph mutations */
97570034214SMatthew G. Knepley   short int     *assignment;              /* Output partition */
97670034214SMatthew G. Knepley   int            fd_stdout, fd_pipe[2];
97770034214SMatthew G. Knepley   PetscInt      *points;
97870034214SMatthew G. Knepley   int            i, v, p;
97970034214SMatthew G. Knepley   PetscErrorCode ierr;
98070034214SMatthew G. Knepley 
98170034214SMatthew G. Knepley   PetscFunctionBegin;
98270034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
98370034214SMatthew G. Knepley   if (!numVertices) {
98477623264SMatthew G. Knepley     ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
98577623264SMatthew G. Knepley     ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
98670034214SMatthew G. Knepley     ierr = ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
98770034214SMatthew G. Knepley     PetscFunctionReturn(0);
98870034214SMatthew G. Knepley   }
98970034214SMatthew G. Knepley   FREE_GRAPH = 0;                         /* Do not let Chaco free my memory */
99070034214SMatthew G. Knepley   for (i = 0; i < start[numVertices]; ++i) ++adjacency[i];
99170034214SMatthew G. Knepley 
99270034214SMatthew G. Knepley   if (global_method == INERTIAL_METHOD) {
99370034214SMatthew G. Knepley     /* manager.createCellCoordinates(nvtxs, &x, &y, &z); */
99470034214SMatthew G. Knepley     SETERRQ(comm, PETSC_ERR_SUP, "Inertial partitioning not yet supported");
99570034214SMatthew G. Knepley   }
99677623264SMatthew G. Knepley   mesh_dims[0] = nparts;
99770034214SMatthew G. Knepley   mesh_dims[1] = 1;
99870034214SMatthew G. Knepley   mesh_dims[2] = 1;
99970034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &assignment);CHKERRQ(ierr);
100070034214SMatthew G. Knepley   /* Chaco outputs to stdout. We redirect this to a buffer. */
100170034214SMatthew G. Knepley   /* TODO: check error codes for UNIX calls */
100270034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
100370034214SMatthew G. Knepley   {
100470034214SMatthew G. Knepley     int piperet;
100570034214SMatthew G. Knepley     piperet = pipe(fd_pipe);
100670034214SMatthew G. Knepley     if (piperet) SETERRQ(comm,PETSC_ERR_SYS,"Could not create pipe");
100770034214SMatthew G. Knepley     fd_stdout = dup(1);
100870034214SMatthew G. Knepley     close(1);
100970034214SMatthew G. Knepley     dup2(fd_pipe[1], 1);
101070034214SMatthew G. Knepley   }
101170034214SMatthew G. Knepley #endif
101270034214SMatthew G. Knepley   ierr = interface(nvtxs, (int*) start, (int*) adjacency, vwgts, ewgts, x, y, z, outassignname, outfilename,
101370034214SMatthew G. Knepley                    assignment, architecture, ndims_tot, mesh_dims, goal, global_method, local_method, rqi_flag,
101470034214SMatthew G. Knepley                    vmax, ndims, eigtol, seed);
101570034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
101670034214SMatthew G. Knepley   {
101770034214SMatthew G. Knepley     char msgLog[10000];
101870034214SMatthew G. Knepley     int  count;
101970034214SMatthew G. Knepley 
102070034214SMatthew G. Knepley     fflush(stdout);
102170034214SMatthew G. Knepley     count = read(fd_pipe[0], msgLog, (10000-1)*sizeof(char));
102270034214SMatthew G. Knepley     if (count < 0) count = 0;
102370034214SMatthew G. Knepley     msgLog[count] = 0;
102470034214SMatthew G. Knepley     close(1);
102570034214SMatthew G. Knepley     dup2(fd_stdout, 1);
102670034214SMatthew G. Knepley     close(fd_stdout);
102770034214SMatthew G. Knepley     close(fd_pipe[0]);
102870034214SMatthew G. Knepley     close(fd_pipe[1]);
102970034214SMatthew G. Knepley     if (ierr) SETERRQ1(comm, PETSC_ERR_LIB, "Error in Chaco library: %s", msgLog);
103070034214SMatthew G. Knepley   }
103170034214SMatthew G. Knepley #endif
103270034214SMatthew G. Knepley   /* Convert to PetscSection+IS */
103377623264SMatthew G. Knepley   ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
103470034214SMatthew G. Knepley   for (v = 0; v < nvtxs; ++v) {
103577623264SMatthew G. Knepley     ierr = PetscSectionAddDof(partSection, assignment[v], 1);CHKERRQ(ierr);
103670034214SMatthew G. Knepley   }
103777623264SMatthew G. Knepley   ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
103870034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr);
103977623264SMatthew G. Knepley   for (p = 0, i = 0; p < nparts; ++p) {
104070034214SMatthew G. Knepley     for (v = 0; v < nvtxs; ++v) {
104170034214SMatthew G. Knepley       if (assignment[v] == p) points[i++] = v;
104270034214SMatthew G. Knepley     }
104370034214SMatthew G. Knepley   }
104470034214SMatthew G. Knepley   if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs);
104570034214SMatthew G. Knepley   ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
104670034214SMatthew G. Knepley   if (global_method == INERTIAL_METHOD) {
104770034214SMatthew G. Knepley     /* manager.destroyCellCoordinates(nvtxs, &x, &y, &z); */
104870034214SMatthew G. Knepley   }
104970034214SMatthew G. Knepley   ierr = PetscFree(assignment);CHKERRQ(ierr);
105070034214SMatthew G. Knepley   for (i = 0; i < start[numVertices]; ++i) --adjacency[i];
105170034214SMatthew G. Knepley   PetscFunctionReturn(0);
105277623264SMatthew G. Knepley #else
105377623264SMatthew G. Knepley   SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-chaco.");
105470034214SMatthew G. Knepley #endif
105577623264SMatthew G. Knepley }
105677623264SMatthew G. Knepley 
105777623264SMatthew G. Knepley #undef __FUNCT__
105877623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerInitialize_Chaco"
105977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerInitialize_Chaco(PetscPartitioner part)
106077623264SMatthew G. Knepley {
106177623264SMatthew G. Knepley   PetscFunctionBegin;
106277623264SMatthew G. Knepley   part->ops->view      = PetscPartitionerView_Chaco;
106377623264SMatthew G. Knepley   part->ops->destroy   = PetscPartitionerDestroy_Chaco;
106477623264SMatthew G. Knepley   part->ops->partition = PetscPartitionerPartition_Chaco;
106577623264SMatthew G. Knepley   PetscFunctionReturn(0);
106677623264SMatthew G. Knepley }
106777623264SMatthew G. Knepley 
106877623264SMatthew G. Knepley /*MC
106977623264SMatthew G. Knepley   PETSCPARTITIONERCHACO = "chaco" - A PetscPartitioner object using the Chaco library
107077623264SMatthew G. Knepley 
107177623264SMatthew G. Knepley   Level: intermediate
107277623264SMatthew G. Knepley 
107377623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
107477623264SMatthew G. Knepley M*/
107577623264SMatthew G. Knepley 
107677623264SMatthew G. Knepley #undef __FUNCT__
107777623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerCreate_Chaco"
107877623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Chaco(PetscPartitioner part)
107977623264SMatthew G. Knepley {
108077623264SMatthew G. Knepley   PetscPartitioner_Chaco *p;
108177623264SMatthew G. Knepley   PetscErrorCode          ierr;
108277623264SMatthew G. Knepley 
108377623264SMatthew G. Knepley   PetscFunctionBegin;
108477623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
108577623264SMatthew G. Knepley   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
108677623264SMatthew G. Knepley   part->data = p;
108777623264SMatthew G. Knepley 
108877623264SMatthew G. Knepley   ierr = PetscPartitionerInitialize_Chaco(part);CHKERRQ(ierr);
108977623264SMatthew G. Knepley   ierr = PetscCitationsRegister(ChacoPartitionerCitation, &ChacoPartitionercite);CHKERRQ(ierr);
109077623264SMatthew G. Knepley   PetscFunctionReturn(0);
109177623264SMatthew G. Knepley }
109277623264SMatthew G. Knepley 
109377623264SMatthew G. Knepley #undef __FUNCT__
109477623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerDestroy_ParMetis"
109577623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy_ParMetis(PetscPartitioner part)
109677623264SMatthew G. Knepley {
109777623264SMatthew G. Knepley   PetscPartitioner_ParMetis *p = (PetscPartitioner_ParMetis *) part->data;
109877623264SMatthew G. Knepley   PetscErrorCode             ierr;
109977623264SMatthew G. Knepley 
110077623264SMatthew G. Knepley   PetscFunctionBegin;
110177623264SMatthew G. Knepley   ierr = PetscFree(p);CHKERRQ(ierr);
110277623264SMatthew G. Knepley   PetscFunctionReturn(0);
110377623264SMatthew G. Knepley }
110477623264SMatthew G. Knepley 
110577623264SMatthew G. Knepley #undef __FUNCT__
110677623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_ParMetis_Ascii"
110777623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_ParMetis_Ascii(PetscPartitioner part, PetscViewer viewer)
110877623264SMatthew G. Knepley {
110977623264SMatthew G. Knepley   PetscViewerFormat format;
111077623264SMatthew G. Knepley   PetscErrorCode    ierr;
111177623264SMatthew G. Knepley 
111277623264SMatthew G. Knepley   PetscFunctionBegin;
111377623264SMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
111477623264SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "ParMetis Graph Partitioner:\n");CHKERRQ(ierr);
111577623264SMatthew G. Knepley   PetscFunctionReturn(0);
111677623264SMatthew G. Knepley }
111777623264SMatthew G. Knepley 
111877623264SMatthew G. Knepley #undef __FUNCT__
111977623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_ParMetis"
112077623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_ParMetis(PetscPartitioner part, PetscViewer viewer)
112177623264SMatthew G. Knepley {
112277623264SMatthew G. Knepley   PetscBool      iascii;
112377623264SMatthew G. Knepley   PetscErrorCode ierr;
112477623264SMatthew G. Knepley 
112577623264SMatthew G. Knepley   PetscFunctionBegin;
112677623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
112777623264SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
112877623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
112977623264SMatthew G. Knepley   if (iascii) {ierr = PetscPartitionerView_ParMetis_Ascii(part, viewer);CHKERRQ(ierr);}
113077623264SMatthew G. Knepley   PetscFunctionReturn(0);
113177623264SMatthew G. Knepley }
113270034214SMatthew G. Knepley 
113370034214SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS)
113470034214SMatthew G. Knepley #include <parmetis.h>
113577623264SMatthew G. Knepley #endif
113670034214SMatthew G. Knepley 
113770034214SMatthew G. Knepley #undef __FUNCT__
113877623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerPartition_ParMetis"
113977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerPartition_ParMetis(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
114070034214SMatthew G. Knepley {
114177623264SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS)
114270034214SMatthew G. Knepley   MPI_Comm       comm;
114370034214SMatthew G. Knepley   PetscInt       nvtxs      = numVertices; /* The number of vertices in full graph */
114470034214SMatthew G. Knepley   PetscInt      *vtxdist;                  /* Distribution of vertices across processes */
114570034214SMatthew G. Knepley   PetscInt      *xadj       = start;       /* Start of edge list for each vertex */
114670034214SMatthew G. Knepley   PetscInt      *adjncy     = adjacency;   /* Edge lists for all vertices */
114770034214SMatthew G. Knepley   PetscInt      *vwgt       = NULL;        /* Vertex weights */
114870034214SMatthew G. Knepley   PetscInt      *adjwgt     = NULL;        /* Edge weights */
114970034214SMatthew G. Knepley   PetscInt       wgtflag    = 0;           /* Indicates which weights are present */
115070034214SMatthew G. Knepley   PetscInt       numflag    = 0;           /* Indicates initial offset (0 or 1) */
115170034214SMatthew G. Knepley   PetscInt       ncon       = 1;           /* The number of weights per vertex */
115270034214SMatthew G. Knepley   PetscReal     *tpwgts;                   /* The fraction of vertex weights assigned to each partition */
115370034214SMatthew G. Knepley   PetscReal     *ubvec;                    /* The balance intolerance for vertex weights */
115470034214SMatthew G. Knepley   PetscInt       options[5];               /* Options */
115570034214SMatthew G. Knepley   /* Outputs */
115670034214SMatthew G. Knepley   PetscInt       edgeCut;                  /* The number of edges cut by the partition */
115770034214SMatthew G. Knepley   PetscInt      *assignment, *points;
115877623264SMatthew G. Knepley   PetscMPIInt    rank, p, v, i;
115970034214SMatthew G. Knepley   PetscErrorCode ierr;
116070034214SMatthew G. Knepley 
116170034214SMatthew G. Knepley   PetscFunctionBegin;
116277623264SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) part, &comm);CHKERRQ(ierr);
116370034214SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
116470034214SMatthew G. Knepley   options[0] = 0; /* Use all defaults */
116570034214SMatthew G. Knepley   /* Calculate vertex distribution */
116670034214SMatthew G. Knepley   ierr = PetscMalloc4(nparts+1,&vtxdist,nparts*ncon,&tpwgts,ncon,&ubvec,nvtxs,&assignment);CHKERRQ(ierr);
116770034214SMatthew G. Knepley   vtxdist[0] = 0;
116870034214SMatthew G. Knepley   ierr = MPI_Allgather(&nvtxs, 1, MPIU_INT, &vtxdist[1], 1, MPIU_INT, comm);CHKERRQ(ierr);
116970034214SMatthew G. Knepley   for (p = 2; p <= nparts; ++p) {
117070034214SMatthew G. Knepley     vtxdist[p] += vtxdist[p-1];
117170034214SMatthew G. Knepley   }
117270034214SMatthew G. Knepley   /* Calculate weights */
117370034214SMatthew G. Knepley   for (p = 0; p < nparts; ++p) {
117470034214SMatthew G. Knepley     tpwgts[p] = 1.0/nparts;
117570034214SMatthew G. Knepley   }
117670034214SMatthew G. Knepley   ubvec[0] = 1.05;
117770034214SMatthew G. Knepley 
117870034214SMatthew G. Knepley   if (nparts == 1) {
117970034214SMatthew G. Knepley     ierr = PetscMemzero(assignment, nvtxs * sizeof(PetscInt));
118070034214SMatthew G. Knepley   } else {
118170034214SMatthew G. Knepley     if (vtxdist[1] == vtxdist[nparts]) {
118270034214SMatthew G. Knepley       if (!rank) {
118370034214SMatthew G. Knepley         PetscStackPush("METIS_PartGraphKway");
118470034214SMatthew G. Knepley         ierr = METIS_PartGraphKway(&nvtxs, &ncon, xadj, adjncy, vwgt, NULL, adjwgt, &nparts, tpwgts, ubvec, NULL, &edgeCut, assignment);
118570034214SMatthew G. Knepley         PetscStackPop;
118670034214SMatthew G. Knepley         if (ierr != METIS_OK) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in METIS_PartGraphKway()");
118770034214SMatthew G. Knepley       }
118870034214SMatthew G. Knepley     } else {
118970034214SMatthew G. Knepley       PetscStackPush("ParMETIS_V3_PartKway");
119070034214SMatthew G. Knepley       ierr = ParMETIS_V3_PartKway(vtxdist, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgeCut, assignment, &comm);
119170034214SMatthew G. Knepley       PetscStackPop;
119270034214SMatthew G. Knepley       if (ierr != METIS_OK) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in ParMETIS_V3_PartKway()");
119370034214SMatthew G. Knepley     }
119470034214SMatthew G. Knepley   }
119570034214SMatthew G. Knepley   /* Convert to PetscSection+IS */
119677623264SMatthew G. Knepley   ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
119777623264SMatthew G. Knepley   for (v = 0; v < nvtxs; ++v) {ierr = PetscSectionAddDof(partSection, assignment[v], 1);CHKERRQ(ierr);}
119877623264SMatthew G. Knepley   ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
119970034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr);
120077623264SMatthew G. Knepley   for (p = 0, i = 0; p < nparts; ++p) {
120170034214SMatthew G. Knepley     for (v = 0; v < nvtxs; ++v) {
120270034214SMatthew G. Knepley       if (assignment[v] == p) points[i++] = v;
120370034214SMatthew G. Knepley     }
120470034214SMatthew G. Knepley   }
120570034214SMatthew G. Knepley   if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs);
120670034214SMatthew G. Knepley   ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
120770034214SMatthew G. Knepley   ierr = PetscFree4(vtxdist,tpwgts,ubvec,assignment);CHKERRQ(ierr);
120870034214SMatthew G. Knepley #else
120977623264SMatthew G. Knepley   SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-parmetis.");
121070034214SMatthew G. Knepley #endif
121177623264SMatthew G. Knepley   PetscFunctionReturn(0);
121270034214SMatthew G. Knepley }
121370034214SMatthew G. Knepley 
121477623264SMatthew G. Knepley #undef __FUNCT__
121577623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerInitialize_ParMetis"
121677623264SMatthew G. Knepley PetscErrorCode PetscPartitionerInitialize_ParMetis(PetscPartitioner part)
121777623264SMatthew G. Knepley {
121877623264SMatthew G. Knepley   PetscFunctionBegin;
121977623264SMatthew G. Knepley   part->ops->view      = PetscPartitionerView_ParMetis;
122077623264SMatthew G. Knepley   part->ops->destroy   = PetscPartitionerDestroy_ParMetis;
122177623264SMatthew G. Knepley   part->ops->partition = PetscPartitionerPartition_ParMetis;
122277623264SMatthew G. Knepley   PetscFunctionReturn(0);
122377623264SMatthew G. Knepley }
122477623264SMatthew G. Knepley 
122577623264SMatthew G. Knepley /*MC
122677623264SMatthew G. Knepley   PETSCPARTITIONERPARMETIS = "parmetis" - A PetscPartitioner object using the ParMetis library
122777623264SMatthew G. Knepley 
122877623264SMatthew G. Knepley   Level: intermediate
122977623264SMatthew G. Knepley 
123077623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
123177623264SMatthew G. Knepley M*/
123277623264SMatthew G. Knepley 
123377623264SMatthew G. Knepley #undef __FUNCT__
123477623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerCreate_ParMetis"
123577623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_ParMetis(PetscPartitioner part)
123677623264SMatthew G. Knepley {
123777623264SMatthew G. Knepley   PetscPartitioner_ParMetis *p;
123877623264SMatthew G. Knepley   PetscErrorCode          ierr;
123977623264SMatthew G. Knepley 
124077623264SMatthew G. Knepley   PetscFunctionBegin;
124177623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
124277623264SMatthew G. Knepley   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
124377623264SMatthew G. Knepley   part->data = p;
124477623264SMatthew G. Knepley 
124577623264SMatthew G. Knepley   ierr = PetscPartitionerInitialize_ParMetis(part);CHKERRQ(ierr);
124677623264SMatthew G. Knepley   ierr = PetscCitationsRegister(ParMetisPartitionerCitation, &ParMetisPartitionercite);CHKERRQ(ierr);
124770034214SMatthew G. Knepley   PetscFunctionReturn(0);
124870034214SMatthew G. Knepley }
124970034214SMatthew G. Knepley 
125070034214SMatthew G. Knepley #undef __FUNCT__
12515680f57bSMatthew G. Knepley #define __FUNCT__ "DMPlexGetPartitioner"
12525680f57bSMatthew G. Knepley /*@
12535680f57bSMatthew G. Knepley   DMPlexGetPartitioner - Get the mesh partitioner
12545680f57bSMatthew G. Knepley 
12555680f57bSMatthew G. Knepley   Not collective
12565680f57bSMatthew G. Knepley 
12575680f57bSMatthew G. Knepley   Input Parameter:
12585680f57bSMatthew G. Knepley . dm - The DM
12595680f57bSMatthew G. Knepley 
12605680f57bSMatthew G. Knepley   Output Parameter:
12615680f57bSMatthew G. Knepley . part - The PetscPartitioner
12625680f57bSMatthew G. Knepley 
12635680f57bSMatthew G. Knepley   Level: developer
12645680f57bSMatthew G. Knepley 
12655680f57bSMatthew G. Knepley .seealso DMPlexDistribute(), PetscPartitionerCreate()
12665680f57bSMatthew G. Knepley @*/
12675680f57bSMatthew G. Knepley PetscErrorCode DMPlexGetPartitioner(DM dm, PetscPartitioner *part)
12685680f57bSMatthew G. Knepley {
12695680f57bSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
12705680f57bSMatthew G. Knepley 
12715680f57bSMatthew G. Knepley   PetscFunctionBegin;
12725680f57bSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12735680f57bSMatthew G. Knepley   PetscValidPointer(part, 2);
12745680f57bSMatthew G. Knepley   *part = mesh->partitioner;
12755680f57bSMatthew G. Knepley   PetscFunctionReturn(0);
12765680f57bSMatthew G. Knepley }
12775680f57bSMatthew G. Knepley 
12785680f57bSMatthew G. Knepley #undef __FUNCT__
12799ffc88e4SToby Isaac #define __FUNCT__ "DMPlexMarkTreeClosure"
12809ffc88e4SToby Isaac static PetscErrorCode DMPlexMarkTreeClosure(DM dm, PetscSegBuffer segpart, PetscBT bt, PetscInt point, PetscInt *partSize)
12819ffc88e4SToby Isaac {
12829ffc88e4SToby Isaac   PetscInt       parent, closureSize, *closure = NULL, i, pStart, pEnd;
12839ffc88e4SToby Isaac   PetscErrorCode ierr;
12849ffc88e4SToby Isaac 
12859ffc88e4SToby Isaac   PetscFunctionBegin;
12869ffc88e4SToby Isaac   ierr = DMPlexGetTreeParent(dm,point,&parent,NULL);CHKERRQ(ierr);
12879ffc88e4SToby Isaac   if (parent == point) PetscFunctionReturn(0);
12889ffc88e4SToby Isaac   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
12899ffc88e4SToby Isaac   ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
12909ffc88e4SToby Isaac   for (i = 0; i < closureSize; i++) {
12919ffc88e4SToby Isaac     PetscInt cpoint = closure[2*i];
12929ffc88e4SToby Isaac 
12939ffc88e4SToby Isaac     if (!PetscBTLookupSet(bt,cpoint-pStart)) {
12949ffc88e4SToby Isaac       PetscInt *PETSC_RESTRICT pt;
12959ffc88e4SToby Isaac       (*partSize)++;
12969ffc88e4SToby Isaac       ierr = PetscSegBufferGetInts(segpart,1,&pt);CHKERRQ(ierr);
12979ffc88e4SToby Isaac       *pt = cpoint;
12989ffc88e4SToby Isaac     }
12999ffc88e4SToby Isaac     ierr = DMPlexMarkTreeClosure(dm,segpart,bt,cpoint,partSize);CHKERRQ(ierr);
13009ffc88e4SToby Isaac   }
13019ffc88e4SToby Isaac   ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
13029ffc88e4SToby Isaac   PetscFunctionReturn(0);
13039ffc88e4SToby Isaac }
13049ffc88e4SToby Isaac 
13059ffc88e4SToby Isaac #undef __FUNCT__
130670034214SMatthew G. Knepley #define __FUNCT__ "DMPlexCreatePartitionClosure"
130770034214SMatthew G. Knepley PetscErrorCode DMPlexCreatePartitionClosure(DM dm, PetscSection pointSection, IS pointPartition, PetscSection *section, IS *partition)
130870034214SMatthew G. Knepley {
130970034214SMatthew G. Knepley   /* const PetscInt  height = 0; */
131070034214SMatthew G. Knepley   const PetscInt *partArray;
131170034214SMatthew G. Knepley   PetscInt       *allPoints, *packPoints;
131270034214SMatthew G. Knepley   PetscInt        rStart, rEnd, rank, pStart, pEnd, newSize;
131370034214SMatthew G. Knepley   PetscErrorCode  ierr;
131470034214SMatthew G. Knepley   PetscBT         bt;
131570034214SMatthew G. Knepley   PetscSegBuffer  segpack,segpart;
131670034214SMatthew G. Knepley 
131770034214SMatthew G. Knepley   PetscFunctionBegin;
131870034214SMatthew G. Knepley   ierr = PetscSectionGetChart(pointSection, &rStart, &rEnd);CHKERRQ(ierr);
131970034214SMatthew G. Knepley   ierr = ISGetIndices(pointPartition, &partArray);CHKERRQ(ierr);
132070034214SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), section);CHKERRQ(ierr);
132170034214SMatthew G. Knepley   ierr = PetscSectionSetChart(*section, rStart, rEnd);CHKERRQ(ierr);
132270034214SMatthew G. Knepley   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
132370034214SMatthew G. Knepley   ierr = PetscBTCreate(pEnd-pStart,&bt);CHKERRQ(ierr);
132470034214SMatthew G. Knepley   ierr = PetscSegBufferCreate(sizeof(PetscInt),1000,&segpack);CHKERRQ(ierr);
132570034214SMatthew G. Knepley   ierr = PetscSegBufferCreate(sizeof(PetscInt),1000,&segpart);CHKERRQ(ierr);
132670034214SMatthew G. Knepley   for (rank = rStart; rank < rEnd; ++rank) {
132770034214SMatthew G. Knepley     PetscInt partSize = 0, numPoints, offset, p, *PETSC_RESTRICT placePoints;
132870034214SMatthew G. Knepley 
132970034214SMatthew G. Knepley     ierr = PetscSectionGetDof(pointSection, rank, &numPoints);CHKERRQ(ierr);
133070034214SMatthew G. Knepley     ierr = PetscSectionGetOffset(pointSection, rank, &offset);CHKERRQ(ierr);
133170034214SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
133270034214SMatthew G. Knepley       PetscInt  point   = partArray[offset+p], closureSize, c;
133370034214SMatthew G. Knepley       PetscInt *closure = NULL;
133470034214SMatthew G. Knepley 
133570034214SMatthew G. Knepley       /* TODO Include support for height > 0 case */
133670034214SMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
133770034214SMatthew G. Knepley       for (c=0; c<closureSize; c++) {
133870034214SMatthew G. Knepley         PetscInt cpoint = closure[c*2];
133970034214SMatthew G. Knepley         if (!PetscBTLookupSet(bt,cpoint-pStart)) {
134070034214SMatthew G. Knepley           PetscInt *PETSC_RESTRICT pt;
134170034214SMatthew G. Knepley           partSize++;
134270034214SMatthew G. Knepley           ierr = PetscSegBufferGetInts(segpart,1,&pt);CHKERRQ(ierr);
134370034214SMatthew G. Knepley           *pt = cpoint;
134470034214SMatthew G. Knepley         }
13459ffc88e4SToby Isaac         ierr = DMPlexMarkTreeClosure(dm,segpart,bt,cpoint,&partSize);CHKERRQ(ierr);
134670034214SMatthew G. Knepley       }
134770034214SMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
134870034214SMatthew G. Knepley     }
134970034214SMatthew G. Knepley     ierr = PetscSectionSetDof(*section, rank, partSize);CHKERRQ(ierr);
135070034214SMatthew G. Knepley     ierr = PetscSegBufferGetInts(segpack,partSize,&placePoints);CHKERRQ(ierr);
135170034214SMatthew G. Knepley     ierr = PetscSegBufferExtractTo(segpart,placePoints);CHKERRQ(ierr);
135270034214SMatthew G. Knepley     ierr = PetscSortInt(partSize,placePoints);CHKERRQ(ierr);
135370034214SMatthew G. Knepley     for (p=0; p<partSize; p++) {ierr = PetscBTClear(bt,placePoints[p]-pStart);CHKERRQ(ierr);}
135470034214SMatthew G. Knepley   }
135570034214SMatthew G. Knepley   ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
135670034214SMatthew G. Knepley   ierr = PetscSegBufferDestroy(&segpart);CHKERRQ(ierr);
135770034214SMatthew G. Knepley 
135870034214SMatthew G. Knepley   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
135970034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(*section, &newSize);CHKERRQ(ierr);
136070034214SMatthew G. Knepley   ierr = PetscMalloc1(newSize, &allPoints);CHKERRQ(ierr);
136170034214SMatthew G. Knepley 
136270034214SMatthew G. Knepley   ierr = PetscSegBufferExtractInPlace(segpack,&packPoints);CHKERRQ(ierr);
136370034214SMatthew G. Knepley   for (rank = rStart; rank < rEnd; ++rank) {
136470034214SMatthew G. Knepley     PetscInt numPoints, offset;
136570034214SMatthew G. Knepley 
136670034214SMatthew G. Knepley     ierr = PetscSectionGetDof(*section, rank, &numPoints);CHKERRQ(ierr);
136770034214SMatthew G. Knepley     ierr = PetscSectionGetOffset(*section, rank, &offset);CHKERRQ(ierr);
136870034214SMatthew G. Knepley     ierr = PetscMemcpy(&allPoints[offset], packPoints, numPoints * sizeof(PetscInt));CHKERRQ(ierr);
136970034214SMatthew G. Knepley     packPoints += numPoints;
137070034214SMatthew G. Knepley   }
137170034214SMatthew G. Knepley 
137270034214SMatthew G. Knepley   ierr = PetscSegBufferDestroy(&segpack);CHKERRQ(ierr);
137370034214SMatthew G. Knepley   ierr = ISRestoreIndices(pointPartition, &partArray);CHKERRQ(ierr);
137470034214SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), newSize, allPoints, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
137570034214SMatthew G. Knepley   PetscFunctionReturn(0);
137670034214SMatthew G. Knepley }
1377aa3148a8SMichael Lange 
1378aa3148a8SMichael Lange #undef __FUNCT__
13795abbe4feSMichael Lange #define __FUNCT__ "DMPlexPartitionLabelClosure"
13805abbe4feSMichael Lange /*@
13815abbe4feSMichael Lange   DMPlexPartitionLabelClosure - Add the closure of all points to the partition label
13825abbe4feSMichael Lange 
13835abbe4feSMichael Lange   Input Parameters:
13845abbe4feSMichael Lange + dm     - The DM
13855abbe4feSMichael Lange - label  - DMLabel assinging ranks to remote roots
13865abbe4feSMichael Lange 
13875abbe4feSMichael Lange   Level: developer
13885abbe4feSMichael Lange 
13895abbe4feSMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap
13905abbe4feSMichael Lange @*/
13915abbe4feSMichael Lange PetscErrorCode DMPlexPartitionLabelClosure(DM dm, DMLabel label)
13925abbe4feSMichael Lange {
13935abbe4feSMichael Lange   IS              rankIS,   pointIS;
13945abbe4feSMichael Lange   const PetscInt *ranks,   *points;
13955abbe4feSMichael Lange   PetscInt        numRanks, numPoints, r, p, c, closureSize;
13965abbe4feSMichael Lange   PetscInt       *closure = NULL;
13975abbe4feSMichael Lange   PetscErrorCode  ierr;
13985abbe4feSMichael Lange 
13995abbe4feSMichael Lange   PetscFunctionBegin;
14005abbe4feSMichael Lange   ierr = DMLabelGetValueIS(label, &rankIS);CHKERRQ(ierr);
14015abbe4feSMichael Lange   ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr);
14025abbe4feSMichael Lange   ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr);
14035abbe4feSMichael Lange   for (r = 0; r < numRanks; ++r) {
14045abbe4feSMichael Lange     const PetscInt rank = ranks[r];
14055abbe4feSMichael Lange 
14065abbe4feSMichael Lange     ierr = DMLabelGetStratumIS(label, rank, &pointIS);CHKERRQ(ierr);
14075abbe4feSMichael Lange     ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
14085abbe4feSMichael Lange     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
14095abbe4feSMichael Lange     for (p = 0; p < numPoints; ++p) {
14105abbe4feSMichael Lange       ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
14115abbe4feSMichael Lange       for (c = 0; c < closureSize*2; c += 2) {ierr = DMLabelSetValue(label, closure[c], rank);CHKERRQ(ierr);}
14125abbe4feSMichael Lange     }
14135abbe4feSMichael Lange     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
14145abbe4feSMichael Lange     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
14155abbe4feSMichael Lange   }
14165abbe4feSMichael Lange   if (closure) {ierr = DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);}
14175abbe4feSMichael Lange   ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr);
14185abbe4feSMichael Lange   ierr = ISDestroy(&rankIS);CHKERRQ(ierr);
14195abbe4feSMichael Lange   PetscFunctionReturn(0);
14205abbe4feSMichael Lange }
14215abbe4feSMichael Lange 
14225abbe4feSMichael Lange #undef __FUNCT__
142324d039d7SMichael Lange #define __FUNCT__ "DMPlexPartitionLabelAdjacency"
142424d039d7SMichael Lange /*@
142524d039d7SMichael Lange   DMPlexPartitionLabelAdjacency - Add one level of adjacent points to the partition label
142624d039d7SMichael Lange 
142724d039d7SMichael Lange   Input Parameters:
142824d039d7SMichael Lange + dm     - The DM
142924d039d7SMichael Lange - label  - DMLabel assinging ranks to remote roots
143024d039d7SMichael Lange 
143124d039d7SMichael Lange   Level: developer
143224d039d7SMichael Lange 
143324d039d7SMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap
143424d039d7SMichael Lange @*/
143524d039d7SMichael Lange PetscErrorCode DMPlexPartitionLabelAdjacency(DM dm, DMLabel label)
143624d039d7SMichael Lange {
143724d039d7SMichael Lange   IS              rankIS,   pointIS;
143824d039d7SMichael Lange   const PetscInt *ranks,   *points;
143924d039d7SMichael Lange   PetscInt        numRanks, numPoints, r, p, a, adjSize;
144024d039d7SMichael Lange   PetscInt       *adj = NULL;
144124d039d7SMichael Lange   PetscErrorCode  ierr;
144224d039d7SMichael Lange 
144324d039d7SMichael Lange   PetscFunctionBegin;
144424d039d7SMichael Lange   ierr = DMLabelGetValueIS(label, &rankIS);CHKERRQ(ierr);
144524d039d7SMichael Lange   ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr);
144624d039d7SMichael Lange   ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr);
144724d039d7SMichael Lange   for (r = 0; r < numRanks; ++r) {
144824d039d7SMichael Lange     const PetscInt rank = ranks[r];
144924d039d7SMichael Lange 
145024d039d7SMichael Lange     ierr = DMLabelGetStratumIS(label, rank, &pointIS);CHKERRQ(ierr);
145124d039d7SMichael Lange     ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
145224d039d7SMichael Lange     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
145324d039d7SMichael Lange     for (p = 0; p < numPoints; ++p) {
145424d039d7SMichael Lange       adjSize = PETSC_DETERMINE;
145524d039d7SMichael Lange       ierr = DMPlexGetAdjacency(dm, points[p], &adjSize, &adj);CHKERRQ(ierr);
145624d039d7SMichael Lange       for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(label, adj[a], rank);CHKERRQ(ierr);}
145724d039d7SMichael Lange     }
145824d039d7SMichael Lange     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
145924d039d7SMichael Lange     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
146024d039d7SMichael Lange   }
146124d039d7SMichael Lange   ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr);
146224d039d7SMichael Lange   ierr = ISDestroy(&rankIS);CHKERRQ(ierr);
146324d039d7SMichael Lange   ierr = PetscFree(adj);CHKERRQ(ierr);
146424d039d7SMichael Lange   PetscFunctionReturn(0);
146524d039d7SMichael Lange }
146624d039d7SMichael Lange 
146724d039d7SMichael Lange #undef __FUNCT__
14681fd9873aSMichael Lange #define __FUNCT__ "DMPlexPartitionLabelInvert"
14691fd9873aSMichael Lange /*@
14701fd9873aSMichael Lange   DMPlexPartitionLabelInvert - Create a partition label of remote roots from a local root label
14711fd9873aSMichael Lange 
14721fd9873aSMichael Lange   Input Parameters:
14731fd9873aSMichael Lange + dm        - The DM
14741fd9873aSMichael Lange . rootLabel - DMLabel assinging ranks to local roots
14751fd9873aSMichael Lange . processSF - A star forest mapping into the local index on each remote rank
14761fd9873aSMichael Lange 
14771fd9873aSMichael Lange   Output Parameter:
14781fd9873aSMichael Lange - leafLabel - DMLabel assinging ranks to remote roots
14791fd9873aSMichael Lange 
14801fd9873aSMichael Lange   Note: The rootLabel defines a send pattern by mapping local points to remote target ranks. The
14811fd9873aSMichael Lange   resulting leafLabel is a receiver mapping of remote roots to their parent rank.
14821fd9873aSMichael Lange 
14831fd9873aSMichael Lange   Level: developer
14841fd9873aSMichael Lange 
14851fd9873aSMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap
14861fd9873aSMichael Lange @*/
14871fd9873aSMichael Lange PetscErrorCode DMPlexPartitionLabelInvert(DM dm, DMLabel rootLabel, PetscSF processSF, DMLabel leafLabel)
14881fd9873aSMichael Lange {
14891fd9873aSMichael Lange   MPI_Comm           comm;
14901fd9873aSMichael Lange   PetscMPIInt        rank, numProcs;
14911fd9873aSMichael Lange   PetscInt           p, n, numNeighbors, size, l, nleaves;
14921fd9873aSMichael Lange   PetscSF            sfPoint;
14931fd9873aSMichael Lange   PetscSFNode       *rootPoints, *leafPoints;
14941fd9873aSMichael Lange   PetscSection       rootSection, leafSection;
14951fd9873aSMichael Lange   const PetscSFNode *remote;
14961fd9873aSMichael Lange   const PetscInt    *local, *neighbors;
14971fd9873aSMichael Lange   IS                 valueIS;
14981fd9873aSMichael Lange   PetscErrorCode     ierr;
14991fd9873aSMichael Lange 
15001fd9873aSMichael Lange   PetscFunctionBegin;
15011fd9873aSMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
15021fd9873aSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
15031fd9873aSMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
15041fd9873aSMichael Lange   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
15051fd9873aSMichael Lange 
15061fd9873aSMichael Lange   /* Convert to (point, rank) and use actual owners */
15071fd9873aSMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
15081fd9873aSMichael Lange   ierr = PetscSectionSetChart(rootSection, 0, numProcs);CHKERRQ(ierr);
15091fd9873aSMichael Lange   ierr = DMLabelGetValueIS(rootLabel, &valueIS);CHKERRQ(ierr);
15101fd9873aSMichael Lange   ierr = ISGetLocalSize(valueIS, &numNeighbors);CHKERRQ(ierr);
15111fd9873aSMichael Lange   ierr = ISGetIndices(valueIS, &neighbors);CHKERRQ(ierr);
15121fd9873aSMichael Lange   for (n = 0; n < numNeighbors; ++n) {
15131fd9873aSMichael Lange     PetscInt numPoints;
15141fd9873aSMichael Lange 
15151fd9873aSMichael Lange     ierr = DMLabelGetStratumSize(rootLabel, neighbors[n], &numPoints);CHKERRQ(ierr);
15161fd9873aSMichael Lange     ierr = PetscSectionAddDof(rootSection, neighbors[n], numPoints);CHKERRQ(ierr);
15171fd9873aSMichael Lange   }
15181fd9873aSMichael Lange   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
15191fd9873aSMichael Lange   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
15201fd9873aSMichael Lange   ierr = PetscMalloc1(size, &rootPoints);CHKERRQ(ierr);
15211fd9873aSMichael Lange   ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr);
15221fd9873aSMichael Lange   for (n = 0; n < numNeighbors; ++n) {
15231fd9873aSMichael Lange     IS              pointIS;
15241fd9873aSMichael Lange     const PetscInt *points;
15251fd9873aSMichael Lange     PetscInt        off, numPoints, p;
15261fd9873aSMichael Lange 
15271fd9873aSMichael Lange     ierr = PetscSectionGetOffset(rootSection, neighbors[n], &off);CHKERRQ(ierr);
15281fd9873aSMichael Lange     ierr = DMLabelGetStratumIS(rootLabel, neighbors[n], &pointIS);CHKERRQ(ierr);
15291fd9873aSMichael Lange     ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
15301fd9873aSMichael Lange     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
15311fd9873aSMichael Lange     for (p = 0; p < numPoints; ++p) {
1532f8987ae8SMichael Lange       if (local) {ierr = PetscFindInt(points[p], nleaves, local, &l);CHKERRQ(ierr);}
1533f8987ae8SMichael Lange       else       {l = -1;}
15341fd9873aSMichael Lange       if (l >= 0) {rootPoints[off+p] = remote[l];}
15351fd9873aSMichael Lange       else        {rootPoints[off+p].index = points[p]; rootPoints[off+p].rank = rank;}
15361fd9873aSMichael Lange     }
15371fd9873aSMichael Lange     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
15381fd9873aSMichael Lange     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
15391fd9873aSMichael Lange   }
15401fd9873aSMichael Lange   ierr = ISRestoreIndices(valueIS, &neighbors);CHKERRQ(ierr);
15411fd9873aSMichael Lange   ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
15421fd9873aSMichael Lange   /* Communicate overlap */
15431fd9873aSMichael Lange   ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr);
15441fd9873aSMichael Lange   ierr = DMPlexDistributeData(dm, processSF, rootSection, MPIU_2INT, rootPoints, leafSection, (void**) &leafPoints);CHKERRQ(ierr);
15451fd9873aSMichael Lange   /* Filter remote contributions (ovLeafPoints) into the overlapSF */
15461fd9873aSMichael Lange   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
15471fd9873aSMichael Lange   for (p = 0; p < size; p++) {
15481fd9873aSMichael Lange     ierr = DMLabelSetValue(leafLabel, leafPoints[p].index, leafPoints[p].rank);CHKERRQ(ierr);
15491fd9873aSMichael Lange   }
15501fd9873aSMichael Lange   ierr = PetscFree(rootPoints);CHKERRQ(ierr);
15511fd9873aSMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
15521fd9873aSMichael Lange   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
15531fd9873aSMichael Lange   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
15541fd9873aSMichael Lange   PetscFunctionReturn(0);
15551fd9873aSMichael Lange }
15561fd9873aSMichael Lange 
15571fd9873aSMichael Lange #undef __FUNCT__
1558aa3148a8SMichael Lange #define __FUNCT__ "DMPlexPartitionLabelCreateSF"
1559aa3148a8SMichael Lange /*@
1560aa3148a8SMichael Lange   DMPlexPartitionLabelCreateSF - Create a star forest from a label that assigns ranks to points
1561aa3148a8SMichael Lange 
1562aa3148a8SMichael Lange   Input Parameters:
1563aa3148a8SMichael Lange + dm    - The DM
1564aa3148a8SMichael Lange . label - DMLabel assinging ranks to remote roots
1565aa3148a8SMichael Lange 
1566aa3148a8SMichael Lange   Output Parameter:
1567aa3148a8SMichael Lange - sf    - The star forest communication context encapsulating the defined mapping
1568aa3148a8SMichael Lange 
1569aa3148a8SMichael Lange   Note: The incoming label is a receiver mapping of remote points to their parent rank.
1570aa3148a8SMichael Lange 
1571aa3148a8SMichael Lange   Level: developer
1572aa3148a8SMichael Lange 
1573aa3148a8SMichael Lange .seealso: DMPlexDistribute(), DMPlexCreateOverlap
1574aa3148a8SMichael Lange @*/
1575aa3148a8SMichael Lange PetscErrorCode DMPlexPartitionLabelCreateSF(DM dm, DMLabel label, PetscSF *sf)
1576aa3148a8SMichael Lange {
1577aa3148a8SMichael Lange   PetscMPIInt    numProcs;
1578aa3148a8SMichael Lange   PetscInt       n, idx, numRemote, p, numPoints, pStart, pEnd;
1579aa3148a8SMichael Lange   PetscSFNode   *remotePoints;
1580aa3148a8SMichael Lange   PetscErrorCode ierr;
1581aa3148a8SMichael Lange 
1582aa3148a8SMichael Lange   PetscFunctionBegin;
1583aa3148a8SMichael Lange   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr);
1584aa3148a8SMichael Lange 
1585aa3148a8SMichael Lange   for (numRemote = 0, n = 0; n < numProcs; ++n) {
1586aa3148a8SMichael Lange     ierr = DMLabelGetStratumSize(label, n, &numPoints);CHKERRQ(ierr);
1587aa3148a8SMichael Lange     numRemote += numPoints;
1588aa3148a8SMichael Lange   }
1589aa3148a8SMichael Lange   ierr = PetscMalloc1(numRemote, &remotePoints);CHKERRQ(ierr);
1590aa3148a8SMichael Lange   for (idx = 0, n = 0; n < numProcs; ++n) {
1591aa3148a8SMichael Lange     IS remoteRootIS;
1592aa3148a8SMichael Lange     const PetscInt *remoteRoots;
1593aa3148a8SMichael Lange     ierr = DMLabelGetStratumSize(label, n, &numPoints);CHKERRQ(ierr);
1594aa3148a8SMichael Lange     if (numPoints <= 0) continue;
1595aa3148a8SMichael Lange     ierr = DMLabelGetStratumIS(label, n, &remoteRootIS);CHKERRQ(ierr);
1596aa3148a8SMichael Lange     ierr = ISGetIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
1597aa3148a8SMichael Lange     for (p = 0; p < numPoints; p++) {
1598aa3148a8SMichael Lange       remotePoints[idx].index = remoteRoots[p];
1599aa3148a8SMichael Lange       remotePoints[idx].rank = n;
1600aa3148a8SMichael Lange       idx++;
1601aa3148a8SMichael Lange     }
1602aa3148a8SMichael Lange     ierr = ISRestoreIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
1603aa3148a8SMichael Lange     ierr = ISDestroy(&remoteRootIS);CHKERRQ(ierr);
1604aa3148a8SMichael Lange   }
1605aa3148a8SMichael Lange   ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), sf);CHKERRQ(ierr);
1606aa3148a8SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
1607aa3148a8SMichael Lange   ierr = PetscSFSetGraph(*sf, pEnd-pStart, numRemote, NULL, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
1608aa3148a8SMichael Lange   PetscFunctionReturn(0);
1609aa3148a8SMichael Lange }
1610