xref: /petsc/src/dm/impls/plex/plexpartition.c (revision ac9a96f1dae736d3e81a2c900ff2094167083dee)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
270034214SMatthew G. Knepley 
377623264SMatthew G. Knepley PetscClassId PETSCPARTITIONER_CLASSID = 0;
477623264SMatthew G. Knepley 
577623264SMatthew G. Knepley PetscFunctionList PetscPartitionerList              = NULL;
677623264SMatthew G. Knepley PetscBool         PetscPartitionerRegisterAllCalled = PETSC_FALSE;
777623264SMatthew G. Knepley 
877623264SMatthew G. Knepley PetscBool ChacoPartitionercite = PETSC_FALSE;
977623264SMatthew G. Knepley const char ChacoPartitionerCitation[] = "@inproceedings{Chaco95,\n"
1077623264SMatthew G. Knepley                                "  author    = {Bruce Hendrickson and Robert Leland},\n"
1177623264SMatthew G. Knepley                                "  title     = {A multilevel algorithm for partitioning graphs},\n"
1277623264SMatthew G. Knepley                                "  booktitle = {Supercomputing '95: Proceedings of the 1995 ACM/IEEE Conference on Supercomputing (CDROM)},"
1377623264SMatthew G. Knepley                                "  isbn      = {0-89791-816-9},\n"
1477623264SMatthew G. Knepley                                "  pages     = {28},\n"
1577623264SMatthew G. Knepley                                "  doi       = {http://doi.acm.org/10.1145/224170.224228},\n"
1677623264SMatthew G. Knepley                                "  publisher = {ACM Press},\n"
1777623264SMatthew G. Knepley                                "  address   = {New York},\n"
1877623264SMatthew G. Knepley                                "  year      = {1995}\n}\n";
1977623264SMatthew G. Knepley 
2077623264SMatthew G. Knepley PetscBool ParMetisPartitionercite = PETSC_FALSE;
2177623264SMatthew G. Knepley const char ParMetisPartitionerCitation[] = "@article{KarypisKumar98,\n"
2277623264SMatthew G. Knepley                                "  author  = {George Karypis and Vipin Kumar},\n"
2377623264SMatthew G. Knepley                                "  title   = {A Parallel Algorithm for Multilevel Graph Partitioning and Sparse Matrix Ordering},\n"
2477623264SMatthew G. Knepley                                "  journal = {Journal of Parallel and Distributed Computing},\n"
2577623264SMatthew G. Knepley                                "  volume  = {48},\n"
2677623264SMatthew G. Knepley                                "  pages   = {71--85},\n"
2777623264SMatthew G. Knepley                                "  year    = {1998}\n}\n";
2877623264SMatthew G. Knepley 
29532c4e7dSMichael Lange /*@C
30532c4e7dSMichael Lange   DMPlexCreatePartitionerGraph - Create a CSR graph of point connections for the partitioner
31532c4e7dSMichael Lange 
32532c4e7dSMichael Lange   Input Parameters:
33532c4e7dSMichael Lange + dm      - The mesh DM dm
34532c4e7dSMichael Lange - height  - Height of the strata from which to construct the graph
35532c4e7dSMichael Lange 
36532c4e7dSMichael Lange   Output Parameter:
37532c4e7dSMichael Lange + numVertices     - Number of vertices in the graph
383fa7752dSToby Isaac . offsets         - Point offsets in the graph
393fa7752dSToby Isaac . adjacency       - Point connectivity in the graph
403fa7752dSToby Isaac - globalNumbering - A map from the local cell numbering to the global numbering used in "adjacency".  Negative indicates that the cell is a duplicate from another process.
41532c4e7dSMichael Lange 
42532c4e7dSMichael Lange   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
43532c4e7dSMichael Lange   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
44532c4e7dSMichael Lange   representation on the mesh.
45532c4e7dSMichael Lange 
46532c4e7dSMichael Lange   Level: developer
47532c4e7dSMichael Lange 
48532c4e7dSMichael Lange .seealso: PetscPartitionerGetType(), PetscPartitionerCreate(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
49532c4e7dSMichael Lange @*/
503fa7752dSToby Isaac PetscErrorCode DMPlexCreatePartitionerGraph(DM dm, PetscInt height, PetscInt *numVertices, PetscInt **offsets, PetscInt **adjacency, IS *globalNumbering)
51532c4e7dSMichael Lange {
5243f7d02bSMichael Lange   PetscInt       p, pStart, pEnd, a, adjSize, idx, size, nroots;
53389e55d8SMichael Lange   PetscInt      *adj = NULL, *vOffsets = NULL, *graph = NULL;
548cfe4c1fSMichael Lange   IS             cellNumbering;
558cfe4c1fSMichael Lange   const PetscInt *cellNum;
56532c4e7dSMichael Lange   PetscBool      useCone, useClosure;
57532c4e7dSMichael Lange   PetscSection   section;
58532c4e7dSMichael Lange   PetscSegBuffer adjBuffer;
598cfe4c1fSMichael Lange   PetscSF        sfPoint;
60532c4e7dSMichael Lange   PetscMPIInt    rank;
61532c4e7dSMichael Lange   PetscErrorCode ierr;
62532c4e7dSMichael Lange 
63532c4e7dSMichael Lange   PetscFunctionBegin;
64532c4e7dSMichael Lange   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
65532c4e7dSMichael Lange   ierr = DMPlexGetHeightStratum(dm, height, &pStart, &pEnd);CHKERRQ(ierr);
668cfe4c1fSMichael Lange   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
678cfe4c1fSMichael Lange   ierr = PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
68532c4e7dSMichael Lange   /* Build adjacency graph via a section/segbuffer */
69532c4e7dSMichael Lange   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) dm), &section);CHKERRQ(ierr);
70532c4e7dSMichael Lange   ierr = PetscSectionSetChart(section, pStart, pEnd);CHKERRQ(ierr);
71532c4e7dSMichael Lange   ierr = PetscSegBufferCreate(sizeof(PetscInt),1000,&adjBuffer);CHKERRQ(ierr);
72532c4e7dSMichael Lange   /* Always use FVM adjacency to create partitioner graph */
73532c4e7dSMichael Lange   ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
74532c4e7dSMichael Lange   ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr);
75532c4e7dSMichael Lange   ierr = DMPlexSetAdjacencyUseCone(dm, PETSC_TRUE);CHKERRQ(ierr);
76532c4e7dSMichael Lange   ierr = DMPlexSetAdjacencyUseClosure(dm, PETSC_FALSE);CHKERRQ(ierr);
77f0927f4eSMatthew G. Knepley   ierr = DMPlexCreateCellNumbering_Internal(dm, PETSC_TRUE, &cellNumbering);CHKERRQ(ierr);
783fa7752dSToby Isaac   if (globalNumbering) {
793fa7752dSToby Isaac     ierr = PetscObjectReference((PetscObject)cellNumbering);CHKERRQ(ierr);
803fa7752dSToby Isaac     *globalNumbering = cellNumbering;
813fa7752dSToby Isaac   }
828cfe4c1fSMichael Lange   ierr = ISGetIndices(cellNumbering, &cellNum);CHKERRQ(ierr);
838cfe4c1fSMichael Lange   for (*numVertices = 0, p = pStart; p < pEnd; p++) {
848cfe4c1fSMichael Lange     /* Skip non-owned cells in parallel (ParMetis expects no overlap) */
858cfe4c1fSMichael Lange     if (nroots > 0) {if (cellNum[p] < 0) continue;}
86532c4e7dSMichael Lange     adjSize = PETSC_DETERMINE;
87532c4e7dSMichael Lange     ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
88532c4e7dSMichael Lange     for (a = 0; a < adjSize; ++a) {
89532c4e7dSMichael Lange       const PetscInt point = adj[a];
90532c4e7dSMichael Lange       if (point != p && pStart <= point && point < pEnd) {
91532c4e7dSMichael Lange         PetscInt *PETSC_RESTRICT pBuf;
92532c4e7dSMichael Lange         ierr = PetscSectionAddDof(section, p, 1);CHKERRQ(ierr);
93532c4e7dSMichael Lange         ierr = PetscSegBufferGetInts(adjBuffer, 1, &pBuf);CHKERRQ(ierr);
94532c4e7dSMichael Lange         *pBuf = point;
95532c4e7dSMichael Lange       }
96532c4e7dSMichael Lange     }
978cfe4c1fSMichael Lange     (*numVertices)++;
98532c4e7dSMichael Lange   }
99532c4e7dSMichael Lange   ierr = DMPlexSetAdjacencyUseCone(dm, useCone);CHKERRQ(ierr);
100532c4e7dSMichael Lange   ierr = DMPlexSetAdjacencyUseClosure(dm, useClosure);CHKERRQ(ierr);
101532c4e7dSMichael Lange   /* Derive CSR graph from section/segbuffer */
102532c4e7dSMichael Lange   ierr = PetscSectionSetUp(section);CHKERRQ(ierr);
103532c4e7dSMichael Lange   ierr = PetscSectionGetStorageSize(section, &size);CHKERRQ(ierr);
104389e55d8SMichael Lange   ierr = PetscMalloc1(*numVertices+1, &vOffsets);CHKERRQ(ierr);
10543f7d02bSMichael Lange   for (idx = 0, p = pStart; p < pEnd; p++) {
10643f7d02bSMichael Lange     if (nroots > 0) {if (cellNum[p] < 0) continue;}
10743f7d02bSMichael Lange     ierr = PetscSectionGetOffset(section, p, &(vOffsets[idx++]));CHKERRQ(ierr);
10843f7d02bSMichael Lange   }
109532c4e7dSMichael Lange   vOffsets[*numVertices] = size;
110532c4e7dSMichael Lange   if (offsets) *offsets = vOffsets;
111389e55d8SMichael Lange   ierr = PetscSegBufferExtractAlloc(adjBuffer, &graph);CHKERRQ(ierr);
112bf4602e4SToby Isaac   {
1138cfe4c1fSMichael Lange     ISLocalToGlobalMapping ltogCells;
1148cfe4c1fSMichael Lange     PetscInt n, size, *cells_arr;
1158cfe4c1fSMichael Lange     /* In parallel, apply a global cell numbering to the graph */
1168cfe4c1fSMichael Lange     ierr = ISRestoreIndices(cellNumbering, &cellNum);CHKERRQ(ierr);
1178cfe4c1fSMichael Lange     ierr = ISLocalToGlobalMappingCreateIS(cellNumbering, &ltogCells);CHKERRQ(ierr);
1188cfe4c1fSMichael Lange     ierr = ISLocalToGlobalMappingGetSize(ltogCells, &size);CHKERRQ(ierr);
1198cfe4c1fSMichael Lange     ierr = ISLocalToGlobalMappingGetIndices(ltogCells, (const PetscInt**)&cells_arr);CHKERRQ(ierr);
1208cfe4c1fSMichael Lange     /* Convert to positive global cell numbers */
1218cfe4c1fSMichael Lange     for (n=0; n<size; n++) {if (cells_arr[n] < 0) cells_arr[n] = -(cells_arr[n]+1);}
1228cfe4c1fSMichael Lange     ierr = ISLocalToGlobalMappingRestoreIndices(ltogCells, (const PetscInt**)&cells_arr);CHKERRQ(ierr);
123389e55d8SMichael Lange     ierr = ISLocalToGlobalMappingApplyBlock(ltogCells, vOffsets[*numVertices], graph, graph);CHKERRQ(ierr);
1248cfe4c1fSMichael Lange     ierr = ISLocalToGlobalMappingDestroy(&ltogCells);CHKERRQ(ierr);
125f0927f4eSMatthew G. Knepley     ierr = ISDestroy(&cellNumbering);CHKERRQ(ierr);
1268cfe4c1fSMichael Lange   }
127389e55d8SMichael Lange   if (adjacency) *adjacency = graph;
128532c4e7dSMichael Lange   /* Clean up */
129532c4e7dSMichael Lange   ierr = PetscSegBufferDestroy(&adjBuffer);CHKERRQ(ierr);
130532c4e7dSMichael Lange   ierr = PetscSectionDestroy(&section);CHKERRQ(ierr);
131532c4e7dSMichael Lange   ierr = PetscFree(adj);CHKERRQ(ierr);
132532c4e7dSMichael Lange   PetscFunctionReturn(0);
133532c4e7dSMichael Lange }
134532c4e7dSMichael Lange 
135d5577e40SMatthew G. Knepley /*@C
136d5577e40SMatthew G. Knepley   DMPlexCreateNeighborCSR - Create a mesh graph (cell-cell adjacency) in parallel CSR format.
137d5577e40SMatthew G. Knepley 
138d5577e40SMatthew G. Knepley   Collective
139d5577e40SMatthew G. Knepley 
140d5577e40SMatthew G. Knepley   Input Arguments:
141d5577e40SMatthew G. Knepley + dm - The DMPlex
142d5577e40SMatthew G. Knepley - cellHeight - The height of mesh points to treat as cells (default should be 0)
143d5577e40SMatthew G. Knepley 
144d5577e40SMatthew G. Knepley   Output Arguments:
145d5577e40SMatthew G. Knepley + numVertices - The number of local vertices in the graph, or cells in the mesh.
146d5577e40SMatthew G. Knepley . offsets     - The offset to the adjacency list for each cell
147d5577e40SMatthew G. Knepley - adjacency   - The adjacency list for all cells
148d5577e40SMatthew G. Knepley 
149d5577e40SMatthew G. Knepley   Note: This is suitable for input to a mesh partitioner like ParMetis.
150d5577e40SMatthew G. Knepley 
151d5577e40SMatthew G. Knepley   Level: advanced
152d5577e40SMatthew G. Knepley 
153d5577e40SMatthew G. Knepley .seealso: DMPlexCreate()
154d5577e40SMatthew G. Knepley @*/
15570034214SMatthew G. Knepley PetscErrorCode DMPlexCreateNeighborCSR(DM dm, PetscInt cellHeight, PetscInt *numVertices, PetscInt **offsets, PetscInt **adjacency)
15670034214SMatthew G. Knepley {
15770034214SMatthew G. Knepley   const PetscInt maxFaceCases = 30;
15870034214SMatthew G. Knepley   PetscInt       numFaceCases = 0;
15970034214SMatthew G. Knepley   PetscInt       numFaceVertices[30]; /* maxFaceCases, C89 sucks sucks sucks */
16070034214SMatthew G. Knepley   PetscInt      *off, *adj;
16170034214SMatthew G. Knepley   PetscInt      *neighborCells = NULL;
16270034214SMatthew G. Knepley   PetscInt       dim, cellDim, depth = 0, faceDepth, cStart, cEnd, c, numCells, cell;
16370034214SMatthew G. Knepley   PetscErrorCode ierr;
16470034214SMatthew G. Knepley 
16570034214SMatthew G. Knepley   PetscFunctionBegin;
16670034214SMatthew G. Knepley   /* For parallel partitioning, I think you have to communicate supports */
167c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
16870034214SMatthew G. Knepley   cellDim = dim - cellHeight;
16970034214SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
17070034214SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
17170034214SMatthew G. Knepley   if (cEnd - cStart == 0) {
17270034214SMatthew G. Knepley     if (numVertices) *numVertices = 0;
17370034214SMatthew G. Knepley     if (offsets)   *offsets   = NULL;
17470034214SMatthew G. Knepley     if (adjacency) *adjacency = NULL;
17570034214SMatthew G. Knepley     PetscFunctionReturn(0);
17670034214SMatthew G. Knepley   }
17770034214SMatthew G. Knepley   numCells  = cEnd - cStart;
17870034214SMatthew G. Knepley   faceDepth = depth - cellHeight;
17970034214SMatthew G. Knepley   if (dim == depth) {
18070034214SMatthew G. Knepley     PetscInt f, fStart, fEnd;
18170034214SMatthew G. Knepley 
18270034214SMatthew G. Knepley     ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr);
18370034214SMatthew G. Knepley     /* Count neighboring cells */
18470034214SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, cellHeight+1, &fStart, &fEnd);CHKERRQ(ierr);
18570034214SMatthew G. Knepley     for (f = fStart; f < fEnd; ++f) {
18670034214SMatthew G. Knepley       const PetscInt *support;
18770034214SMatthew G. Knepley       PetscInt        supportSize;
18870034214SMatthew G. Knepley       ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr);
18970034214SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr);
19070034214SMatthew G. Knepley       if (supportSize == 2) {
1919ffc88e4SToby Isaac         PetscInt numChildren;
1929ffc88e4SToby Isaac 
1939ffc88e4SToby Isaac         ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr);
1949ffc88e4SToby Isaac         if (!numChildren) {
19570034214SMatthew G. Knepley           ++off[support[0]-cStart+1];
19670034214SMatthew G. Knepley           ++off[support[1]-cStart+1];
19770034214SMatthew G. Knepley         }
19870034214SMatthew G. Knepley       }
1999ffc88e4SToby Isaac     }
20070034214SMatthew G. Knepley     /* Prefix sum */
20170034214SMatthew G. Knepley     for (c = 1; c <= numCells; ++c) off[c] += off[c-1];
20270034214SMatthew G. Knepley     if (adjacency) {
20370034214SMatthew G. Knepley       PetscInt *tmp;
20470034214SMatthew G. Knepley 
20570034214SMatthew G. Knepley       ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr);
206854ce69bSBarry Smith       ierr = PetscMalloc1(numCells+1, &tmp);CHKERRQ(ierr);
20770034214SMatthew G. Knepley       ierr = PetscMemcpy(tmp, off, (numCells+1) * sizeof(PetscInt));CHKERRQ(ierr);
20870034214SMatthew G. Knepley       /* Get neighboring cells */
20970034214SMatthew G. Knepley       for (f = fStart; f < fEnd; ++f) {
21070034214SMatthew G. Knepley         const PetscInt *support;
21170034214SMatthew G. Knepley         PetscInt        supportSize;
21270034214SMatthew G. Knepley         ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr);
21370034214SMatthew G. Knepley         ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr);
21470034214SMatthew G. Knepley         if (supportSize == 2) {
2159ffc88e4SToby Isaac           PetscInt numChildren;
2169ffc88e4SToby Isaac 
2179ffc88e4SToby Isaac           ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr);
2189ffc88e4SToby Isaac           if (!numChildren) {
21970034214SMatthew G. Knepley             adj[tmp[support[0]-cStart]++] = support[1];
22070034214SMatthew G. Knepley             adj[tmp[support[1]-cStart]++] = support[0];
22170034214SMatthew G. Knepley           }
22270034214SMatthew G. Knepley         }
2239ffc88e4SToby Isaac       }
22470034214SMatthew G. Knepley #if defined(PETSC_USE_DEBUG)
22570034214SMatthew 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);
22670034214SMatthew G. Knepley #endif
22770034214SMatthew G. Knepley       ierr = PetscFree(tmp);CHKERRQ(ierr);
22870034214SMatthew G. Knepley     }
22970034214SMatthew G. Knepley     if (numVertices) *numVertices = numCells;
23070034214SMatthew G. Knepley     if (offsets)   *offsets   = off;
23170034214SMatthew G. Knepley     if (adjacency) *adjacency = adj;
23270034214SMatthew G. Knepley     PetscFunctionReturn(0);
23370034214SMatthew G. Knepley   }
23470034214SMatthew G. Knepley   /* Setup face recognition */
23570034214SMatthew G. Knepley   if (faceDepth == 1) {
23670034214SMatthew 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 */
23770034214SMatthew G. Knepley 
23870034214SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
23970034214SMatthew G. Knepley       PetscInt corners;
24070034214SMatthew G. Knepley 
24170034214SMatthew G. Knepley       ierr = DMPlexGetConeSize(dm, c, &corners);CHKERRQ(ierr);
24270034214SMatthew G. Knepley       if (!cornersSeen[corners]) {
24370034214SMatthew G. Knepley         PetscInt nFV;
24470034214SMatthew G. Knepley 
24570034214SMatthew G. Knepley         if (numFaceCases >= maxFaceCases) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Exceeded maximum number of face recognition cases");
24670034214SMatthew G. Knepley         cornersSeen[corners] = 1;
24770034214SMatthew G. Knepley 
24870034214SMatthew G. Knepley         ierr = DMPlexGetNumFaceVertices(dm, cellDim, corners, &nFV);CHKERRQ(ierr);
24970034214SMatthew G. Knepley 
25070034214SMatthew G. Knepley         numFaceVertices[numFaceCases++] = nFV;
25170034214SMatthew G. Knepley       }
25270034214SMatthew G. Knepley     }
25370034214SMatthew G. Knepley   }
25470034214SMatthew G. Knepley   ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr);
25570034214SMatthew G. Knepley   /* Count neighboring cells */
25670034214SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) {
25770034214SMatthew G. Knepley     PetscInt numNeighbors = PETSC_DETERMINE, n;
25870034214SMatthew G. Knepley 
2598b0b4c70SToby Isaac     ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &numNeighbors, &neighborCells);CHKERRQ(ierr);
26070034214SMatthew G. Knepley     /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */
26170034214SMatthew G. Knepley     for (n = 0; n < numNeighbors; ++n) {
26270034214SMatthew G. Knepley       PetscInt        cellPair[2];
26370034214SMatthew G. Knepley       PetscBool       found    = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE;
26470034214SMatthew G. Knepley       PetscInt        meetSize = 0;
26570034214SMatthew G. Knepley       const PetscInt *meet    = NULL;
26670034214SMatthew G. Knepley 
26770034214SMatthew G. Knepley       cellPair[0] = cell; cellPair[1] = neighborCells[n];
26870034214SMatthew G. Knepley       if (cellPair[0] == cellPair[1]) continue;
26970034214SMatthew G. Knepley       if (!found) {
27070034214SMatthew G. Knepley         ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
27170034214SMatthew G. Knepley         if (meetSize) {
27270034214SMatthew G. Knepley           PetscInt f;
27370034214SMatthew G. Knepley 
27470034214SMatthew G. Knepley           for (f = 0; f < numFaceCases; ++f) {
27570034214SMatthew G. Knepley             if (numFaceVertices[f] == meetSize) {
27670034214SMatthew G. Knepley               found = PETSC_TRUE;
27770034214SMatthew G. Knepley               break;
27870034214SMatthew G. Knepley             }
27970034214SMatthew G. Knepley           }
28070034214SMatthew G. Knepley         }
28170034214SMatthew G. Knepley         ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
28270034214SMatthew G. Knepley       }
28370034214SMatthew G. Knepley       if (found) ++off[cell-cStart+1];
28470034214SMatthew G. Knepley     }
28570034214SMatthew G. Knepley   }
28670034214SMatthew G. Knepley   /* Prefix sum */
28770034214SMatthew G. Knepley   for (cell = 1; cell <= numCells; ++cell) off[cell] += off[cell-1];
28870034214SMatthew G. Knepley 
28970034214SMatthew G. Knepley   if (adjacency) {
29070034214SMatthew G. Knepley     ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr);
29170034214SMatthew G. Knepley     /* Get neighboring cells */
29270034214SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
29370034214SMatthew G. Knepley       PetscInt numNeighbors = PETSC_DETERMINE, n;
29470034214SMatthew G. Knepley       PetscInt cellOffset   = 0;
29570034214SMatthew G. Knepley 
2968b0b4c70SToby Isaac       ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &numNeighbors, &neighborCells);CHKERRQ(ierr);
29770034214SMatthew G. Knepley       /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */
29870034214SMatthew G. Knepley       for (n = 0; n < numNeighbors; ++n) {
29970034214SMatthew G. Knepley         PetscInt        cellPair[2];
30070034214SMatthew G. Knepley         PetscBool       found    = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE;
30170034214SMatthew G. Knepley         PetscInt        meetSize = 0;
30270034214SMatthew G. Knepley         const PetscInt *meet    = NULL;
30370034214SMatthew G. Knepley 
30470034214SMatthew G. Knepley         cellPair[0] = cell; cellPair[1] = neighborCells[n];
30570034214SMatthew G. Knepley         if (cellPair[0] == cellPair[1]) continue;
30670034214SMatthew G. Knepley         if (!found) {
30770034214SMatthew G. Knepley           ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
30870034214SMatthew G. Knepley           if (meetSize) {
30970034214SMatthew G. Knepley             PetscInt f;
31070034214SMatthew G. Knepley 
31170034214SMatthew G. Knepley             for (f = 0; f < numFaceCases; ++f) {
31270034214SMatthew G. Knepley               if (numFaceVertices[f] == meetSize) {
31370034214SMatthew G. Knepley                 found = PETSC_TRUE;
31470034214SMatthew G. Knepley                 break;
31570034214SMatthew G. Knepley               }
31670034214SMatthew G. Knepley             }
31770034214SMatthew G. Knepley           }
31870034214SMatthew G. Knepley           ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
31970034214SMatthew G. Knepley         }
32070034214SMatthew G. Knepley         if (found) {
32170034214SMatthew G. Knepley           adj[off[cell-cStart]+cellOffset] = neighborCells[n];
32270034214SMatthew G. Knepley           ++cellOffset;
32370034214SMatthew G. Knepley         }
32470034214SMatthew G. Knepley       }
32570034214SMatthew G. Knepley     }
32670034214SMatthew G. Knepley   }
32770034214SMatthew G. Knepley   ierr = PetscFree(neighborCells);CHKERRQ(ierr);
32870034214SMatthew G. Knepley   if (numVertices) *numVertices = numCells;
32970034214SMatthew G. Knepley   if (offsets)   *offsets   = off;
33070034214SMatthew G. Knepley   if (adjacency) *adjacency = adj;
33170034214SMatthew G. Knepley   PetscFunctionReturn(0);
33270034214SMatthew G. Knepley }
33370034214SMatthew G. Knepley 
33477623264SMatthew G. Knepley /*@C
33577623264SMatthew G. Knepley   PetscPartitionerRegister - Adds a new PetscPartitioner implementation
33677623264SMatthew G. Knepley 
33777623264SMatthew G. Knepley   Not Collective
33877623264SMatthew G. Knepley 
33977623264SMatthew G. Knepley   Input Parameters:
34077623264SMatthew G. Knepley + name        - The name of a new user-defined creation routine
34177623264SMatthew G. Knepley - create_func - The creation routine itself
34277623264SMatthew G. Knepley 
34377623264SMatthew G. Knepley   Notes:
34477623264SMatthew G. Knepley   PetscPartitionerRegister() may be called multiple times to add several user-defined PetscPartitioners
34577623264SMatthew G. Knepley 
34677623264SMatthew G. Knepley   Sample usage:
34777623264SMatthew G. Knepley .vb
34877623264SMatthew G. Knepley     PetscPartitionerRegister("my_part", MyPetscPartitionerCreate);
34977623264SMatthew G. Knepley .ve
35077623264SMatthew G. Knepley 
35177623264SMatthew G. Knepley   Then, your PetscPartitioner type can be chosen with the procedural interface via
35277623264SMatthew G. Knepley .vb
35377623264SMatthew G. Knepley     PetscPartitionerCreate(MPI_Comm, PetscPartitioner *);
35477623264SMatthew G. Knepley     PetscPartitionerSetType(PetscPartitioner, "my_part");
35577623264SMatthew G. Knepley .ve
35677623264SMatthew G. Knepley    or at runtime via the option
35777623264SMatthew G. Knepley .vb
35877623264SMatthew G. Knepley     -petscpartitioner_type my_part
35977623264SMatthew G. Knepley .ve
36077623264SMatthew G. Knepley 
36177623264SMatthew G. Knepley   Level: advanced
36277623264SMatthew G. Knepley 
36377623264SMatthew G. Knepley .keywords: PetscPartitioner, register
36477623264SMatthew G. Knepley .seealso: PetscPartitionerRegisterAll(), PetscPartitionerRegisterDestroy()
36577623264SMatthew G. Knepley 
36677623264SMatthew G. Knepley @*/
36777623264SMatthew G. Knepley PetscErrorCode PetscPartitionerRegister(const char sname[], PetscErrorCode (*function)(PetscPartitioner))
36877623264SMatthew G. Knepley {
36977623264SMatthew G. Knepley   PetscErrorCode ierr;
37077623264SMatthew G. Knepley 
37177623264SMatthew G. Knepley   PetscFunctionBegin;
37277623264SMatthew G. Knepley   ierr = PetscFunctionListAdd(&PetscPartitionerList, sname, function);CHKERRQ(ierr);
37377623264SMatthew G. Knepley   PetscFunctionReturn(0);
37477623264SMatthew G. Knepley }
37577623264SMatthew G. Knepley 
37677623264SMatthew G. Knepley /*@C
37777623264SMatthew G. Knepley   PetscPartitionerSetType - Builds a particular PetscPartitioner
37877623264SMatthew G. Knepley 
37977623264SMatthew G. Knepley   Collective on PetscPartitioner
38077623264SMatthew G. Knepley 
38177623264SMatthew G. Knepley   Input Parameters:
38277623264SMatthew G. Knepley + part - The PetscPartitioner object
38377623264SMatthew G. Knepley - name - The kind of partitioner
38477623264SMatthew G. Knepley 
38577623264SMatthew G. Knepley   Options Database Key:
38677623264SMatthew G. Knepley . -petscpartitioner_type <type> - Sets the PetscPartitioner type; use -help for a list of available types
38777623264SMatthew G. Knepley 
38877623264SMatthew G. Knepley   Level: intermediate
38977623264SMatthew G. Knepley 
39077623264SMatthew G. Knepley .keywords: PetscPartitioner, set, type
39177623264SMatthew G. Knepley .seealso: PetscPartitionerGetType(), PetscPartitionerCreate()
39277623264SMatthew G. Knepley @*/
39377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetType(PetscPartitioner part, PetscPartitionerType name)
39477623264SMatthew G. Knepley {
39577623264SMatthew G. Knepley   PetscErrorCode (*r)(PetscPartitioner);
39677623264SMatthew G. Knepley   PetscBool      match;
39777623264SMatthew G. Knepley   PetscErrorCode ierr;
39877623264SMatthew G. Knepley 
39977623264SMatthew G. Knepley   PetscFunctionBegin;
40077623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
40177623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) part, name, &match);CHKERRQ(ierr);
40277623264SMatthew G. Knepley   if (match) PetscFunctionReturn(0);
40377623264SMatthew G. Knepley 
4040f51fdf8SToby Isaac   ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr);
40577623264SMatthew G. Knepley   ierr = PetscFunctionListFind(PetscPartitionerList, name, &r);CHKERRQ(ierr);
40677623264SMatthew G. Knepley   if (!r) SETERRQ1(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscPartitioner type: %s", name);
40777623264SMatthew G. Knepley 
40877623264SMatthew G. Knepley   if (part->ops->destroy) {
40977623264SMatthew G. Knepley     ierr              = (*part->ops->destroy)(part);CHKERRQ(ierr);
41077623264SMatthew G. Knepley     part->ops->destroy = NULL;
41177623264SMatthew G. Knepley   }
41277623264SMatthew G. Knepley   ierr = (*r)(part);CHKERRQ(ierr);
41377623264SMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) part, name);CHKERRQ(ierr);
41477623264SMatthew G. Knepley   PetscFunctionReturn(0);
41577623264SMatthew G. Knepley }
41677623264SMatthew G. Knepley 
41777623264SMatthew G. Knepley /*@C
41877623264SMatthew G. Knepley   PetscPartitionerGetType - Gets the PetscPartitioner type name (as a string) from the object.
41977623264SMatthew G. Knepley 
42077623264SMatthew G. Knepley   Not Collective
42177623264SMatthew G. Knepley 
42277623264SMatthew G. Knepley   Input Parameter:
42377623264SMatthew G. Knepley . part - The PetscPartitioner
42477623264SMatthew G. Knepley 
42577623264SMatthew G. Knepley   Output Parameter:
42677623264SMatthew G. Knepley . name - The PetscPartitioner type name
42777623264SMatthew G. Knepley 
42877623264SMatthew G. Knepley   Level: intermediate
42977623264SMatthew G. Knepley 
43077623264SMatthew G. Knepley .keywords: PetscPartitioner, get, type, name
43177623264SMatthew G. Knepley .seealso: PetscPartitionerSetType(), PetscPartitionerCreate()
43277623264SMatthew G. Knepley @*/
43377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerGetType(PetscPartitioner part, PetscPartitionerType *name)
43477623264SMatthew G. Knepley {
43577623264SMatthew G. Knepley   PetscErrorCode ierr;
43677623264SMatthew G. Knepley 
43777623264SMatthew G. Knepley   PetscFunctionBegin;
43877623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
43977623264SMatthew G. Knepley   PetscValidPointer(name, 2);
4400f51fdf8SToby Isaac   ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr);
44177623264SMatthew G. Knepley   *name = ((PetscObject) part)->type_name;
44277623264SMatthew G. Knepley   PetscFunctionReturn(0);
44377623264SMatthew G. Knepley }
44477623264SMatthew G. Knepley 
44577623264SMatthew G. Knepley /*@C
44677623264SMatthew G. Knepley   PetscPartitionerView - Views a PetscPartitioner
44777623264SMatthew G. Knepley 
44877623264SMatthew G. Knepley   Collective on PetscPartitioner
44977623264SMatthew G. Knepley 
45077623264SMatthew G. Knepley   Input Parameter:
45177623264SMatthew G. Knepley + part - the PetscPartitioner object to view
45277623264SMatthew G. Knepley - v    - the viewer
45377623264SMatthew G. Knepley 
45477623264SMatthew G. Knepley   Level: developer
45577623264SMatthew G. Knepley 
45677623264SMatthew G. Knepley .seealso: PetscPartitionerDestroy()
45777623264SMatthew G. Knepley @*/
45877623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView(PetscPartitioner part, PetscViewer v)
45977623264SMatthew G. Knepley {
46077623264SMatthew G. Knepley   PetscErrorCode ierr;
46177623264SMatthew G. Knepley 
46277623264SMatthew G. Knepley   PetscFunctionBegin;
46377623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
46477623264SMatthew G. Knepley   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) part), &v);CHKERRQ(ierr);}
46577623264SMatthew G. Knepley   if (part->ops->view) {ierr = (*part->ops->view)(part, v);CHKERRQ(ierr);}
46677623264SMatthew G. Knepley   PetscFunctionReturn(0);
46777623264SMatthew G. Knepley }
46877623264SMatthew G. Knepley 
46977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetTypeFromOptions_Internal(PetscPartitioner part)
47077623264SMatthew G. Knepley {
47177623264SMatthew G. Knepley   const char    *defaultType;
47277623264SMatthew G. Knepley   char           name[256];
47377623264SMatthew G. Knepley   PetscBool      flg;
47477623264SMatthew G. Knepley   PetscErrorCode ierr;
47577623264SMatthew G. Knepley 
47677623264SMatthew G. Knepley   PetscFunctionBegin;
47777623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
4787cf71cfdSLisandro Dalcin   if (!((PetscObject) part)->type_name)
4797cf71cfdSLisandro Dalcin #if defined(PETSC_HAVE_CHACO)
4807cf71cfdSLisandro Dalcin     defaultType = PETSCPARTITIONERCHACO;
4817cf71cfdSLisandro Dalcin #elif defined(PETSC_HAVE_PARMETIS)
4827cf71cfdSLisandro Dalcin     defaultType = PETSCPARTITIONERPARMETIS;
4837cf71cfdSLisandro Dalcin #else
4847cf71cfdSLisandro Dalcin     defaultType = PETSCPARTITIONERSIMPLE;
4857cf71cfdSLisandro Dalcin #endif
4867cf71cfdSLisandro Dalcin   else
4877cf71cfdSLisandro Dalcin     defaultType = ((PetscObject) part)->type_name;
4880f51fdf8SToby Isaac   ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr);
48977623264SMatthew G. Knepley 
49077623264SMatthew G. Knepley   ierr = PetscObjectOptionsBegin((PetscObject) part);CHKERRQ(ierr);
49177623264SMatthew G. Knepley   ierr = PetscOptionsFList("-petscpartitioner_type", "Graph partitioner", "PetscPartitionerSetType", PetscPartitionerList, defaultType, name, 256, &flg);CHKERRQ(ierr);
49277623264SMatthew G. Knepley   if (flg) {
49377623264SMatthew G. Knepley     ierr = PetscPartitionerSetType(part, name);CHKERRQ(ierr);
49477623264SMatthew G. Knepley   } else if (!((PetscObject) part)->type_name) {
49577623264SMatthew G. Knepley     ierr = PetscPartitionerSetType(part, defaultType);CHKERRQ(ierr);
49677623264SMatthew G. Knepley   }
49777623264SMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
49877623264SMatthew G. Knepley   PetscFunctionReturn(0);
49977623264SMatthew G. Knepley }
50077623264SMatthew G. Knepley 
50177623264SMatthew G. Knepley /*@
50277623264SMatthew G. Knepley   PetscPartitionerSetFromOptions - sets parameters in a PetscPartitioner from the options database
50377623264SMatthew G. Knepley 
50477623264SMatthew G. Knepley   Collective on PetscPartitioner
50577623264SMatthew G. Knepley 
50677623264SMatthew G. Knepley   Input Parameter:
50777623264SMatthew G. Knepley . part - the PetscPartitioner object to set options for
50877623264SMatthew G. Knepley 
50977623264SMatthew G. Knepley   Level: developer
51077623264SMatthew G. Knepley 
51177623264SMatthew G. Knepley .seealso: PetscPartitionerView()
51277623264SMatthew G. Knepley @*/
51377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetFromOptions(PetscPartitioner part)
51477623264SMatthew G. Knepley {
51577623264SMatthew G. Knepley   PetscErrorCode ierr;
51677623264SMatthew G. Knepley 
51777623264SMatthew G. Knepley   PetscFunctionBegin;
51877623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
51977623264SMatthew G. Knepley   ierr = PetscPartitionerSetTypeFromOptions_Internal(part);CHKERRQ(ierr);
52077623264SMatthew G. Knepley 
52177623264SMatthew G. Knepley   ierr = PetscObjectOptionsBegin((PetscObject) part);CHKERRQ(ierr);
522077101c0SMatthew G. Knepley   if (part->ops->setfromoptions) {ierr = (*part->ops->setfromoptions)(PetscOptionsObject,part);CHKERRQ(ierr);}
52377623264SMatthew G. Knepley   /* process any options handlers added with PetscObjectAddOptionsHandler() */
5240633abcbSJed Brown   ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) part);CHKERRQ(ierr);
52577623264SMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
52677623264SMatthew G. Knepley   ierr = PetscPartitionerViewFromOptions(part, NULL, "-petscpartitioner_view");CHKERRQ(ierr);
52777623264SMatthew G. Knepley   PetscFunctionReturn(0);
52877623264SMatthew G. Knepley }
52977623264SMatthew G. Knepley 
53077623264SMatthew G. Knepley /*@C
53177623264SMatthew G. Knepley   PetscPartitionerSetUp - Construct data structures for the PetscPartitioner
53277623264SMatthew G. Knepley 
53377623264SMatthew G. Knepley   Collective on PetscPartitioner
53477623264SMatthew G. Knepley 
53577623264SMatthew G. Knepley   Input Parameter:
53677623264SMatthew G. Knepley . part - the PetscPartitioner object to setup
53777623264SMatthew G. Knepley 
53877623264SMatthew G. Knepley   Level: developer
53977623264SMatthew G. Knepley 
54077623264SMatthew G. Knepley .seealso: PetscPartitionerView(), PetscPartitionerDestroy()
54177623264SMatthew G. Knepley @*/
54277623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetUp(PetscPartitioner part)
54377623264SMatthew G. Knepley {
54477623264SMatthew G. Knepley   PetscErrorCode ierr;
54577623264SMatthew G. Knepley 
54677623264SMatthew G. Knepley   PetscFunctionBegin;
54777623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
54877623264SMatthew G. Knepley   if (part->ops->setup) {ierr = (*part->ops->setup)(part);CHKERRQ(ierr);}
54977623264SMatthew G. Knepley   PetscFunctionReturn(0);
55077623264SMatthew G. Knepley }
55177623264SMatthew G. Knepley 
55277623264SMatthew G. Knepley /*@
55377623264SMatthew G. Knepley   PetscPartitionerDestroy - Destroys a PetscPartitioner object
55477623264SMatthew G. Knepley 
55577623264SMatthew G. Knepley   Collective on PetscPartitioner
55677623264SMatthew G. Knepley 
55777623264SMatthew G. Knepley   Input Parameter:
55877623264SMatthew G. Knepley . part - the PetscPartitioner object to destroy
55977623264SMatthew G. Knepley 
56077623264SMatthew G. Knepley   Level: developer
56177623264SMatthew G. Knepley 
56277623264SMatthew G. Knepley .seealso: PetscPartitionerView()
56377623264SMatthew G. Knepley @*/
56477623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy(PetscPartitioner *part)
56577623264SMatthew G. Knepley {
56677623264SMatthew G. Knepley   PetscErrorCode ierr;
56777623264SMatthew G. Knepley 
56877623264SMatthew G. Knepley   PetscFunctionBegin;
56977623264SMatthew G. Knepley   if (!*part) PetscFunctionReturn(0);
57077623264SMatthew G. Knepley   PetscValidHeaderSpecific((*part), PETSCPARTITIONER_CLASSID, 1);
57177623264SMatthew G. Knepley 
57277623264SMatthew G. Knepley   if (--((PetscObject)(*part))->refct > 0) {*part = 0; PetscFunctionReturn(0);}
57377623264SMatthew G. Knepley   ((PetscObject) (*part))->refct = 0;
57477623264SMatthew G. Knepley 
57577623264SMatthew G. Knepley   if ((*part)->ops->destroy) {ierr = (*(*part)->ops->destroy)(*part);CHKERRQ(ierr);}
57677623264SMatthew G. Knepley   ierr = PetscHeaderDestroy(part);CHKERRQ(ierr);
57777623264SMatthew G. Knepley   PetscFunctionReturn(0);
57877623264SMatthew G. Knepley }
57977623264SMatthew G. Knepley 
58077623264SMatthew G. Knepley /*@
58177623264SMatthew G. Knepley   PetscPartitionerCreate - Creates an empty PetscPartitioner object. The type can then be set with PetscPartitionerSetType().
58277623264SMatthew G. Knepley 
58377623264SMatthew G. Knepley   Collective on MPI_Comm
58477623264SMatthew G. Knepley 
58577623264SMatthew G. Knepley   Input Parameter:
58677623264SMatthew G. Knepley . comm - The communicator for the PetscPartitioner object
58777623264SMatthew G. Knepley 
58877623264SMatthew G. Knepley   Output Parameter:
58977623264SMatthew G. Knepley . part - The PetscPartitioner object
59077623264SMatthew G. Knepley 
59177623264SMatthew G. Knepley   Level: beginner
59277623264SMatthew G. Knepley 
593dae52e14SToby Isaac .seealso: PetscPartitionerSetType(), PETSCPARTITIONERCHACO, PETSCPARTITIONERPARMETIS, PETSCPARTITIONERSHELL, PETSCPARTITIONERSIMPLE, PETSCPARTITIONERGATHER
59477623264SMatthew G. Knepley @*/
59577623264SMatthew G. Knepley PetscErrorCode PetscPartitionerCreate(MPI_Comm comm, PetscPartitioner *part)
59677623264SMatthew G. Knepley {
59777623264SMatthew G. Knepley   PetscPartitioner p;
59877623264SMatthew G. Knepley   PetscErrorCode   ierr;
59977623264SMatthew G. Knepley 
60077623264SMatthew G. Knepley   PetscFunctionBegin;
60177623264SMatthew G. Knepley   PetscValidPointer(part, 2);
60277623264SMatthew G. Knepley   *part = NULL;
60383cde681SMatthew G. Knepley   ierr = DMInitializePackage();CHKERRQ(ierr);
60477623264SMatthew G. Knepley 
60573107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(p, PETSCPARTITIONER_CLASSID, "PetscPartitioner", "Graph Partitioner", "PetscPartitioner", comm, PetscPartitionerDestroy, PetscPartitionerView);CHKERRQ(ierr);
60677623264SMatthew G. Knepley 
60777623264SMatthew G. Knepley   *part = p;
60877623264SMatthew G. Knepley   PetscFunctionReturn(0);
60977623264SMatthew G. Knepley }
61077623264SMatthew G. Knepley 
61177623264SMatthew G. Knepley /*@
61277623264SMatthew G. Knepley   PetscPartitionerPartition - Create a non-overlapping partition of the cells in the mesh
61377623264SMatthew G. Knepley 
61477623264SMatthew G. Knepley   Collective on DM
61577623264SMatthew G. Knepley 
61677623264SMatthew G. Knepley   Input Parameters:
61777623264SMatthew G. Knepley + part    - The PetscPartitioner
618f8987ae8SMichael Lange - dm      - The mesh DM
61977623264SMatthew G. Knepley 
62077623264SMatthew G. Knepley   Output Parameters:
62177623264SMatthew G. Knepley + partSection     - The PetscSection giving the division of points by partition
622f8987ae8SMichael Lange - partition       - The list of points by partition
62377623264SMatthew G. Knepley 
62477623264SMatthew G. Knepley   Note: Instead of cells, points at a given height can be partitioned by calling PetscPartitionerSetPointHeight()
62577623264SMatthew G. Knepley 
62677623264SMatthew G. Knepley   Level: developer
62777623264SMatthew G. Knepley 
62877623264SMatthew G. Knepley .seealso DMPlexDistribute(), PetscPartitionerSetPointHeight(), PetscPartitionerCreate()
6294b15ede2SMatthew G. Knepley @*/
630f8987ae8SMichael Lange PetscErrorCode PetscPartitionerPartition(PetscPartitioner part, DM dm, PetscSection partSection, IS *partition)
63177623264SMatthew G. Knepley {
63277623264SMatthew G. Knepley   PetscMPIInt    size;
63377623264SMatthew G. Knepley   PetscErrorCode ierr;
63477623264SMatthew G. Knepley 
63577623264SMatthew G. Knepley   PetscFunctionBegin;
63677623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
63777623264SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
63877623264SMatthew G. Knepley   PetscValidHeaderSpecific(partSection, PETSC_SECTION_CLASSID, 4);
63977623264SMatthew G. Knepley   PetscValidPointer(partition, 5);
64077623264SMatthew G. Knepley   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) part), &size);CHKERRQ(ierr);
64177623264SMatthew G. Knepley   if (size == 1) {
64277623264SMatthew G. Knepley     PetscInt *points;
64377623264SMatthew G. Knepley     PetscInt  cStart, cEnd, c;
64477623264SMatthew G. Knepley 
64577623264SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, part->height, &cStart, &cEnd);CHKERRQ(ierr);
64677623264SMatthew G. Knepley     ierr = PetscSectionSetChart(partSection, 0, size);CHKERRQ(ierr);
64777623264SMatthew G. Knepley     ierr = PetscSectionSetDof(partSection, 0, cEnd-cStart);CHKERRQ(ierr);
64877623264SMatthew G. Knepley     ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
64977623264SMatthew G. Knepley     ierr = PetscMalloc1(cEnd-cStart, &points);CHKERRQ(ierr);
65077623264SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) points[c] = c;
65177623264SMatthew G. Knepley     ierr = ISCreateGeneral(PetscObjectComm((PetscObject) part), cEnd-cStart, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
65277623264SMatthew G. Knepley     PetscFunctionReturn(0);
65377623264SMatthew G. Knepley   }
65413fef8f3SMatthew G. Knepley   /* Obvious cheating, but I cannot think of a better way to do this. The DMSetFromOptions() has refinement in it,
65513fef8f3SMatthew G. Knepley      so we cannot call it and have it feed down to the partitioner before partitioning */
65613fef8f3SMatthew G. Knepley   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
65777623264SMatthew G. Knepley   if (part->height == 0) {
65877623264SMatthew G. Knepley     PetscInt numVertices;
65977623264SMatthew G. Knepley     PetscInt *start     = NULL;
66077623264SMatthew G. Knepley     PetscInt *adjacency = NULL;
6613fa7752dSToby Isaac     IS       globalNumbering;
66277623264SMatthew G. Knepley 
6633fa7752dSToby Isaac     ierr = DMPlexCreatePartitionerGraph(dm, 0, &numVertices, &start, &adjacency, &globalNumbering);CHKERRQ(ierr);
66477623264SMatthew G. Knepley     if (!part->ops->partition) SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_WRONGSTATE, "PetscPartitioner has no type");
66577623264SMatthew G. Knepley     ierr = (*part->ops->partition)(part, dm, size, numVertices, start, adjacency, partSection, partition);CHKERRQ(ierr);
66677623264SMatthew G. Knepley     ierr = PetscFree(start);CHKERRQ(ierr);
66777623264SMatthew G. Knepley     ierr = PetscFree(adjacency);CHKERRQ(ierr);
6683fa7752dSToby Isaac     if (globalNumbering) { /* partition is wrt global unique numbering: change this to be wrt local numbering */
6693fa7752dSToby Isaac       const PetscInt *globalNum;
6703fa7752dSToby Isaac       const PetscInt *partIdx;
6713fa7752dSToby Isaac       PetscInt *map, cStart, cEnd;
6723fa7752dSToby Isaac       PetscInt *adjusted, i, localSize, offset;
6733fa7752dSToby Isaac       IS    newPartition;
6743fa7752dSToby Isaac 
6753fa7752dSToby Isaac       ierr = ISGetLocalSize(*partition,&localSize);CHKERRQ(ierr);
6763fa7752dSToby Isaac       ierr = PetscMalloc1(localSize,&adjusted);CHKERRQ(ierr);
6773fa7752dSToby Isaac       ierr = ISGetIndices(globalNumbering,&globalNum);CHKERRQ(ierr);
6783fa7752dSToby Isaac       ierr = ISGetIndices(*partition,&partIdx);CHKERRQ(ierr);
6793fa7752dSToby Isaac       ierr = PetscMalloc1(localSize,&map);CHKERRQ(ierr);
6803fa7752dSToby Isaac       ierr = DMPlexGetHeightStratum(dm, part->height, &cStart, &cEnd);CHKERRQ(ierr);
6813fa7752dSToby Isaac       for (i = cStart, offset = 0; i < cEnd; i++) {
6823fa7752dSToby Isaac         if (globalNum[i - cStart] >= 0) map[offset++] = i;
6833fa7752dSToby Isaac       }
6843fa7752dSToby Isaac       for (i = 0; i < localSize; i++) {
6853fa7752dSToby Isaac         adjusted[i] = map[partIdx[i]];
6863fa7752dSToby Isaac       }
6873fa7752dSToby Isaac       ierr = PetscFree(map);CHKERRQ(ierr);
6883fa7752dSToby Isaac       ierr = ISRestoreIndices(*partition,&partIdx);CHKERRQ(ierr);
6893fa7752dSToby Isaac       ierr = ISRestoreIndices(globalNumbering,&globalNum);CHKERRQ(ierr);
6903fa7752dSToby Isaac       ierr = ISCreateGeneral(PETSC_COMM_SELF,localSize,adjusted,PETSC_OWN_POINTER,&newPartition);CHKERRQ(ierr);
6913fa7752dSToby Isaac       ierr = ISDestroy(&globalNumbering);CHKERRQ(ierr);
6923fa7752dSToby Isaac       ierr = ISDestroy(partition);CHKERRQ(ierr);
6933fa7752dSToby Isaac       *partition = newPartition;
6943fa7752dSToby Isaac     }
69577623264SMatthew G. Knepley   } else SETERRQ1(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_OUTOFRANGE, "Invalid height %D for points to partition", part->height);
69677623264SMatthew G. Knepley   PetscFunctionReturn(0);
6973fa7752dSToby Isaac 
69877623264SMatthew G. Knepley }
69977623264SMatthew G. Knepley 
700d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerDestroy_Shell(PetscPartitioner part)
70177623264SMatthew G. Knepley {
70277623264SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
70377623264SMatthew G. Knepley   PetscErrorCode          ierr;
70477623264SMatthew G. Knepley 
70577623264SMatthew G. Knepley   PetscFunctionBegin;
70677623264SMatthew G. Knepley   ierr = PetscSectionDestroy(&p->section);CHKERRQ(ierr);
70777623264SMatthew G. Knepley   ierr = ISDestroy(&p->partition);CHKERRQ(ierr);
70877623264SMatthew G. Knepley   ierr = PetscFree(p);CHKERRQ(ierr);
70977623264SMatthew G. Knepley   PetscFunctionReturn(0);
71077623264SMatthew G. Knepley }
71177623264SMatthew G. Knepley 
712d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Shell_Ascii(PetscPartitioner part, PetscViewer viewer)
71377623264SMatthew G. Knepley {
714077101c0SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
71577623264SMatthew G. Knepley   PetscViewerFormat       format;
71677623264SMatthew G. Knepley   PetscErrorCode          ierr;
71777623264SMatthew G. Knepley 
71877623264SMatthew G. Knepley   PetscFunctionBegin;
71977623264SMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
72077623264SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Shell Graph Partitioner:\n");CHKERRQ(ierr);
721077101c0SMatthew G. Knepley   if (p->random) {
722077101c0SMatthew G. Knepley     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
723077101c0SMatthew G. Knepley     ierr = PetscViewerASCIIPrintf(viewer, "using random partition\n");CHKERRQ(ierr);
724077101c0SMatthew G. Knepley     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
725077101c0SMatthew G. Knepley   }
72677623264SMatthew G. Knepley   PetscFunctionReturn(0);
72777623264SMatthew G. Knepley }
72877623264SMatthew G. Knepley 
729d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Shell(PetscPartitioner part, PetscViewer viewer)
73077623264SMatthew G. Knepley {
73177623264SMatthew G. Knepley   PetscBool      iascii;
73277623264SMatthew G. Knepley   PetscErrorCode ierr;
73377623264SMatthew G. Knepley 
73477623264SMatthew G. Knepley   PetscFunctionBegin;
73577623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
73677623264SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
73777623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
73877623264SMatthew G. Knepley   if (iascii) {ierr = PetscPartitionerView_Shell_Ascii(part, viewer);CHKERRQ(ierr);}
73977623264SMatthew G. Knepley   PetscFunctionReturn(0);
74077623264SMatthew G. Knepley }
74177623264SMatthew G. Knepley 
742d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerSetFromOptions_Shell(PetscOptionItems *PetscOptionsObject, PetscPartitioner part)
743077101c0SMatthew G. Knepley {
744077101c0SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
745077101c0SMatthew G. Knepley   PetscErrorCode          ierr;
746077101c0SMatthew G. Knepley 
747077101c0SMatthew G. Knepley   PetscFunctionBegin;
748077101c0SMatthew G. Knepley   ierr = PetscOptionsHead(PetscOptionsObject, "PetscPartitionerShell Options");CHKERRQ(ierr);
749077101c0SMatthew G. Knepley   ierr = PetscOptionsBool("-petscpartitioner_shell_random", "Use a random partition", "PetscPartitionerView", PETSC_FALSE, &p->random, NULL);CHKERRQ(ierr);
750077101c0SMatthew G. Knepley   ierr = PetscOptionsTail();CHKERRQ(ierr);
751077101c0SMatthew G. Knepley   PetscFunctionReturn(0);
752077101c0SMatthew G. Knepley }
753077101c0SMatthew G. Knepley 
754d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_Shell(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
75577623264SMatthew G. Knepley {
75677623264SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
75777623264SMatthew G. Knepley   PetscInt                np;
75877623264SMatthew G. Knepley   PetscErrorCode          ierr;
75977623264SMatthew G. Knepley 
76077623264SMatthew G. Knepley   PetscFunctionBegin;
761077101c0SMatthew G. Knepley   if (p->random) {
762077101c0SMatthew G. Knepley     PetscRandom r;
763aa1d5631SMatthew G. Knepley     PetscInt   *sizes, *points, v, p;
764aa1d5631SMatthew G. Knepley     PetscMPIInt rank;
765077101c0SMatthew G. Knepley 
766aa1d5631SMatthew G. Knepley     ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
767077101c0SMatthew G. Knepley     ierr = PetscRandomCreate(PETSC_COMM_SELF, &r);CHKERRQ(ierr);
768c717d290SMatthew G. Knepley     ierr = PetscRandomSetInterval(r, 0.0, (PetscScalar) nparts);CHKERRQ(ierr);
769077101c0SMatthew G. Knepley     ierr = PetscRandomSetFromOptions(r);CHKERRQ(ierr);
770077101c0SMatthew G. Knepley     ierr = PetscCalloc2(nparts, &sizes, numVertices, &points);CHKERRQ(ierr);
771aa1d5631SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {points[v] = v;}
772*ac9a96f1SMichael Lange     for (p = 0; p < nparts; ++p) {sizes[p] = numVertices/nparts + (PetscInt) (p < numVertices % nparts);}
773aa1d5631SMatthew G. Knepley     for (v = numVertices-1; v > 0; --v) {
774077101c0SMatthew G. Knepley       PetscReal val;
775aa1d5631SMatthew G. Knepley       PetscInt  w, tmp;
776077101c0SMatthew G. Knepley 
777aa1d5631SMatthew G. Knepley       ierr = PetscRandomSetInterval(r, 0.0, (PetscScalar) (v+1));CHKERRQ(ierr);
778077101c0SMatthew G. Knepley       ierr = PetscRandomGetValueReal(r, &val);CHKERRQ(ierr);
779aa1d5631SMatthew G. Knepley       w    = PetscFloorReal(val);
780aa1d5631SMatthew G. Knepley       tmp       = points[v];
781aa1d5631SMatthew G. Knepley       points[v] = points[w];
782aa1d5631SMatthew G. Knepley       points[w] = tmp;
783077101c0SMatthew G. Knepley     }
784077101c0SMatthew G. Knepley     ierr = PetscRandomDestroy(&r);CHKERRQ(ierr);
785077101c0SMatthew G. Knepley     ierr = PetscPartitionerShellSetPartition(part, nparts, sizes, points);CHKERRQ(ierr);
786077101c0SMatthew G. Knepley     ierr = PetscFree2(sizes, points);CHKERRQ(ierr);
787077101c0SMatthew G. Knepley   }
788c717d290SMatthew G. Knepley   if (!p->section) SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_WRONG, "Shell partitioner information not provided. Please call PetscPartitionerShellSetPartition()");
78977623264SMatthew G. Knepley   ierr = PetscSectionGetChart(p->section, NULL, &np);CHKERRQ(ierr);
79077623264SMatthew G. Knepley   if (nparts != np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of requested partitions %d != configured partitions %d", nparts, np);
79177623264SMatthew G. Knepley   ierr = ISGetLocalSize(p->partition, &np);CHKERRQ(ierr);
79277623264SMatthew G. Knepley   if (numVertices != np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of input vertices %d != configured vertices %d", numVertices, np);
7935680f57bSMatthew G. Knepley   ierr = PetscSectionCopy(p->section, partSection);CHKERRQ(ierr);
79477623264SMatthew G. Knepley   *partition = p->partition;
79577623264SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) p->partition);CHKERRQ(ierr);
79677623264SMatthew G. Knepley   PetscFunctionReturn(0);
79777623264SMatthew G. Knepley }
79877623264SMatthew G. Knepley 
799d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerInitialize_Shell(PetscPartitioner part)
80077623264SMatthew G. Knepley {
80177623264SMatthew G. Knepley   PetscFunctionBegin;
80277623264SMatthew G. Knepley   part->ops->view           = PetscPartitionerView_Shell;
803077101c0SMatthew G. Knepley   part->ops->setfromoptions = PetscPartitionerSetFromOptions_Shell;
80477623264SMatthew G. Knepley   part->ops->destroy        = PetscPartitionerDestroy_Shell;
80577623264SMatthew G. Knepley   part->ops->partition      = PetscPartitionerPartition_Shell;
80677623264SMatthew G. Knepley   PetscFunctionReturn(0);
80777623264SMatthew G. Knepley }
80877623264SMatthew G. Knepley 
80977623264SMatthew G. Knepley /*MC
81077623264SMatthew G. Knepley   PETSCPARTITIONERSHELL = "shell" - A PetscPartitioner object
81177623264SMatthew G. Knepley 
81277623264SMatthew G. Knepley   Level: intermediate
81377623264SMatthew G. Knepley 
81477623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
81577623264SMatthew G. Knepley M*/
81677623264SMatthew G. Knepley 
81777623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Shell(PetscPartitioner part)
81877623264SMatthew G. Knepley {
81977623264SMatthew G. Knepley   PetscPartitioner_Shell *p;
82077623264SMatthew G. Knepley   PetscErrorCode          ierr;
82177623264SMatthew G. Knepley 
82277623264SMatthew G. Knepley   PetscFunctionBegin;
82377623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
82477623264SMatthew G. Knepley   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
82577623264SMatthew G. Knepley   part->data = p;
82677623264SMatthew G. Knepley 
82777623264SMatthew G. Knepley   ierr = PetscPartitionerInitialize_Shell(part);CHKERRQ(ierr);
828077101c0SMatthew G. Knepley   p->random = PETSC_FALSE;
82977623264SMatthew G. Knepley   PetscFunctionReturn(0);
83077623264SMatthew G. Knepley }
83177623264SMatthew G. Knepley 
8325680f57bSMatthew G. Knepley /*@C
8335680f57bSMatthew G. Knepley   PetscPartitionerShellSetPartition - Set an artifical partition for a mesh
8345680f57bSMatthew G. Knepley 
835077101c0SMatthew G. Knepley   Collective on Part
8365680f57bSMatthew G. Knepley 
8375680f57bSMatthew G. Knepley   Input Parameters:
8385680f57bSMatthew G. Knepley + part     - The PetscPartitioner
8399852e123SBarry Smith . size - The number of partitions
8409852e123SBarry Smith . sizes    - array of size size (or NULL) providing the number of points in each partition
8419758bf69SToby Isaac - points   - array of size sum(sizes) (may be NULL iff sizes is NULL), a permutation of the points that groups those assigned to each partition in order (i.e., partition 0 first, partition 1 next, etc.)
8425680f57bSMatthew G. Knepley 
8435680f57bSMatthew G. Knepley   Level: developer
8445680f57bSMatthew G. Knepley 
845b7e49471SLawrence Mitchell   Notes:
846b7e49471SLawrence Mitchell 
847b7e49471SLawrence Mitchell     It is safe to free the sizes and points arrays after use in this routine.
848b7e49471SLawrence Mitchell 
8495680f57bSMatthew G. Knepley .seealso DMPlexDistribute(), PetscPartitionerCreate()
8505680f57bSMatthew G. Knepley @*/
8519852e123SBarry Smith PetscErrorCode PetscPartitionerShellSetPartition(PetscPartitioner part, PetscInt size, const PetscInt sizes[], const PetscInt points[])
8525680f57bSMatthew G. Knepley {
8535680f57bSMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
8545680f57bSMatthew G. Knepley   PetscInt                proc, numPoints;
8555680f57bSMatthew G. Knepley   PetscErrorCode          ierr;
8565680f57bSMatthew G. Knepley 
8575680f57bSMatthew G. Knepley   PetscFunctionBegin;
8585680f57bSMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
8595680f57bSMatthew G. Knepley   if (sizes)  {PetscValidPointer(sizes, 3);}
860c717d290SMatthew G. Knepley   if (points) {PetscValidPointer(points, 4);}
8615680f57bSMatthew G. Knepley   ierr = PetscSectionDestroy(&p->section);CHKERRQ(ierr);
8625680f57bSMatthew G. Knepley   ierr = ISDestroy(&p->partition);CHKERRQ(ierr);
8635680f57bSMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) part), &p->section);CHKERRQ(ierr);
8649852e123SBarry Smith   ierr = PetscSectionSetChart(p->section, 0, size);CHKERRQ(ierr);
8655680f57bSMatthew G. Knepley   if (sizes) {
8669852e123SBarry Smith     for (proc = 0; proc < size; ++proc) {
8675680f57bSMatthew G. Knepley       ierr = PetscSectionSetDof(p->section, proc, sizes[proc]);CHKERRQ(ierr);
8685680f57bSMatthew G. Knepley     }
8695680f57bSMatthew G. Knepley   }
8705680f57bSMatthew G. Knepley   ierr = PetscSectionSetUp(p->section);CHKERRQ(ierr);
8715680f57bSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(p->section, &numPoints);CHKERRQ(ierr);
8725680f57bSMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) part), numPoints, points, PETSC_COPY_VALUES, &p->partition);CHKERRQ(ierr);
8735680f57bSMatthew G. Knepley   PetscFunctionReturn(0);
8745680f57bSMatthew G. Knepley }
8755680f57bSMatthew G. Knepley 
876077101c0SMatthew G. Knepley /*@
877077101c0SMatthew G. Knepley   PetscPartitionerShellSetRandom - Set the flag to use a random partition
878077101c0SMatthew G. Knepley 
879077101c0SMatthew G. Knepley   Collective on Part
880077101c0SMatthew G. Knepley 
881077101c0SMatthew G. Knepley   Input Parameters:
882077101c0SMatthew G. Knepley + part   - The PetscPartitioner
883077101c0SMatthew G. Knepley - random - The flag to use a random partition
884077101c0SMatthew G. Knepley 
885077101c0SMatthew G. Knepley   Level: intermediate
886077101c0SMatthew G. Knepley 
887077101c0SMatthew G. Knepley .seealso PetscPartitionerShellGetRandom(), PetscPartitionerCreate()
888077101c0SMatthew G. Knepley @*/
889077101c0SMatthew G. Knepley PetscErrorCode PetscPartitionerShellSetRandom(PetscPartitioner part, PetscBool random)
890077101c0SMatthew G. Knepley {
891077101c0SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
892077101c0SMatthew G. Knepley 
893077101c0SMatthew G. Knepley   PetscFunctionBegin;
894077101c0SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
895077101c0SMatthew G. Knepley   p->random = random;
896077101c0SMatthew G. Knepley   PetscFunctionReturn(0);
897077101c0SMatthew G. Knepley }
898077101c0SMatthew G. Knepley 
899077101c0SMatthew G. Knepley /*@
900077101c0SMatthew G. Knepley   PetscPartitionerShellGetRandom - get the flag to use a random partition
901077101c0SMatthew G. Knepley 
902077101c0SMatthew G. Knepley   Collective on Part
903077101c0SMatthew G. Knepley 
904077101c0SMatthew G. Knepley   Input Parameter:
905077101c0SMatthew G. Knepley . part   - The PetscPartitioner
906077101c0SMatthew G. Knepley 
907077101c0SMatthew G. Knepley   Output Parameter
908077101c0SMatthew G. Knepley . random - The flag to use a random partition
909077101c0SMatthew G. Knepley 
910077101c0SMatthew G. Knepley   Level: intermediate
911077101c0SMatthew G. Knepley 
912077101c0SMatthew G. Knepley .seealso PetscPartitionerShellSetRandom(), PetscPartitionerCreate()
913077101c0SMatthew G. Knepley @*/
914077101c0SMatthew G. Knepley PetscErrorCode PetscPartitionerShellGetRandom(PetscPartitioner part, PetscBool *random)
915077101c0SMatthew G. Knepley {
916077101c0SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
917077101c0SMatthew G. Knepley 
918077101c0SMatthew G. Knepley   PetscFunctionBegin;
919077101c0SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
920077101c0SMatthew G. Knepley   PetscValidPointer(random, 2);
921077101c0SMatthew G. Knepley   *random = p->random;
922077101c0SMatthew G. Knepley   PetscFunctionReturn(0);
923077101c0SMatthew G. Knepley }
924077101c0SMatthew G. Knepley 
925d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerDestroy_Simple(PetscPartitioner part)
926555a9cf8SMatthew G. Knepley {
927555a9cf8SMatthew G. Knepley   PetscPartitioner_Simple *p = (PetscPartitioner_Simple *) part->data;
928555a9cf8SMatthew G. Knepley   PetscErrorCode          ierr;
929555a9cf8SMatthew G. Knepley 
930555a9cf8SMatthew G. Knepley   PetscFunctionBegin;
931555a9cf8SMatthew G. Knepley   ierr = PetscFree(p);CHKERRQ(ierr);
932555a9cf8SMatthew G. Knepley   PetscFunctionReturn(0);
933555a9cf8SMatthew G. Knepley }
934555a9cf8SMatthew G. Knepley 
935d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Simple_Ascii(PetscPartitioner part, PetscViewer viewer)
936555a9cf8SMatthew G. Knepley {
937555a9cf8SMatthew G. Knepley   PetscViewerFormat format;
938555a9cf8SMatthew G. Knepley   PetscErrorCode    ierr;
939555a9cf8SMatthew G. Knepley 
940555a9cf8SMatthew G. Knepley   PetscFunctionBegin;
941555a9cf8SMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
942555a9cf8SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Simple Graph Partitioner:\n");CHKERRQ(ierr);
943555a9cf8SMatthew G. Knepley   PetscFunctionReturn(0);
944555a9cf8SMatthew G. Knepley }
945555a9cf8SMatthew G. Knepley 
946d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Simple(PetscPartitioner part, PetscViewer viewer)
947555a9cf8SMatthew G. Knepley {
948555a9cf8SMatthew G. Knepley   PetscBool      iascii;
949555a9cf8SMatthew G. Knepley   PetscErrorCode ierr;
950555a9cf8SMatthew G. Knepley 
951555a9cf8SMatthew G. Knepley   PetscFunctionBegin;
952555a9cf8SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
953555a9cf8SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
954555a9cf8SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
955555a9cf8SMatthew G. Knepley   if (iascii) {ierr = PetscPartitionerView_Simple_Ascii(part, viewer);CHKERRQ(ierr);}
956555a9cf8SMatthew G. Knepley   PetscFunctionReturn(0);
957555a9cf8SMatthew G. Knepley }
958555a9cf8SMatthew G. Knepley 
959d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_Simple(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
960555a9cf8SMatthew G. Knepley {
961cead94edSToby Isaac   MPI_Comm       comm;
962555a9cf8SMatthew G. Knepley   PetscInt       np;
963cead94edSToby Isaac   PetscMPIInt    size;
964555a9cf8SMatthew G. Knepley   PetscErrorCode ierr;
965555a9cf8SMatthew G. Knepley 
966555a9cf8SMatthew G. Knepley   PetscFunctionBegin;
967cead94edSToby Isaac   comm = PetscObjectComm((PetscObject)dm);
968cead94edSToby Isaac   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
969555a9cf8SMatthew G. Knepley   ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
970555a9cf8SMatthew G. Knepley   ierr = ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition);CHKERRQ(ierr);
971cead94edSToby Isaac   if (size == 1) {
972cead94edSToby Isaac     for (np = 0; np < nparts; ++np) {ierr = PetscSectionSetDof(partSection, np, numVertices/nparts + ((numVertices % nparts) > np));CHKERRQ(ierr);}
973cead94edSToby Isaac   }
974cead94edSToby Isaac   else {
975cead94edSToby Isaac     PetscMPIInt rank;
976cead94edSToby Isaac     PetscInt nvGlobal, *offsets, myFirst, myLast;
977cead94edSToby Isaac 
978a679a563SToby Isaac     ierr = PetscMalloc1(size+1,&offsets);CHKERRQ(ierr);
979cead94edSToby Isaac     offsets[0] = 0;
980cead94edSToby Isaac     ierr = MPI_Allgather(&numVertices,1,MPIU_INT,&offsets[1],1,MPIU_INT,comm);CHKERRQ(ierr);
981cead94edSToby Isaac     for (np = 2; np <= size; np++) {
982cead94edSToby Isaac       offsets[np] += offsets[np-1];
983cead94edSToby Isaac     }
984cead94edSToby Isaac     nvGlobal = offsets[size];
985cead94edSToby Isaac     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
986cead94edSToby Isaac     myFirst = offsets[rank];
987cead94edSToby Isaac     myLast  = offsets[rank + 1] - 1;
988cead94edSToby Isaac     ierr = PetscFree(offsets);CHKERRQ(ierr);
989cead94edSToby Isaac     if (numVertices) {
990cead94edSToby Isaac       PetscInt firstPart = 0, firstLargePart = 0;
991cead94edSToby Isaac       PetscInt lastPart = 0, lastLargePart = 0;
992cead94edSToby Isaac       PetscInt rem = nvGlobal % nparts;
993cead94edSToby Isaac       PetscInt pSmall = nvGlobal/nparts;
994cead94edSToby Isaac       PetscInt pBig = nvGlobal/nparts + 1;
995cead94edSToby Isaac 
996cead94edSToby Isaac 
997cead94edSToby Isaac       if (rem) {
998cead94edSToby Isaac         firstLargePart = myFirst / pBig;
999cead94edSToby Isaac         lastLargePart  = myLast  / pBig;
1000cead94edSToby Isaac 
1001cead94edSToby Isaac         if (firstLargePart < rem) {
1002cead94edSToby Isaac           firstPart = firstLargePart;
1003cead94edSToby Isaac         }
1004cead94edSToby Isaac         else {
1005cead94edSToby Isaac           firstPart = rem + (myFirst - (rem * pBig)) / pSmall;
1006cead94edSToby Isaac         }
1007cead94edSToby Isaac         if (lastLargePart < rem) {
1008cead94edSToby Isaac           lastPart = lastLargePart;
1009cead94edSToby Isaac         }
1010cead94edSToby Isaac         else {
1011cead94edSToby Isaac           lastPart = rem + (myLast - (rem * pBig)) / pSmall;
1012cead94edSToby Isaac         }
1013cead94edSToby Isaac       }
1014cead94edSToby Isaac       else {
1015cead94edSToby Isaac         firstPart = myFirst / (nvGlobal/nparts);
1016cead94edSToby Isaac         lastPart  = myLast  / (nvGlobal/nparts);
1017cead94edSToby Isaac       }
1018cead94edSToby Isaac 
1019cead94edSToby Isaac       for (np = firstPart; np <= lastPart; np++) {
1020cead94edSToby Isaac         PetscInt PartStart =  np    * (nvGlobal/nparts) + PetscMin(nvGlobal % nparts,np);
1021cead94edSToby Isaac         PetscInt PartEnd   = (np+1) * (nvGlobal/nparts) + PetscMin(nvGlobal % nparts,np+1);
1022cead94edSToby Isaac 
1023cead94edSToby Isaac         PartStart = PetscMax(PartStart,myFirst);
1024cead94edSToby Isaac         PartEnd   = PetscMin(PartEnd,myLast+1);
1025cead94edSToby Isaac         ierr = PetscSectionSetDof(partSection,np,PartEnd-PartStart);CHKERRQ(ierr);
1026cead94edSToby Isaac       }
1027cead94edSToby Isaac     }
1028cead94edSToby Isaac   }
1029cead94edSToby Isaac   ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
1030555a9cf8SMatthew G. Knepley   PetscFunctionReturn(0);
1031555a9cf8SMatthew G. Knepley }
1032555a9cf8SMatthew G. Knepley 
1033d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerInitialize_Simple(PetscPartitioner part)
1034555a9cf8SMatthew G. Knepley {
1035555a9cf8SMatthew G. Knepley   PetscFunctionBegin;
1036555a9cf8SMatthew G. Knepley   part->ops->view      = PetscPartitionerView_Simple;
1037555a9cf8SMatthew G. Knepley   part->ops->destroy   = PetscPartitionerDestroy_Simple;
1038555a9cf8SMatthew G. Knepley   part->ops->partition = PetscPartitionerPartition_Simple;
1039555a9cf8SMatthew G. Knepley   PetscFunctionReturn(0);
1040555a9cf8SMatthew G. Knepley }
1041555a9cf8SMatthew G. Knepley 
1042555a9cf8SMatthew G. Knepley /*MC
1043555a9cf8SMatthew G. Knepley   PETSCPARTITIONERSIMPLE = "simple" - A PetscPartitioner object
1044555a9cf8SMatthew G. Knepley 
1045555a9cf8SMatthew G. Knepley   Level: intermediate
1046555a9cf8SMatthew G. Knepley 
1047555a9cf8SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
1048555a9cf8SMatthew G. Knepley M*/
1049555a9cf8SMatthew G. Knepley 
1050555a9cf8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Simple(PetscPartitioner part)
1051555a9cf8SMatthew G. Knepley {
1052555a9cf8SMatthew G. Knepley   PetscPartitioner_Simple *p;
1053555a9cf8SMatthew G. Knepley   PetscErrorCode           ierr;
1054555a9cf8SMatthew G. Knepley 
1055555a9cf8SMatthew G. Knepley   PetscFunctionBegin;
1056555a9cf8SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
1057555a9cf8SMatthew G. Knepley   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
1058555a9cf8SMatthew G. Knepley   part->data = p;
1059555a9cf8SMatthew G. Knepley 
1060555a9cf8SMatthew G. Knepley   ierr = PetscPartitionerInitialize_Simple(part);CHKERRQ(ierr);
1061555a9cf8SMatthew G. Knepley   PetscFunctionReturn(0);
1062555a9cf8SMatthew G. Knepley }
1063555a9cf8SMatthew G. Knepley 
1064d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerDestroy_Gather(PetscPartitioner part)
1065dae52e14SToby Isaac {
1066dae52e14SToby Isaac   PetscPartitioner_Gather *p = (PetscPartitioner_Gather *) part->data;
1067dae52e14SToby Isaac   PetscErrorCode          ierr;
1068dae52e14SToby Isaac 
1069dae52e14SToby Isaac   PetscFunctionBegin;
1070dae52e14SToby Isaac   ierr = PetscFree(p);CHKERRQ(ierr);
1071dae52e14SToby Isaac   PetscFunctionReturn(0);
1072dae52e14SToby Isaac }
1073dae52e14SToby Isaac 
1074d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Gather_Ascii(PetscPartitioner part, PetscViewer viewer)
1075dae52e14SToby Isaac {
1076dae52e14SToby Isaac   PetscViewerFormat format;
1077dae52e14SToby Isaac   PetscErrorCode    ierr;
1078dae52e14SToby Isaac 
1079dae52e14SToby Isaac   PetscFunctionBegin;
1080dae52e14SToby Isaac   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
1081dae52e14SToby Isaac   ierr = PetscViewerASCIIPrintf(viewer, "Gather Graph Partitioner:\n");CHKERRQ(ierr);
1082dae52e14SToby Isaac   PetscFunctionReturn(0);
1083dae52e14SToby Isaac }
1084dae52e14SToby Isaac 
1085d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Gather(PetscPartitioner part, PetscViewer viewer)
1086dae52e14SToby Isaac {
1087dae52e14SToby Isaac   PetscBool      iascii;
1088dae52e14SToby Isaac   PetscErrorCode ierr;
1089dae52e14SToby Isaac 
1090dae52e14SToby Isaac   PetscFunctionBegin;
1091dae52e14SToby Isaac   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
1092dae52e14SToby Isaac   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
1093dae52e14SToby Isaac   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
1094dae52e14SToby Isaac   if (iascii) {ierr = PetscPartitionerView_Gather_Ascii(part, viewer);CHKERRQ(ierr);}
1095dae52e14SToby Isaac   PetscFunctionReturn(0);
1096dae52e14SToby Isaac }
1097dae52e14SToby Isaac 
1098d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_Gather(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
1099dae52e14SToby Isaac {
1100dae52e14SToby Isaac   PetscInt       np;
1101dae52e14SToby Isaac   PetscErrorCode ierr;
1102dae52e14SToby Isaac 
1103dae52e14SToby Isaac   PetscFunctionBegin;
1104dae52e14SToby Isaac   ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
1105dae52e14SToby Isaac   ierr = ISCreateStride(PETSC_COMM_SELF, numVertices, 0, 1, partition);CHKERRQ(ierr);
1106dae52e14SToby Isaac   ierr = PetscSectionSetDof(partSection,0,numVertices);CHKERRQ(ierr);
1107dae52e14SToby Isaac   for (np = 1; np < nparts; ++np) {ierr = PetscSectionSetDof(partSection, np, 0);CHKERRQ(ierr);}
1108dae52e14SToby Isaac   ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
1109dae52e14SToby Isaac   PetscFunctionReturn(0);
1110dae52e14SToby Isaac }
1111dae52e14SToby Isaac 
1112d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerInitialize_Gather(PetscPartitioner part)
1113dae52e14SToby Isaac {
1114dae52e14SToby Isaac   PetscFunctionBegin;
1115dae52e14SToby Isaac   part->ops->view      = PetscPartitionerView_Gather;
1116dae52e14SToby Isaac   part->ops->destroy   = PetscPartitionerDestroy_Gather;
1117dae52e14SToby Isaac   part->ops->partition = PetscPartitionerPartition_Gather;
1118dae52e14SToby Isaac   PetscFunctionReturn(0);
1119dae52e14SToby Isaac }
1120dae52e14SToby Isaac 
1121dae52e14SToby Isaac /*MC
1122dae52e14SToby Isaac   PETSCPARTITIONERGATHER = "gather" - A PetscPartitioner object
1123dae52e14SToby Isaac 
1124dae52e14SToby Isaac   Level: intermediate
1125dae52e14SToby Isaac 
1126dae52e14SToby Isaac .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
1127dae52e14SToby Isaac M*/
1128dae52e14SToby Isaac 
1129dae52e14SToby Isaac PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Gather(PetscPartitioner part)
1130dae52e14SToby Isaac {
1131dae52e14SToby Isaac   PetscPartitioner_Gather *p;
1132dae52e14SToby Isaac   PetscErrorCode           ierr;
1133dae52e14SToby Isaac 
1134dae52e14SToby Isaac   PetscFunctionBegin;
1135dae52e14SToby Isaac   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
1136dae52e14SToby Isaac   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
1137dae52e14SToby Isaac   part->data = p;
1138dae52e14SToby Isaac 
1139dae52e14SToby Isaac   ierr = PetscPartitionerInitialize_Gather(part);CHKERRQ(ierr);
1140dae52e14SToby Isaac   PetscFunctionReturn(0);
1141dae52e14SToby Isaac }
1142dae52e14SToby Isaac 
1143dae52e14SToby Isaac 
1144d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerDestroy_Chaco(PetscPartitioner part)
114577623264SMatthew G. Knepley {
114677623264SMatthew G. Knepley   PetscPartitioner_Chaco *p = (PetscPartitioner_Chaco *) part->data;
114777623264SMatthew G. Knepley   PetscErrorCode          ierr;
114877623264SMatthew G. Knepley 
114977623264SMatthew G. Knepley   PetscFunctionBegin;
115077623264SMatthew G. Knepley   ierr = PetscFree(p);CHKERRQ(ierr);
115177623264SMatthew G. Knepley   PetscFunctionReturn(0);
115277623264SMatthew G. Knepley }
115377623264SMatthew G. Knepley 
1154d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Chaco_Ascii(PetscPartitioner part, PetscViewer viewer)
115577623264SMatthew G. Knepley {
115677623264SMatthew G. Knepley   PetscViewerFormat format;
115777623264SMatthew G. Knepley   PetscErrorCode    ierr;
115877623264SMatthew G. Knepley 
115977623264SMatthew G. Knepley   PetscFunctionBegin;
116077623264SMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
116177623264SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Chaco Graph Partitioner:\n");CHKERRQ(ierr);
116277623264SMatthew G. Knepley   PetscFunctionReturn(0);
116377623264SMatthew G. Knepley }
116477623264SMatthew G. Knepley 
1165d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_Chaco(PetscPartitioner part, PetscViewer viewer)
116677623264SMatthew G. Knepley {
116777623264SMatthew G. Knepley   PetscBool      iascii;
116877623264SMatthew G. Knepley   PetscErrorCode ierr;
116977623264SMatthew G. Knepley 
117077623264SMatthew G. Knepley   PetscFunctionBegin;
117177623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
117277623264SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
117377623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
117477623264SMatthew G. Knepley   if (iascii) {ierr = PetscPartitionerView_Chaco_Ascii(part, viewer);CHKERRQ(ierr);}
117577623264SMatthew G. Knepley   PetscFunctionReturn(0);
117677623264SMatthew G. Knepley }
117777623264SMatthew G. Knepley 
117870034214SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO)
117970034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
118070034214SMatthew G. Knepley #include <unistd.h>
118170034214SMatthew G. Knepley #endif
118211d1e910SBarry Smith #if defined(PETSC_HAVE_CHACO_INT_ASSIGNMENT)
118311d1e910SBarry Smith #include <chaco.h>
118411d1e910SBarry Smith #else
118511d1e910SBarry Smith /* Older versions of Chaco do not have an include file */
118670034214SMatthew G. Knepley PETSC_EXTERN int interface(int nvtxs, int *start, int *adjacency, int *vwgts,
118770034214SMatthew G. Knepley                        float *ewgts, float *x, float *y, float *z, char *outassignname,
118870034214SMatthew G. Knepley                        char *outfilename, short *assignment, int architecture, int ndims_tot,
118970034214SMatthew G. Knepley                        int mesh_dims[3], double *goal, int global_method, int local_method,
119070034214SMatthew G. Knepley                        int rqi_flag, int vmax, int ndims, double eigtol, long seed);
119111d1e910SBarry Smith #endif
119270034214SMatthew G. Knepley extern int FREE_GRAPH;
119377623264SMatthew G. Knepley #endif
119470034214SMatthew G. Knepley 
1195d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_Chaco(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
119670034214SMatthew G. Knepley {
119777623264SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO)
119870034214SMatthew G. Knepley   enum {DEFAULT_METHOD = 1, INERTIAL_METHOD = 3};
119970034214SMatthew G. Knepley   MPI_Comm       comm;
120070034214SMatthew G. Knepley   int            nvtxs          = numVertices; /* number of vertices in full graph */
120170034214SMatthew G. Knepley   int           *vwgts          = NULL;   /* weights for all vertices */
120270034214SMatthew G. Knepley   float         *ewgts          = NULL;   /* weights for all edges */
120370034214SMatthew G. Knepley   float         *x              = NULL, *y = NULL, *z = NULL; /* coordinates for inertial method */
120470034214SMatthew G. Knepley   char          *outassignname  = NULL;   /*  name of assignment output file */
120570034214SMatthew G. Knepley   char          *outfilename    = NULL;   /* output file name */
120670034214SMatthew G. Knepley   int            architecture   = 1;      /* 0 => hypercube, d => d-dimensional mesh */
120770034214SMatthew G. Knepley   int            ndims_tot      = 0;      /* total number of cube dimensions to divide */
120870034214SMatthew G. Knepley   int            mesh_dims[3];            /* dimensions of mesh of processors */
120970034214SMatthew G. Knepley   double        *goal          = NULL;    /* desired set sizes for each set */
121070034214SMatthew G. Knepley   int            global_method = 1;       /* global partitioning algorithm */
121170034214SMatthew G. Knepley   int            local_method  = 1;       /* local partitioning algorithm */
121270034214SMatthew G. Knepley   int            rqi_flag      = 0;       /* should I use RQI/Symmlq eigensolver? */
121370034214SMatthew G. Knepley   int            vmax          = 200;     /* how many vertices to coarsen down to? */
121470034214SMatthew G. Knepley   int            ndims         = 1;       /* number of eigenvectors (2^d sets) */
121570034214SMatthew G. Knepley   double         eigtol        = 0.001;   /* tolerance on eigenvectors */
121670034214SMatthew G. Knepley   long           seed          = 123636512; /* for random graph mutations */
121711d1e910SBarry Smith #if defined(PETSC_HAVE_CHACO_INT_ASSIGNMENT)
121811d1e910SBarry Smith   int           *assignment;              /* Output partition */
121911d1e910SBarry Smith #else
122070034214SMatthew G. Knepley   short int     *assignment;              /* Output partition */
122111d1e910SBarry Smith #endif
122270034214SMatthew G. Knepley   int            fd_stdout, fd_pipe[2];
122370034214SMatthew G. Knepley   PetscInt      *points;
122470034214SMatthew G. Knepley   int            i, v, p;
122570034214SMatthew G. Knepley   PetscErrorCode ierr;
122670034214SMatthew G. Knepley 
122770034214SMatthew G. Knepley   PetscFunctionBegin;
122870034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
122970034214SMatthew G. Knepley   if (!numVertices) {
123077623264SMatthew G. Knepley     ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
123177623264SMatthew G. Knepley     ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
123270034214SMatthew G. Knepley     ierr = ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
123370034214SMatthew G. Knepley     PetscFunctionReturn(0);
123470034214SMatthew G. Knepley   }
123570034214SMatthew G. Knepley   FREE_GRAPH = 0;                         /* Do not let Chaco free my memory */
123670034214SMatthew G. Knepley   for (i = 0; i < start[numVertices]; ++i) ++adjacency[i];
123770034214SMatthew G. Knepley 
123870034214SMatthew G. Knepley   if (global_method == INERTIAL_METHOD) {
123970034214SMatthew G. Knepley     /* manager.createCellCoordinates(nvtxs, &x, &y, &z); */
124070034214SMatthew G. Knepley     SETERRQ(comm, PETSC_ERR_SUP, "Inertial partitioning not yet supported");
124170034214SMatthew G. Knepley   }
124277623264SMatthew G. Knepley   mesh_dims[0] = nparts;
124370034214SMatthew G. Knepley   mesh_dims[1] = 1;
124470034214SMatthew G. Knepley   mesh_dims[2] = 1;
124570034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &assignment);CHKERRQ(ierr);
124670034214SMatthew G. Knepley   /* Chaco outputs to stdout. We redirect this to a buffer. */
124770034214SMatthew G. Knepley   /* TODO: check error codes for UNIX calls */
124870034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
124970034214SMatthew G. Knepley   {
125070034214SMatthew G. Knepley     int piperet;
125170034214SMatthew G. Knepley     piperet = pipe(fd_pipe);
125270034214SMatthew G. Knepley     if (piperet) SETERRQ(comm,PETSC_ERR_SYS,"Could not create pipe");
125370034214SMatthew G. Knepley     fd_stdout = dup(1);
125470034214SMatthew G. Knepley     close(1);
125570034214SMatthew G. Knepley     dup2(fd_pipe[1], 1);
125670034214SMatthew G. Knepley   }
125770034214SMatthew G. Knepley #endif
125870034214SMatthew G. Knepley   ierr = interface(nvtxs, (int*) start, (int*) adjacency, vwgts, ewgts, x, y, z, outassignname, outfilename,
125970034214SMatthew G. Knepley                    assignment, architecture, ndims_tot, mesh_dims, goal, global_method, local_method, rqi_flag,
126070034214SMatthew G. Knepley                    vmax, ndims, eigtol, seed);
126170034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
126270034214SMatthew G. Knepley   {
126370034214SMatthew G. Knepley     char msgLog[10000];
126470034214SMatthew G. Knepley     int  count;
126570034214SMatthew G. Knepley 
126670034214SMatthew G. Knepley     fflush(stdout);
126770034214SMatthew G. Knepley     count = read(fd_pipe[0], msgLog, (10000-1)*sizeof(char));
126870034214SMatthew G. Knepley     if (count < 0) count = 0;
126970034214SMatthew G. Knepley     msgLog[count] = 0;
127070034214SMatthew G. Knepley     close(1);
127170034214SMatthew G. Knepley     dup2(fd_stdout, 1);
127270034214SMatthew G. Knepley     close(fd_stdout);
127370034214SMatthew G. Knepley     close(fd_pipe[0]);
127470034214SMatthew G. Knepley     close(fd_pipe[1]);
127570034214SMatthew G. Knepley     if (ierr) SETERRQ1(comm, PETSC_ERR_LIB, "Error in Chaco library: %s", msgLog);
127670034214SMatthew G. Knepley   }
127770034214SMatthew G. Knepley #endif
127870034214SMatthew G. Knepley   /* Convert to PetscSection+IS */
127977623264SMatthew G. Knepley   ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
128070034214SMatthew G. Knepley   for (v = 0; v < nvtxs; ++v) {
128177623264SMatthew G. Knepley     ierr = PetscSectionAddDof(partSection, assignment[v], 1);CHKERRQ(ierr);
128270034214SMatthew G. Knepley   }
128377623264SMatthew G. Knepley   ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
128470034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr);
128577623264SMatthew G. Knepley   for (p = 0, i = 0; p < nparts; ++p) {
128670034214SMatthew G. Knepley     for (v = 0; v < nvtxs; ++v) {
128770034214SMatthew G. Knepley       if (assignment[v] == p) points[i++] = v;
128870034214SMatthew G. Knepley     }
128970034214SMatthew G. Knepley   }
129070034214SMatthew G. Knepley   if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs);
129170034214SMatthew G. Knepley   ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
129270034214SMatthew G. Knepley   if (global_method == INERTIAL_METHOD) {
129370034214SMatthew G. Knepley     /* manager.destroyCellCoordinates(nvtxs, &x, &y, &z); */
129470034214SMatthew G. Knepley   }
129570034214SMatthew G. Knepley   ierr = PetscFree(assignment);CHKERRQ(ierr);
129670034214SMatthew G. Knepley   for (i = 0; i < start[numVertices]; ++i) --adjacency[i];
129770034214SMatthew G. Knepley   PetscFunctionReturn(0);
129877623264SMatthew G. Knepley #else
129977623264SMatthew G. Knepley   SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-chaco.");
130070034214SMatthew G. Knepley #endif
130177623264SMatthew G. Knepley }
130277623264SMatthew G. Knepley 
1303d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerInitialize_Chaco(PetscPartitioner part)
130477623264SMatthew G. Knepley {
130577623264SMatthew G. Knepley   PetscFunctionBegin;
130677623264SMatthew G. Knepley   part->ops->view      = PetscPartitionerView_Chaco;
130777623264SMatthew G. Knepley   part->ops->destroy   = PetscPartitionerDestroy_Chaco;
130877623264SMatthew G. Knepley   part->ops->partition = PetscPartitionerPartition_Chaco;
130977623264SMatthew G. Knepley   PetscFunctionReturn(0);
131077623264SMatthew G. Knepley }
131177623264SMatthew G. Knepley 
131277623264SMatthew G. Knepley /*MC
131377623264SMatthew G. Knepley   PETSCPARTITIONERCHACO = "chaco" - A PetscPartitioner object using the Chaco library
131477623264SMatthew G. Knepley 
131577623264SMatthew G. Knepley   Level: intermediate
131677623264SMatthew G. Knepley 
131777623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
131877623264SMatthew G. Knepley M*/
131977623264SMatthew G. Knepley 
132077623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Chaco(PetscPartitioner part)
132177623264SMatthew G. Knepley {
132277623264SMatthew G. Knepley   PetscPartitioner_Chaco *p;
132377623264SMatthew G. Knepley   PetscErrorCode          ierr;
132477623264SMatthew G. Knepley 
132577623264SMatthew G. Knepley   PetscFunctionBegin;
132677623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
132777623264SMatthew G. Knepley   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
132877623264SMatthew G. Knepley   part->data = p;
132977623264SMatthew G. Knepley 
133077623264SMatthew G. Knepley   ierr = PetscPartitionerInitialize_Chaco(part);CHKERRQ(ierr);
133177623264SMatthew G. Knepley   ierr = PetscCitationsRegister(ChacoPartitionerCitation, &ChacoPartitionercite);CHKERRQ(ierr);
133277623264SMatthew G. Knepley   PetscFunctionReturn(0);
133377623264SMatthew G. Knepley }
133477623264SMatthew G. Knepley 
1335d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerDestroy_ParMetis(PetscPartitioner part)
133677623264SMatthew G. Knepley {
133777623264SMatthew G. Knepley   PetscPartitioner_ParMetis *p = (PetscPartitioner_ParMetis *) part->data;
133877623264SMatthew G. Knepley   PetscErrorCode             ierr;
133977623264SMatthew G. Knepley 
134077623264SMatthew G. Knepley   PetscFunctionBegin;
134177623264SMatthew G. Knepley   ierr = PetscFree(p);CHKERRQ(ierr);
134277623264SMatthew G. Knepley   PetscFunctionReturn(0);
134377623264SMatthew G. Knepley }
134477623264SMatthew G. Knepley 
1345d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_ParMetis_Ascii(PetscPartitioner part, PetscViewer viewer)
134677623264SMatthew G. Knepley {
134777623264SMatthew G. Knepley   PetscViewerFormat format;
134877623264SMatthew G. Knepley   PetscErrorCode    ierr;
134977623264SMatthew G. Knepley 
135077623264SMatthew G. Knepley   PetscFunctionBegin;
135177623264SMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
135277623264SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "ParMetis Graph Partitioner:\n");CHKERRQ(ierr);
135377623264SMatthew G. Knepley   PetscFunctionReturn(0);
135477623264SMatthew G. Knepley }
135577623264SMatthew G. Knepley 
1356d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerView_ParMetis(PetscPartitioner part, PetscViewer viewer)
135777623264SMatthew G. Knepley {
135877623264SMatthew G. Knepley   PetscBool      iascii;
135977623264SMatthew G. Knepley   PetscErrorCode ierr;
136077623264SMatthew G. Knepley 
136177623264SMatthew G. Knepley   PetscFunctionBegin;
136277623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
136377623264SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
136477623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
136577623264SMatthew G. Knepley   if (iascii) {ierr = PetscPartitionerView_ParMetis_Ascii(part, viewer);CHKERRQ(ierr);}
136677623264SMatthew G. Knepley   PetscFunctionReturn(0);
136777623264SMatthew G. Knepley }
136870034214SMatthew G. Knepley 
136970034214SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS)
137070034214SMatthew G. Knepley #include <parmetis.h>
137177623264SMatthew G. Knepley #endif
137270034214SMatthew G. Knepley 
1373d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerPartition_ParMetis(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
137470034214SMatthew G. Knepley {
137577623264SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS)
137670034214SMatthew G. Knepley   MPI_Comm       comm;
1377deea36a5SMatthew G. Knepley   PetscSection   section;
137870034214SMatthew G. Knepley   PetscInt       nvtxs      = numVertices; /* The number of vertices in full graph */
137970034214SMatthew G. Knepley   PetscInt      *vtxdist;                  /* Distribution of vertices across processes */
138070034214SMatthew G. Knepley   PetscInt      *xadj       = start;       /* Start of edge list for each vertex */
138170034214SMatthew G. Knepley   PetscInt      *adjncy     = adjacency;   /* Edge lists for all vertices */
138270034214SMatthew G. Knepley   PetscInt      *vwgt       = NULL;        /* Vertex weights */
138370034214SMatthew G. Knepley   PetscInt      *adjwgt     = NULL;        /* Edge weights */
138470034214SMatthew G. Knepley   PetscInt       wgtflag    = 0;           /* Indicates which weights are present */
138570034214SMatthew G. Knepley   PetscInt       numflag    = 0;           /* Indicates initial offset (0 or 1) */
138670034214SMatthew G. Knepley   PetscInt       ncon       = 1;           /* The number of weights per vertex */
138770034214SMatthew G. Knepley   PetscReal     *tpwgts;                   /* The fraction of vertex weights assigned to each partition */
138870034214SMatthew G. Knepley   PetscReal     *ubvec;                    /* The balance intolerance for vertex weights */
138970034214SMatthew G. Knepley   PetscInt       options[5];               /* Options */
139070034214SMatthew G. Knepley   /* Outputs */
139170034214SMatthew G. Knepley   PetscInt       edgeCut;                  /* The number of edges cut by the partition */
139270034214SMatthew G. Knepley   PetscInt      *assignment, *points;
139377623264SMatthew G. Knepley   PetscMPIInt    rank, p, v, i;
139470034214SMatthew G. Knepley   PetscErrorCode ierr;
139570034214SMatthew G. Knepley 
139670034214SMatthew G. Knepley   PetscFunctionBegin;
139777623264SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) part, &comm);CHKERRQ(ierr);
139870034214SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
139970034214SMatthew G. Knepley   options[0] = 0; /* Use all defaults */
140070034214SMatthew G. Knepley   /* Calculate vertex distribution */
1401cd0de0f2SShri   ierr = PetscMalloc5(nparts+1,&vtxdist,nparts*ncon,&tpwgts,ncon,&ubvec,nvtxs,&assignment,nvtxs,&vwgt);CHKERRQ(ierr);
140270034214SMatthew G. Knepley   vtxdist[0] = 0;
140370034214SMatthew G. Knepley   ierr = MPI_Allgather(&nvtxs, 1, MPIU_INT, &vtxdist[1], 1, MPIU_INT, comm);CHKERRQ(ierr);
140470034214SMatthew G. Knepley   for (p = 2; p <= nparts; ++p) {
140570034214SMatthew G. Knepley     vtxdist[p] += vtxdist[p-1];
140670034214SMatthew G. Knepley   }
140770034214SMatthew G. Knepley   /* Calculate weights */
140870034214SMatthew G. Knepley   for (p = 0; p < nparts; ++p) {
140970034214SMatthew G. Knepley     tpwgts[p] = 1.0/nparts;
141070034214SMatthew G. Knepley   }
141170034214SMatthew G. Knepley   ubvec[0] = 1.05;
1412deea36a5SMatthew G. Knepley   /* Weight cells by dofs on cell by default */
1413deea36a5SMatthew G. Knepley   ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1414deea36a5SMatthew G. Knepley   if (section) {
1415deea36a5SMatthew G. Knepley     PetscInt cStart, cEnd, dof;
141670034214SMatthew G. Knepley 
1417deea36a5SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
1418deea36a5SMatthew G. Knepley     for (v = cStart; v < cEnd; ++v) {
1419deea36a5SMatthew G. Knepley       ierr = PetscSectionGetDof(section, v, &dof);CHKERRQ(ierr);
1420deea36a5SMatthew G. Knepley       vwgt[v-cStart] = PetscMax(dof, 1);
1421deea36a5SMatthew G. Knepley     }
1422deea36a5SMatthew G. Knepley   } else {
1423deea36a5SMatthew G. Knepley     for (v = 0; v < nvtxs; ++v) vwgt[v] = 1;
1424cd0de0f2SShri   }
1425cd0de0f2SShri 
142670034214SMatthew G. Knepley   if (nparts == 1) {
14279fc93327SToby Isaac     ierr = PetscMemzero(assignment, nvtxs * sizeof(PetscInt));CHKERRQ(ierr);
142870034214SMatthew G. Knepley   } else {
142970034214SMatthew G. Knepley     if (vtxdist[1] == vtxdist[nparts]) {
143070034214SMatthew G. Knepley       if (!rank) {
143170034214SMatthew G. Knepley         PetscStackPush("METIS_PartGraphKway");
143270034214SMatthew G. Knepley         ierr = METIS_PartGraphKway(&nvtxs, &ncon, xadj, adjncy, vwgt, NULL, adjwgt, &nparts, tpwgts, ubvec, NULL, &edgeCut, assignment);
143370034214SMatthew G. Knepley         PetscStackPop;
143470034214SMatthew G. Knepley         if (ierr != METIS_OK) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in METIS_PartGraphKway()");
143570034214SMatthew G. Knepley       }
143670034214SMatthew G. Knepley     } else {
143770034214SMatthew G. Knepley       PetscStackPush("ParMETIS_V3_PartKway");
143870034214SMatthew G. Knepley       ierr = ParMETIS_V3_PartKway(vtxdist, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgeCut, assignment, &comm);
143970034214SMatthew G. Knepley       PetscStackPop;
1440c717d290SMatthew G. Knepley       if (ierr != METIS_OK) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error %d in ParMETIS_V3_PartKway()", ierr);
144170034214SMatthew G. Knepley     }
144270034214SMatthew G. Knepley   }
144370034214SMatthew G. Knepley   /* Convert to PetscSection+IS */
144477623264SMatthew G. Knepley   ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
144577623264SMatthew G. Knepley   for (v = 0; v < nvtxs; ++v) {ierr = PetscSectionAddDof(partSection, assignment[v], 1);CHKERRQ(ierr);}
144677623264SMatthew G. Knepley   ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
144770034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr);
144877623264SMatthew G. Knepley   for (p = 0, i = 0; p < nparts; ++p) {
144970034214SMatthew G. Knepley     for (v = 0; v < nvtxs; ++v) {
145070034214SMatthew G. Knepley       if (assignment[v] == p) points[i++] = v;
145170034214SMatthew G. Knepley     }
145270034214SMatthew G. Knepley   }
145370034214SMatthew G. Knepley   if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs);
145470034214SMatthew G. Knepley   ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
1455cd0de0f2SShri   ierr = PetscFree5(vtxdist,tpwgts,ubvec,assignment,vwgt);CHKERRQ(ierr);
14569b80ac48SMatthew G. Knepley   PetscFunctionReturn(0);
145770034214SMatthew G. Knepley #else
145877623264SMatthew G. Knepley   SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-parmetis.");
145970034214SMatthew G. Knepley #endif
146070034214SMatthew G. Knepley }
146170034214SMatthew G. Knepley 
1462d5577e40SMatthew G. Knepley static PetscErrorCode PetscPartitionerInitialize_ParMetis(PetscPartitioner part)
146377623264SMatthew G. Knepley {
146477623264SMatthew G. Knepley   PetscFunctionBegin;
146577623264SMatthew G. Knepley   part->ops->view      = PetscPartitionerView_ParMetis;
146677623264SMatthew G. Knepley   part->ops->destroy   = PetscPartitionerDestroy_ParMetis;
146777623264SMatthew G. Knepley   part->ops->partition = PetscPartitionerPartition_ParMetis;
146877623264SMatthew G. Knepley   PetscFunctionReturn(0);
146977623264SMatthew G. Knepley }
147077623264SMatthew G. Knepley 
147177623264SMatthew G. Knepley /*MC
147277623264SMatthew G. Knepley   PETSCPARTITIONERPARMETIS = "parmetis" - A PetscPartitioner object using the ParMetis library
147377623264SMatthew G. Knepley 
147477623264SMatthew G. Knepley   Level: intermediate
147577623264SMatthew G. Knepley 
147677623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
147777623264SMatthew G. Knepley M*/
147877623264SMatthew G. Knepley 
147977623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_ParMetis(PetscPartitioner part)
148077623264SMatthew G. Knepley {
148177623264SMatthew G. Knepley   PetscPartitioner_ParMetis *p;
148277623264SMatthew G. Knepley   PetscErrorCode          ierr;
148377623264SMatthew G. Knepley 
148477623264SMatthew G. Knepley   PetscFunctionBegin;
148577623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
148677623264SMatthew G. Knepley   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
148777623264SMatthew G. Knepley   part->data = p;
148877623264SMatthew G. Knepley 
148977623264SMatthew G. Knepley   ierr = PetscPartitionerInitialize_ParMetis(part);CHKERRQ(ierr);
149077623264SMatthew G. Knepley   ierr = PetscCitationsRegister(ParMetisPartitionerCitation, &ParMetisPartitionercite);CHKERRQ(ierr);
149170034214SMatthew G. Knepley   PetscFunctionReturn(0);
149270034214SMatthew G. Knepley }
149370034214SMatthew G. Knepley 
14945680f57bSMatthew G. Knepley /*@
14955680f57bSMatthew G. Knepley   DMPlexGetPartitioner - Get the mesh partitioner
14965680f57bSMatthew G. Knepley 
14975680f57bSMatthew G. Knepley   Not collective
14985680f57bSMatthew G. Knepley 
14995680f57bSMatthew G. Knepley   Input Parameter:
15005680f57bSMatthew G. Knepley . dm - The DM
15015680f57bSMatthew G. Knepley 
15025680f57bSMatthew G. Knepley   Output Parameter:
15035680f57bSMatthew G. Knepley . part - The PetscPartitioner
15045680f57bSMatthew G. Knepley 
15055680f57bSMatthew G. Knepley   Level: developer
15065680f57bSMatthew G. Knepley 
150798599a47SLawrence Mitchell   Note: This gets a borrowed reference, so the user should not destroy this PetscPartitioner.
150898599a47SLawrence Mitchell 
150998599a47SLawrence Mitchell .seealso DMPlexDistribute(), DMPlexSetPartitioner(), PetscPartitionerCreate()
15105680f57bSMatthew G. Knepley @*/
15115680f57bSMatthew G. Knepley PetscErrorCode DMPlexGetPartitioner(DM dm, PetscPartitioner *part)
15125680f57bSMatthew G. Knepley {
15135680f57bSMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
15144f3833eaSMatthew G. Knepley   PetscErrorCode ierr;
15155680f57bSMatthew G. Knepley 
15165680f57bSMatthew G. Knepley   PetscFunctionBegin;
15175680f57bSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
15185680f57bSMatthew G. Knepley   PetscValidPointer(part, 2);
15194f3833eaSMatthew G. Knepley   ierr = PetscPartitionerSetTypeFromOptions_Internal(mesh->partitioner);CHKERRQ(ierr);
15205680f57bSMatthew G. Knepley   *part = mesh->partitioner;
15215680f57bSMatthew G. Knepley   PetscFunctionReturn(0);
15225680f57bSMatthew G. Knepley }
15235680f57bSMatthew G. Knepley 
152471bb2955SLawrence Mitchell /*@
152571bb2955SLawrence Mitchell   DMPlexSetPartitioner - Set the mesh partitioner
152671bb2955SLawrence Mitchell 
152771bb2955SLawrence Mitchell   logically collective on dm and part
152871bb2955SLawrence Mitchell 
152971bb2955SLawrence Mitchell   Input Parameters:
153071bb2955SLawrence Mitchell + dm - The DM
153171bb2955SLawrence Mitchell - part - The partitioner
153271bb2955SLawrence Mitchell 
153371bb2955SLawrence Mitchell   Level: developer
153471bb2955SLawrence Mitchell 
153571bb2955SLawrence Mitchell   Note: Any existing PetscPartitioner will be destroyed.
153671bb2955SLawrence Mitchell 
153771bb2955SLawrence Mitchell .seealso DMPlexDistribute(), DMPlexGetPartitioner(), PetscPartitionerCreate()
153871bb2955SLawrence Mitchell @*/
153971bb2955SLawrence Mitchell PetscErrorCode DMPlexSetPartitioner(DM dm, PetscPartitioner part)
154071bb2955SLawrence Mitchell {
154171bb2955SLawrence Mitchell   DM_Plex       *mesh = (DM_Plex *) dm->data;
154271bb2955SLawrence Mitchell   PetscErrorCode ierr;
154371bb2955SLawrence Mitchell 
154471bb2955SLawrence Mitchell   PetscFunctionBegin;
154571bb2955SLawrence Mitchell   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
154671bb2955SLawrence Mitchell   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 2);
154771bb2955SLawrence Mitchell   ierr = PetscObjectReference((PetscObject)part);CHKERRQ(ierr);
154871bb2955SLawrence Mitchell   ierr = PetscPartitionerDestroy(&mesh->partitioner);CHKERRQ(ierr);
154971bb2955SLawrence Mitchell   mesh->partitioner = part;
155071bb2955SLawrence Mitchell   PetscFunctionReturn(0);
155171bb2955SLawrence Mitchell }
155271bb2955SLawrence Mitchell 
15536a5a2ffdSToby Isaac static PetscErrorCode DMPlexPartitionLabelClosure_Tree(DM dm, DMLabel label, PetscInt rank, PetscInt point, PetscBool up, PetscBool down)
1554270bba0cSToby Isaac {
1555270bba0cSToby Isaac   PetscErrorCode ierr;
1556270bba0cSToby Isaac 
1557270bba0cSToby Isaac   PetscFunctionBegin;
15586a5a2ffdSToby Isaac   if (up) {
15596a5a2ffdSToby Isaac     PetscInt parent;
15606a5a2ffdSToby Isaac 
1561270bba0cSToby Isaac     ierr = DMPlexGetTreeParent(dm,point,&parent,NULL);CHKERRQ(ierr);
15626a5a2ffdSToby Isaac     if (parent != point) {
15636a5a2ffdSToby Isaac       PetscInt closureSize, *closure = NULL, i;
15646a5a2ffdSToby Isaac 
1565270bba0cSToby Isaac       ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1566270bba0cSToby Isaac       for (i = 0; i < closureSize; i++) {
1567270bba0cSToby Isaac         PetscInt cpoint = closure[2*i];
1568270bba0cSToby Isaac 
1569270bba0cSToby Isaac         ierr = DMLabelSetValue(label,cpoint,rank);CHKERRQ(ierr);
15706a5a2ffdSToby Isaac         ierr = DMPlexPartitionLabelClosure_Tree(dm,label,rank,cpoint,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
1571270bba0cSToby Isaac       }
1572270bba0cSToby Isaac       ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
15736a5a2ffdSToby Isaac     }
15746a5a2ffdSToby Isaac   }
15756a5a2ffdSToby Isaac   if (down) {
15766a5a2ffdSToby Isaac     PetscInt numChildren;
15776a5a2ffdSToby Isaac     const PetscInt *children;
15786a5a2ffdSToby Isaac 
15796a5a2ffdSToby Isaac     ierr = DMPlexGetTreeChildren(dm,point,&numChildren,&children);CHKERRQ(ierr);
15806a5a2ffdSToby Isaac     if (numChildren) {
15816a5a2ffdSToby Isaac       PetscInt i;
15826a5a2ffdSToby Isaac 
15836a5a2ffdSToby Isaac       for (i = 0; i < numChildren; i++) {
15846a5a2ffdSToby Isaac         PetscInt cpoint = children[i];
15856a5a2ffdSToby Isaac 
15866a5a2ffdSToby Isaac         ierr = DMLabelSetValue(label,cpoint,rank);CHKERRQ(ierr);
15876a5a2ffdSToby Isaac         ierr = DMPlexPartitionLabelClosure_Tree(dm,label,rank,cpoint,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr);
15886a5a2ffdSToby Isaac       }
15896a5a2ffdSToby Isaac     }
15906a5a2ffdSToby Isaac   }
1591270bba0cSToby Isaac   PetscFunctionReturn(0);
1592270bba0cSToby Isaac }
1593270bba0cSToby Isaac 
15945abbe4feSMichael Lange /*@
15955abbe4feSMichael Lange   DMPlexPartitionLabelClosure - Add the closure of all points to the partition label
15965abbe4feSMichael Lange 
15975abbe4feSMichael Lange   Input Parameters:
15985abbe4feSMichael Lange + dm     - The DM
15995abbe4feSMichael Lange - label  - DMLabel assinging ranks to remote roots
16005abbe4feSMichael Lange 
16015abbe4feSMichael Lange   Level: developer
16025abbe4feSMichael Lange 
16035abbe4feSMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap
16045abbe4feSMichael Lange @*/
16055abbe4feSMichael Lange PetscErrorCode DMPlexPartitionLabelClosure(DM dm, DMLabel label)
16069ffc88e4SToby Isaac {
16075abbe4feSMichael Lange   IS              rankIS,   pointIS;
16085abbe4feSMichael Lange   const PetscInt *ranks,   *points;
16095abbe4feSMichael Lange   PetscInt        numRanks, numPoints, r, p, c, closureSize;
16105abbe4feSMichael Lange   PetscInt       *closure = NULL;
16119ffc88e4SToby Isaac   PetscErrorCode  ierr;
16129ffc88e4SToby Isaac 
16139ffc88e4SToby Isaac   PetscFunctionBegin;
16145abbe4feSMichael Lange   ierr = DMLabelGetValueIS(label, &rankIS);CHKERRQ(ierr);
16155abbe4feSMichael Lange   ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr);
16165abbe4feSMichael Lange   ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr);
16175abbe4feSMichael Lange   for (r = 0; r < numRanks; ++r) {
16185abbe4feSMichael Lange     const PetscInt rank = ranks[r];
16199ffc88e4SToby Isaac 
16205abbe4feSMichael Lange     ierr = DMLabelGetStratumIS(label, rank, &pointIS);CHKERRQ(ierr);
16215abbe4feSMichael Lange     ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
16225abbe4feSMichael Lange     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
16235abbe4feSMichael Lange     for (p = 0; p < numPoints; ++p) {
16245abbe4feSMichael Lange       ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1625270bba0cSToby Isaac       for (c = 0; c < closureSize*2; c += 2) {
1626270bba0cSToby Isaac         ierr = DMLabelSetValue(label, closure[c], rank);CHKERRQ(ierr);
16276a5a2ffdSToby Isaac         ierr = DMPlexPartitionLabelClosure_Tree(dm,label,rank,closure[c],PETSC_TRUE,PETSC_TRUE);CHKERRQ(ierr);
1628270bba0cSToby Isaac       }
16299ffc88e4SToby Isaac     }
16305abbe4feSMichael Lange     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
16315abbe4feSMichael Lange     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
16329ffc88e4SToby Isaac   }
16337de78196SMichael Lange   if (closure) {ierr = DMPlexRestoreTransitiveClosure(dm, 0, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);}
16345abbe4feSMichael Lange   ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr);
16355abbe4feSMichael Lange   ierr = ISDestroy(&rankIS);CHKERRQ(ierr);
16369ffc88e4SToby Isaac   PetscFunctionReturn(0);
16379ffc88e4SToby Isaac }
16389ffc88e4SToby Isaac 
163924d039d7SMichael Lange /*@
164024d039d7SMichael Lange   DMPlexPartitionLabelAdjacency - Add one level of adjacent points to the partition label
164124d039d7SMichael Lange 
164224d039d7SMichael Lange   Input Parameters:
164324d039d7SMichael Lange + dm     - The DM
164424d039d7SMichael Lange - label  - DMLabel assinging ranks to remote roots
164524d039d7SMichael Lange 
164624d039d7SMichael Lange   Level: developer
164724d039d7SMichael Lange 
164824d039d7SMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap
164924d039d7SMichael Lange @*/
165024d039d7SMichael Lange PetscErrorCode DMPlexPartitionLabelAdjacency(DM dm, DMLabel label)
165170034214SMatthew G. Knepley {
165224d039d7SMichael Lange   IS              rankIS,   pointIS;
165324d039d7SMichael Lange   const PetscInt *ranks,   *points;
165424d039d7SMichael Lange   PetscInt        numRanks, numPoints, r, p, a, adjSize;
165524d039d7SMichael Lange   PetscInt       *adj = NULL;
165670034214SMatthew G. Knepley   PetscErrorCode  ierr;
165770034214SMatthew G. Knepley 
165870034214SMatthew G. Knepley   PetscFunctionBegin;
165924d039d7SMichael Lange   ierr = DMLabelGetValueIS(label, &rankIS);CHKERRQ(ierr);
166024d039d7SMichael Lange   ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr);
166124d039d7SMichael Lange   ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr);
166224d039d7SMichael Lange   for (r = 0; r < numRanks; ++r) {
166324d039d7SMichael Lange     const PetscInt rank = ranks[r];
166470034214SMatthew G. Knepley 
166524d039d7SMichael Lange     ierr = DMLabelGetStratumIS(label, rank, &pointIS);CHKERRQ(ierr);
166624d039d7SMichael Lange     ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
166724d039d7SMichael Lange     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
166870034214SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
166924d039d7SMichael Lange       adjSize = PETSC_DETERMINE;
167024d039d7SMichael Lange       ierr = DMPlexGetAdjacency(dm, points[p], &adjSize, &adj);CHKERRQ(ierr);
167124d039d7SMichael Lange       for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(label, adj[a], rank);CHKERRQ(ierr);}
167270034214SMatthew G. Knepley     }
167324d039d7SMichael Lange     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
167424d039d7SMichael Lange     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
167570034214SMatthew G. Knepley   }
167624d039d7SMichael Lange   ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr);
167724d039d7SMichael Lange   ierr = ISDestroy(&rankIS);CHKERRQ(ierr);
167824d039d7SMichael Lange   ierr = PetscFree(adj);CHKERRQ(ierr);
167924d039d7SMichael Lange   PetscFunctionReturn(0);
168070034214SMatthew G. Knepley }
168170034214SMatthew G. Knepley 
1682be200f8dSMichael Lange /*@
1683be200f8dSMichael Lange   DMPlexPartitionLabelPropagate - Propagate points in a partition label over the point SF
1684be200f8dSMichael Lange 
1685be200f8dSMichael Lange   Input Parameters:
1686be200f8dSMichael Lange + dm     - The DM
1687be200f8dSMichael Lange - label  - DMLabel assinging ranks to remote roots
1688be200f8dSMichael Lange 
1689be200f8dSMichael Lange   Level: developer
1690be200f8dSMichael Lange 
1691be200f8dSMichael Lange   Note: This is required when generating multi-level overlaps to capture
1692be200f8dSMichael Lange   overlap points from non-neighbouring partitions.
1693be200f8dSMichael Lange 
1694be200f8dSMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap
1695be200f8dSMichael Lange @*/
1696be200f8dSMichael Lange PetscErrorCode DMPlexPartitionLabelPropagate(DM dm, DMLabel label)
1697be200f8dSMichael Lange {
1698be200f8dSMichael Lange   MPI_Comm        comm;
1699be200f8dSMichael Lange   PetscMPIInt     rank;
1700be200f8dSMichael Lange   PetscSF         sfPoint;
17015d04f6ebSMichael Lange   DMLabel         lblRoots, lblLeaves;
1702be200f8dSMichael Lange   IS              rankIS, pointIS;
1703be200f8dSMichael Lange   const PetscInt *ranks;
1704be200f8dSMichael Lange   PetscInt        numRanks, r;
1705be200f8dSMichael Lange   PetscErrorCode  ierr;
1706be200f8dSMichael Lange 
1707be200f8dSMichael Lange   PetscFunctionBegin;
1708be200f8dSMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
1709be200f8dSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1710be200f8dSMichael Lange   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
17115d04f6ebSMichael Lange   /* Pull point contributions from remote leaves into local roots */
17125d04f6ebSMichael Lange   ierr = DMLabelGather(label, sfPoint, &lblLeaves);CHKERRQ(ierr);
17135d04f6ebSMichael Lange   ierr = DMLabelGetValueIS(lblLeaves, &rankIS);CHKERRQ(ierr);
17145d04f6ebSMichael Lange   ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr);
17155d04f6ebSMichael Lange   ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr);
17165d04f6ebSMichael Lange   for (r = 0; r < numRanks; ++r) {
17175d04f6ebSMichael Lange     const PetscInt remoteRank = ranks[r];
17185d04f6ebSMichael Lange     if (remoteRank == rank) continue;
17195d04f6ebSMichael Lange     ierr = DMLabelGetStratumIS(lblLeaves, remoteRank, &pointIS);CHKERRQ(ierr);
17205d04f6ebSMichael Lange     ierr = DMLabelInsertIS(label, pointIS, remoteRank);CHKERRQ(ierr);
17215d04f6ebSMichael Lange     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
17225d04f6ebSMichael Lange   }
17235d04f6ebSMichael Lange   ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr);
17245d04f6ebSMichael Lange   ierr = ISDestroy(&rankIS);CHKERRQ(ierr);
17255d04f6ebSMichael Lange   ierr = DMLabelDestroy(&lblLeaves);CHKERRQ(ierr);
1726be200f8dSMichael Lange   /* Push point contributions from roots into remote leaves */
1727be200f8dSMichael Lange   ierr = DMLabelDistribute(label, sfPoint, &lblRoots);CHKERRQ(ierr);
1728be200f8dSMichael Lange   ierr = DMLabelGetValueIS(lblRoots, &rankIS);CHKERRQ(ierr);
1729be200f8dSMichael Lange   ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr);
1730be200f8dSMichael Lange   ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr);
1731be200f8dSMichael Lange   for (r = 0; r < numRanks; ++r) {
1732be200f8dSMichael Lange     const PetscInt remoteRank = ranks[r];
1733be200f8dSMichael Lange     if (remoteRank == rank) continue;
1734be200f8dSMichael Lange     ierr = DMLabelGetStratumIS(lblRoots, remoteRank, &pointIS);CHKERRQ(ierr);
1735be200f8dSMichael Lange     ierr = DMLabelInsertIS(label, pointIS, remoteRank);CHKERRQ(ierr);
1736be200f8dSMichael Lange     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
1737be200f8dSMichael Lange   }
1738be200f8dSMichael Lange   ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr);
1739be200f8dSMichael Lange   ierr = ISDestroy(&rankIS);CHKERRQ(ierr);
1740be200f8dSMichael Lange   ierr = DMLabelDestroy(&lblRoots);CHKERRQ(ierr);
1741be200f8dSMichael Lange   PetscFunctionReturn(0);
1742be200f8dSMichael Lange }
1743be200f8dSMichael Lange 
17441fd9873aSMichael Lange /*@
17451fd9873aSMichael Lange   DMPlexPartitionLabelInvert - Create a partition label of remote roots from a local root label
17461fd9873aSMichael Lange 
17471fd9873aSMichael Lange   Input Parameters:
17481fd9873aSMichael Lange + dm        - The DM
17491fd9873aSMichael Lange . rootLabel - DMLabel assinging ranks to local roots
17501fd9873aSMichael Lange . processSF - A star forest mapping into the local index on each remote rank
17511fd9873aSMichael Lange 
17521fd9873aSMichael Lange   Output Parameter:
17531fd9873aSMichael Lange - leafLabel - DMLabel assinging ranks to remote roots
17541fd9873aSMichael Lange 
17551fd9873aSMichael Lange   Note: The rootLabel defines a send pattern by mapping local points to remote target ranks. The
17561fd9873aSMichael Lange   resulting leafLabel is a receiver mapping of remote roots to their parent rank.
17571fd9873aSMichael Lange 
17581fd9873aSMichael Lange   Level: developer
17591fd9873aSMichael Lange 
17601fd9873aSMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap
17611fd9873aSMichael Lange @*/
17621fd9873aSMichael Lange PetscErrorCode DMPlexPartitionLabelInvert(DM dm, DMLabel rootLabel, PetscSF processSF, DMLabel leafLabel)
17631fd9873aSMichael Lange {
17641fd9873aSMichael Lange   MPI_Comm           comm;
17659852e123SBarry Smith   PetscMPIInt        rank, size;
17669852e123SBarry Smith   PetscInt           p, n, numNeighbors, ssize, l, nleaves;
17671fd9873aSMichael Lange   PetscSF            sfPoint;
17681fd9873aSMichael Lange   PetscSFNode       *rootPoints, *leafPoints;
17691fd9873aSMichael Lange   PetscSection       rootSection, leafSection;
17701fd9873aSMichael Lange   const PetscSFNode *remote;
17711fd9873aSMichael Lange   const PetscInt    *local, *neighbors;
17721fd9873aSMichael Lange   IS                 valueIS;
17731fd9873aSMichael Lange   PetscErrorCode     ierr;
17741fd9873aSMichael Lange 
17751fd9873aSMichael Lange   PetscFunctionBegin;
17761fd9873aSMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
17771fd9873aSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
17789852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
17791fd9873aSMichael Lange   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
17801fd9873aSMichael Lange 
17811fd9873aSMichael Lange   /* Convert to (point, rank) and use actual owners */
17821fd9873aSMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
17839852e123SBarry Smith   ierr = PetscSectionSetChart(rootSection, 0, size);CHKERRQ(ierr);
17841fd9873aSMichael Lange   ierr = DMLabelGetValueIS(rootLabel, &valueIS);CHKERRQ(ierr);
17851fd9873aSMichael Lange   ierr = ISGetLocalSize(valueIS, &numNeighbors);CHKERRQ(ierr);
17861fd9873aSMichael Lange   ierr = ISGetIndices(valueIS, &neighbors);CHKERRQ(ierr);
17871fd9873aSMichael Lange   for (n = 0; n < numNeighbors; ++n) {
17881fd9873aSMichael Lange     PetscInt numPoints;
17891fd9873aSMichael Lange 
17901fd9873aSMichael Lange     ierr = DMLabelGetStratumSize(rootLabel, neighbors[n], &numPoints);CHKERRQ(ierr);
17911fd9873aSMichael Lange     ierr = PetscSectionAddDof(rootSection, neighbors[n], numPoints);CHKERRQ(ierr);
17921fd9873aSMichael Lange   }
17931fd9873aSMichael Lange   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
17949852e123SBarry Smith   ierr = PetscSectionGetStorageSize(rootSection, &ssize);CHKERRQ(ierr);
17959852e123SBarry Smith   ierr = PetscMalloc1(ssize, &rootPoints);CHKERRQ(ierr);
17961fd9873aSMichael Lange   ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr);
17971fd9873aSMichael Lange   for (n = 0; n < numNeighbors; ++n) {
17981fd9873aSMichael Lange     IS              pointIS;
17991fd9873aSMichael Lange     const PetscInt *points;
18001fd9873aSMichael Lange     PetscInt        off, numPoints, p;
18011fd9873aSMichael Lange 
18021fd9873aSMichael Lange     ierr = PetscSectionGetOffset(rootSection, neighbors[n], &off);CHKERRQ(ierr);
18031fd9873aSMichael Lange     ierr = DMLabelGetStratumIS(rootLabel, neighbors[n], &pointIS);CHKERRQ(ierr);
18041fd9873aSMichael Lange     ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
18051fd9873aSMichael Lange     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
18061fd9873aSMichael Lange     for (p = 0; p < numPoints; ++p) {
1807f8987ae8SMichael Lange       if (local) {ierr = PetscFindInt(points[p], nleaves, local, &l);CHKERRQ(ierr);}
1808f8987ae8SMichael Lange       else       {l = -1;}
18091fd9873aSMichael Lange       if (l >= 0) {rootPoints[off+p] = remote[l];}
18101fd9873aSMichael Lange       else        {rootPoints[off+p].index = points[p]; rootPoints[off+p].rank = rank;}
18111fd9873aSMichael Lange     }
18121fd9873aSMichael Lange     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
18131fd9873aSMichael Lange     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
18141fd9873aSMichael Lange   }
18151fd9873aSMichael Lange   ierr = ISRestoreIndices(valueIS, &neighbors);CHKERRQ(ierr);
18161fd9873aSMichael Lange   ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
18171fd9873aSMichael Lange   /* Communicate overlap */
18181fd9873aSMichael Lange   ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr);
18191fd9873aSMichael Lange   ierr = DMPlexDistributeData(dm, processSF, rootSection, MPIU_2INT, rootPoints, leafSection, (void**) &leafPoints);CHKERRQ(ierr);
18201fd9873aSMichael Lange   /* Filter remote contributions (ovLeafPoints) into the overlapSF */
18219852e123SBarry Smith   ierr = PetscSectionGetStorageSize(leafSection, &ssize);CHKERRQ(ierr);
18229852e123SBarry Smith   for (p = 0; p < ssize; p++) {
18231fd9873aSMichael Lange     ierr = DMLabelSetValue(leafLabel, leafPoints[p].index, leafPoints[p].rank);CHKERRQ(ierr);
18241fd9873aSMichael Lange   }
18251fd9873aSMichael Lange   ierr = PetscFree(rootPoints);CHKERRQ(ierr);
18261fd9873aSMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
18271fd9873aSMichael Lange   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
18281fd9873aSMichael Lange   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
18291fd9873aSMichael Lange   PetscFunctionReturn(0);
18301fd9873aSMichael Lange }
18311fd9873aSMichael Lange 
1832aa3148a8SMichael Lange /*@
1833aa3148a8SMichael Lange   DMPlexPartitionLabelCreateSF - Create a star forest from a label that assigns ranks to points
1834aa3148a8SMichael Lange 
1835aa3148a8SMichael Lange   Input Parameters:
1836aa3148a8SMichael Lange + dm    - The DM
1837aa3148a8SMichael Lange . label - DMLabel assinging ranks to remote roots
1838aa3148a8SMichael Lange 
1839aa3148a8SMichael Lange   Output Parameter:
1840aa3148a8SMichael Lange - sf    - The star forest communication context encapsulating the defined mapping
1841aa3148a8SMichael Lange 
1842aa3148a8SMichael Lange   Note: The incoming label is a receiver mapping of remote points to their parent rank.
1843aa3148a8SMichael Lange 
1844aa3148a8SMichael Lange   Level: developer
1845aa3148a8SMichael Lange 
1846aa3148a8SMichael Lange .seealso: DMPlexDistribute(), DMPlexCreateOverlap
1847aa3148a8SMichael Lange @*/
1848aa3148a8SMichael Lange PetscErrorCode DMPlexPartitionLabelCreateSF(DM dm, DMLabel label, PetscSF *sf)
1849aa3148a8SMichael Lange {
18509852e123SBarry Smith   PetscMPIInt     rank, size;
185143f7d02bSMichael Lange   PetscInt        n, numRemote, p, numPoints, pStart, pEnd, idx = 0;
1852aa3148a8SMichael Lange   PetscSFNode    *remotePoints;
185343f7d02bSMichael Lange   IS              remoteRootIS;
185443f7d02bSMichael Lange   const PetscInt *remoteRoots;
1855aa3148a8SMichael Lange   PetscErrorCode ierr;
1856aa3148a8SMichael Lange 
1857aa3148a8SMichael Lange   PetscFunctionBegin;
185843f7d02bSMichael Lange   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
18599852e123SBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size);CHKERRQ(ierr);
1860aa3148a8SMichael Lange 
18619852e123SBarry Smith   for (numRemote = 0, n = 0; n < size; ++n) {
1862aa3148a8SMichael Lange     ierr = DMLabelGetStratumSize(label, n, &numPoints);CHKERRQ(ierr);
1863aa3148a8SMichael Lange     numRemote += numPoints;
1864aa3148a8SMichael Lange   }
1865aa3148a8SMichael Lange   ierr = PetscMalloc1(numRemote, &remotePoints);CHKERRQ(ierr);
186643f7d02bSMichael Lange   /* Put owned points first */
186743f7d02bSMichael Lange   ierr = DMLabelGetStratumSize(label, rank, &numPoints);CHKERRQ(ierr);
186843f7d02bSMichael Lange   if (numPoints > 0) {
186943f7d02bSMichael Lange     ierr = DMLabelGetStratumIS(label, rank, &remoteRootIS);CHKERRQ(ierr);
187043f7d02bSMichael Lange     ierr = ISGetIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
187143f7d02bSMichael Lange     for (p = 0; p < numPoints; p++) {
187243f7d02bSMichael Lange       remotePoints[idx].index = remoteRoots[p];
187343f7d02bSMichael Lange       remotePoints[idx].rank = rank;
187443f7d02bSMichael Lange       idx++;
187543f7d02bSMichael Lange     }
187643f7d02bSMichael Lange     ierr = ISRestoreIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
187743f7d02bSMichael Lange     ierr = ISDestroy(&remoteRootIS);CHKERRQ(ierr);
187843f7d02bSMichael Lange   }
187943f7d02bSMichael Lange   /* Now add remote points */
18809852e123SBarry Smith   for (n = 0; n < size; ++n) {
1881aa3148a8SMichael Lange     ierr = DMLabelGetStratumSize(label, n, &numPoints);CHKERRQ(ierr);
188243f7d02bSMichael Lange     if (numPoints <= 0 || n == rank) continue;
1883aa3148a8SMichael Lange     ierr = DMLabelGetStratumIS(label, n, &remoteRootIS);CHKERRQ(ierr);
1884aa3148a8SMichael Lange     ierr = ISGetIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
1885aa3148a8SMichael Lange     for (p = 0; p < numPoints; p++) {
1886aa3148a8SMichael Lange       remotePoints[idx].index = remoteRoots[p];
1887aa3148a8SMichael Lange       remotePoints[idx].rank = n;
1888aa3148a8SMichael Lange       idx++;
1889aa3148a8SMichael Lange     }
1890aa3148a8SMichael Lange     ierr = ISRestoreIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
1891aa3148a8SMichael Lange     ierr = ISDestroy(&remoteRootIS);CHKERRQ(ierr);
1892aa3148a8SMichael Lange   }
1893aa3148a8SMichael Lange   ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), sf);CHKERRQ(ierr);
1894aa3148a8SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
1895aa3148a8SMichael Lange   ierr = PetscSFSetGraph(*sf, pEnd-pStart, numRemote, NULL, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
189670034214SMatthew G. Knepley   PetscFunctionReturn(0);
189770034214SMatthew G. Knepley }
1898