xref: /petsc/src/dm/impls/plex/plexpartition.c (revision f8987ae87a5d7b49db17cc3bf536b2c44eec4fd8)
170034214SMatthew G. Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
270034214SMatthew G. Knepley 
377623264SMatthew G. Knepley PetscClassId PETSCPARTITIONER_CLASSID = 0;
477623264SMatthew G. Knepley 
577623264SMatthew G. Knepley PetscFunctionList PetscPartitionerList              = NULL;
677623264SMatthew G. Knepley PetscBool         PetscPartitionerRegisterAllCalled = PETSC_FALSE;
777623264SMatthew G. Knepley 
877623264SMatthew G. Knepley PetscBool ChacoPartitionercite = PETSC_FALSE;
977623264SMatthew G. Knepley const char ChacoPartitionerCitation[] = "@inproceedings{Chaco95,\n"
1077623264SMatthew G. Knepley                                "  author    = {Bruce Hendrickson and Robert Leland},\n"
1177623264SMatthew G. Knepley                                "  title     = {A multilevel algorithm for partitioning graphs},\n"
1277623264SMatthew G. Knepley                                "  booktitle = {Supercomputing '95: Proceedings of the 1995 ACM/IEEE Conference on Supercomputing (CDROM)},"
1377623264SMatthew G. Knepley                                "  isbn      = {0-89791-816-9},\n"
1477623264SMatthew G. Knepley                                "  pages     = {28},\n"
1577623264SMatthew G. Knepley                                "  doi       = {http://doi.acm.org/10.1145/224170.224228},\n"
1677623264SMatthew G. Knepley                                "  publisher = {ACM Press},\n"
1777623264SMatthew G. Knepley                                "  address   = {New York},\n"
1877623264SMatthew G. Knepley                                "  year      = {1995}\n}\n";
1977623264SMatthew G. Knepley 
2077623264SMatthew G. Knepley PetscBool ParMetisPartitionercite = PETSC_FALSE;
2177623264SMatthew G. Knepley const char ParMetisPartitionerCitation[] = "@article{KarypisKumar98,\n"
2277623264SMatthew G. Knepley                                "  author  = {George Karypis and Vipin Kumar},\n"
2377623264SMatthew G. Knepley                                "  title   = {A Parallel Algorithm for Multilevel Graph Partitioning and Sparse Matrix Ordering},\n"
2477623264SMatthew G. Knepley                                "  journal = {Journal of Parallel and Distributed Computing},\n"
2577623264SMatthew G. Knepley                                "  volume  = {48},\n"
2677623264SMatthew G. Knepley                                "  pages   = {71--85},\n"
2777623264SMatthew G. Knepley                                "  year    = {1998}\n}\n";
2877623264SMatthew G. Knepley 
2970034214SMatthew G. Knepley #undef __FUNCT__
3070034214SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateNeighborCSR"
3170034214SMatthew G. Knepley PetscErrorCode DMPlexCreateNeighborCSR(DM dm, PetscInt cellHeight, PetscInt *numVertices, PetscInt **offsets, PetscInt **adjacency)
3270034214SMatthew G. Knepley {
3370034214SMatthew G. Knepley   const PetscInt maxFaceCases = 30;
3470034214SMatthew G. Knepley   PetscInt       numFaceCases = 0;
3570034214SMatthew G. Knepley   PetscInt       numFaceVertices[30]; /* maxFaceCases, C89 sucks sucks sucks */
3670034214SMatthew G. Knepley   PetscInt      *off, *adj;
3770034214SMatthew G. Knepley   PetscInt      *neighborCells = NULL;
3870034214SMatthew G. Knepley   PetscInt       dim, cellDim, depth = 0, faceDepth, cStart, cEnd, c, numCells, cell;
3970034214SMatthew G. Knepley   PetscErrorCode ierr;
4070034214SMatthew G. Knepley 
4170034214SMatthew G. Knepley   PetscFunctionBegin;
4270034214SMatthew G. Knepley   /* For parallel partitioning, I think you have to communicate supports */
43c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
4470034214SMatthew G. Knepley   cellDim = dim - cellHeight;
4570034214SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
4670034214SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
4770034214SMatthew G. Knepley   if (cEnd - cStart == 0) {
4870034214SMatthew G. Knepley     if (numVertices) *numVertices = 0;
4970034214SMatthew G. Knepley     if (offsets)   *offsets   = NULL;
5070034214SMatthew G. Knepley     if (adjacency) *adjacency = NULL;
5170034214SMatthew G. Knepley     PetscFunctionReturn(0);
5270034214SMatthew G. Knepley   }
5370034214SMatthew G. Knepley   numCells  = cEnd - cStart;
5470034214SMatthew G. Knepley   faceDepth = depth - cellHeight;
5570034214SMatthew G. Knepley   if (dim == depth) {
5670034214SMatthew G. Knepley     PetscInt f, fStart, fEnd;
5770034214SMatthew G. Knepley 
5870034214SMatthew G. Knepley     ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr);
5970034214SMatthew G. Knepley     /* Count neighboring cells */
6070034214SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, cellHeight+1, &fStart, &fEnd);CHKERRQ(ierr);
6170034214SMatthew G. Knepley     for (f = fStart; f < fEnd; ++f) {
6270034214SMatthew G. Knepley       const PetscInt *support;
6370034214SMatthew G. Knepley       PetscInt        supportSize;
6470034214SMatthew G. Knepley       ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr);
6570034214SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr);
6670034214SMatthew G. Knepley       if (supportSize == 2) {
679ffc88e4SToby Isaac         PetscInt numChildren;
689ffc88e4SToby Isaac 
699ffc88e4SToby Isaac         ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr);
709ffc88e4SToby Isaac         if (!numChildren) {
7170034214SMatthew G. Knepley           ++off[support[0]-cStart+1];
7270034214SMatthew G. Knepley           ++off[support[1]-cStart+1];
7370034214SMatthew G. Knepley         }
7470034214SMatthew G. Knepley       }
759ffc88e4SToby Isaac     }
7670034214SMatthew G. Knepley     /* Prefix sum */
7770034214SMatthew G. Knepley     for (c = 1; c <= numCells; ++c) off[c] += off[c-1];
7870034214SMatthew G. Knepley     if (adjacency) {
7970034214SMatthew G. Knepley       PetscInt *tmp;
8070034214SMatthew G. Knepley 
8170034214SMatthew G. Knepley       ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr);
8270034214SMatthew G. Knepley       ierr = PetscMalloc1((numCells+1), &tmp);CHKERRQ(ierr);
8370034214SMatthew G. Knepley       ierr = PetscMemcpy(tmp, off, (numCells+1) * sizeof(PetscInt));CHKERRQ(ierr);
8470034214SMatthew G. Knepley       /* Get neighboring cells */
8570034214SMatthew G. Knepley       for (f = fStart; f < fEnd; ++f) {
8670034214SMatthew G. Knepley         const PetscInt *support;
8770034214SMatthew G. Knepley         PetscInt        supportSize;
8870034214SMatthew G. Knepley         ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr);
8970034214SMatthew G. Knepley         ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr);
9070034214SMatthew G. Knepley         if (supportSize == 2) {
919ffc88e4SToby Isaac           PetscInt numChildren;
929ffc88e4SToby Isaac 
939ffc88e4SToby Isaac           ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr);
949ffc88e4SToby Isaac           if (!numChildren) {
9570034214SMatthew G. Knepley             adj[tmp[support[0]-cStart]++] = support[1];
9670034214SMatthew G. Knepley             adj[tmp[support[1]-cStart]++] = support[0];
9770034214SMatthew G. Knepley           }
9870034214SMatthew G. Knepley         }
999ffc88e4SToby Isaac       }
10070034214SMatthew G. Knepley #if defined(PETSC_USE_DEBUG)
10170034214SMatthew 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);
10270034214SMatthew G. Knepley #endif
10370034214SMatthew G. Knepley       ierr = PetscFree(tmp);CHKERRQ(ierr);
10470034214SMatthew G. Knepley     }
10570034214SMatthew G. Knepley     if (numVertices) *numVertices = numCells;
10670034214SMatthew G. Knepley     if (offsets)   *offsets   = off;
10770034214SMatthew G. Knepley     if (adjacency) *adjacency = adj;
10870034214SMatthew G. Knepley     PetscFunctionReturn(0);
10970034214SMatthew G. Knepley   }
11070034214SMatthew G. Knepley   /* Setup face recognition */
11170034214SMatthew G. Knepley   if (faceDepth == 1) {
11270034214SMatthew 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 */
11370034214SMatthew G. Knepley 
11470034214SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
11570034214SMatthew G. Knepley       PetscInt corners;
11670034214SMatthew G. Knepley 
11770034214SMatthew G. Knepley       ierr = DMPlexGetConeSize(dm, c, &corners);CHKERRQ(ierr);
11870034214SMatthew G. Knepley       if (!cornersSeen[corners]) {
11970034214SMatthew G. Knepley         PetscInt nFV;
12070034214SMatthew G. Knepley 
12170034214SMatthew G. Knepley         if (numFaceCases >= maxFaceCases) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Exceeded maximum number of face recognition cases");
12270034214SMatthew G. Knepley         cornersSeen[corners] = 1;
12370034214SMatthew G. Knepley 
12470034214SMatthew G. Knepley         ierr = DMPlexGetNumFaceVertices(dm, cellDim, corners, &nFV);CHKERRQ(ierr);
12570034214SMatthew G. Knepley 
12670034214SMatthew G. Knepley         numFaceVertices[numFaceCases++] = nFV;
12770034214SMatthew G. Knepley       }
12870034214SMatthew G. Knepley     }
12970034214SMatthew G. Knepley   }
13070034214SMatthew G. Knepley   ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr);
13170034214SMatthew G. Knepley   /* Count neighboring cells */
13270034214SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) {
13370034214SMatthew G. Knepley     PetscInt numNeighbors = PETSC_DETERMINE, n;
13470034214SMatthew G. Knepley 
1358b0b4c70SToby Isaac     ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &numNeighbors, &neighborCells);CHKERRQ(ierr);
13670034214SMatthew G. Knepley     /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */
13770034214SMatthew G. Knepley     for (n = 0; n < numNeighbors; ++n) {
13870034214SMatthew G. Knepley       PetscInt        cellPair[2];
13970034214SMatthew G. Knepley       PetscBool       found    = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE;
14070034214SMatthew G. Knepley       PetscInt        meetSize = 0;
14170034214SMatthew G. Knepley       const PetscInt *meet    = NULL;
14270034214SMatthew G. Knepley 
14370034214SMatthew G. Knepley       cellPair[0] = cell; cellPair[1] = neighborCells[n];
14470034214SMatthew G. Knepley       if (cellPair[0] == cellPair[1]) continue;
14570034214SMatthew G. Knepley       if (!found) {
14670034214SMatthew G. Knepley         ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
14770034214SMatthew G. Knepley         if (meetSize) {
14870034214SMatthew G. Knepley           PetscInt f;
14970034214SMatthew G. Knepley 
15070034214SMatthew G. Knepley           for (f = 0; f < numFaceCases; ++f) {
15170034214SMatthew G. Knepley             if (numFaceVertices[f] == meetSize) {
15270034214SMatthew G. Knepley               found = PETSC_TRUE;
15370034214SMatthew G. Knepley               break;
15470034214SMatthew G. Knepley             }
15570034214SMatthew G. Knepley           }
15670034214SMatthew G. Knepley         }
15770034214SMatthew G. Knepley         ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
15870034214SMatthew G. Knepley       }
15970034214SMatthew G. Knepley       if (found) ++off[cell-cStart+1];
16070034214SMatthew G. Knepley     }
16170034214SMatthew G. Knepley   }
16270034214SMatthew G. Knepley   /* Prefix sum */
16370034214SMatthew G. Knepley   for (cell = 1; cell <= numCells; ++cell) off[cell] += off[cell-1];
16470034214SMatthew G. Knepley 
16570034214SMatthew G. Knepley   if (adjacency) {
16670034214SMatthew G. Knepley     ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr);
16770034214SMatthew G. Knepley     /* Get neighboring cells */
16870034214SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
16970034214SMatthew G. Knepley       PetscInt numNeighbors = PETSC_DETERMINE, n;
17070034214SMatthew G. Knepley       PetscInt cellOffset   = 0;
17170034214SMatthew G. Knepley 
1728b0b4c70SToby Isaac       ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &numNeighbors, &neighborCells);CHKERRQ(ierr);
17370034214SMatthew G. Knepley       /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */
17470034214SMatthew G. Knepley       for (n = 0; n < numNeighbors; ++n) {
17570034214SMatthew G. Knepley         PetscInt        cellPair[2];
17670034214SMatthew G. Knepley         PetscBool       found    = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE;
17770034214SMatthew G. Knepley         PetscInt        meetSize = 0;
17870034214SMatthew G. Knepley         const PetscInt *meet    = NULL;
17970034214SMatthew G. Knepley 
18070034214SMatthew G. Knepley         cellPair[0] = cell; cellPair[1] = neighborCells[n];
18170034214SMatthew G. Knepley         if (cellPair[0] == cellPair[1]) continue;
18270034214SMatthew G. Knepley         if (!found) {
18370034214SMatthew G. Knepley           ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
18470034214SMatthew G. Knepley           if (meetSize) {
18570034214SMatthew G. Knepley             PetscInt f;
18670034214SMatthew G. Knepley 
18770034214SMatthew G. Knepley             for (f = 0; f < numFaceCases; ++f) {
18870034214SMatthew G. Knepley               if (numFaceVertices[f] == meetSize) {
18970034214SMatthew G. Knepley                 found = PETSC_TRUE;
19070034214SMatthew G. Knepley                 break;
19170034214SMatthew G. Knepley               }
19270034214SMatthew G. Knepley             }
19370034214SMatthew G. Knepley           }
19470034214SMatthew G. Knepley           ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
19570034214SMatthew G. Knepley         }
19670034214SMatthew G. Knepley         if (found) {
19770034214SMatthew G. Knepley           adj[off[cell-cStart]+cellOffset] = neighborCells[n];
19870034214SMatthew G. Knepley           ++cellOffset;
19970034214SMatthew G. Knepley         }
20070034214SMatthew G. Knepley       }
20170034214SMatthew G. Knepley     }
20270034214SMatthew G. Knepley   }
20370034214SMatthew G. Knepley   ierr = PetscFree(neighborCells);CHKERRQ(ierr);
20470034214SMatthew G. Knepley   if (numVertices) *numVertices = numCells;
20570034214SMatthew G. Knepley   if (offsets)   *offsets   = off;
20670034214SMatthew G. Knepley   if (adjacency) *adjacency = adj;
20770034214SMatthew G. Knepley   PetscFunctionReturn(0);
20870034214SMatthew G. Knepley }
20970034214SMatthew G. Knepley 
21077623264SMatthew G. Knepley #undef __FUNCT__
21177623264SMatthew G. Knepley #define __FUNCT__ "DMPlexEnlargePartition"
21277623264SMatthew G. Knepley /* Expand the partition by BFS on the adjacency graph */
21377623264SMatthew G. Knepley PetscErrorCode DMPlexEnlargePartition(DM dm, const PetscInt start[], const PetscInt adjacency[], PetscSection origPartSection, IS origPartition, PetscSection partSection, IS *partition)
21477623264SMatthew G. Knepley {
21577623264SMatthew G. Knepley   PetscHashI      h;
21677623264SMatthew G. Knepley   const PetscInt *points;
21777623264SMatthew G. Knepley   PetscInt      **tmpPoints, *newPoints, totPoints = 0;
21877623264SMatthew G. Knepley   PetscInt        pStart, pEnd, part, q;
21977623264SMatthew G. Knepley   PetscBool       useCone;
22077623264SMatthew G. Knepley   PetscErrorCode  ierr;
22177623264SMatthew G. Knepley 
22277623264SMatthew G. Knepley   PetscFunctionBegin;
22377623264SMatthew G. Knepley   PetscHashICreate(h);
22477623264SMatthew G. Knepley   ierr = PetscSectionGetChart(origPartSection, &pStart, &pEnd);CHKERRQ(ierr);
22577623264SMatthew G. Knepley   ierr = PetscSectionSetChart(partSection, pStart, pEnd);CHKERRQ(ierr);
22677623264SMatthew G. Knepley   ierr = ISGetIndices(origPartition, &points);CHKERRQ(ierr);
22777623264SMatthew G. Knepley   ierr = PetscMalloc1((pEnd - pStart), &tmpPoints);CHKERRQ(ierr);
22877623264SMatthew G. Knepley   ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
22977623264SMatthew G. Knepley   ierr = DMPlexSetAdjacencyUseCone(dm, PETSC_TRUE);CHKERRQ(ierr);
23077623264SMatthew G. Knepley   for (part = pStart; part < pEnd; ++part) {
23177623264SMatthew G. Knepley     PetscInt *adj = NULL;
23277623264SMatthew G. Knepley     PetscInt  numPoints, nP, numNewPoints, off, p, n = 0;
23377623264SMatthew G. Knepley 
23477623264SMatthew G. Knepley     PetscHashIClear(h);
23577623264SMatthew G. Knepley     ierr = PetscSectionGetDof(origPartSection, part, &numPoints);CHKERRQ(ierr);
23677623264SMatthew G. Knepley     ierr = PetscSectionGetOffset(origPartSection, part, &off);CHKERRQ(ierr);
23777623264SMatthew G. Knepley     /* Add all existing points to h */
23877623264SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
23977623264SMatthew G. Knepley       const PetscInt point = points[off+p];
24077623264SMatthew G. Knepley       PetscHashIAdd(h, point, 1);
24177623264SMatthew G. Knepley     }
24277623264SMatthew G. Knepley     PetscHashISize(h, nP);
24377623264SMatthew G. Knepley     if (nP != numPoints) SETERRQ2(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONG, "Invalid partition has %d points, but only %d were unique", numPoints, nP);
24477623264SMatthew G. Knepley     /* Add all points in next BFS level */
24577623264SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
24677623264SMatthew G. Knepley       const PetscInt point   = points[off+p];
24777623264SMatthew G. Knepley       PetscInt       adjSize = PETSC_DETERMINE, a;
24877623264SMatthew G. Knepley 
24977623264SMatthew G. Knepley       ierr = DMPlexGetAdjacency(dm, point, &adjSize, &adj);CHKERRQ(ierr);
25077623264SMatthew G. Knepley       for (a = 0; a < adjSize; ++a) PetscHashIAdd(h, adj[a], 1);
25177623264SMatthew G. Knepley     }
25277623264SMatthew G. Knepley     PetscHashISize(h, numNewPoints);
25377623264SMatthew G. Knepley     ierr = PetscSectionSetDof(partSection, part, numNewPoints);CHKERRQ(ierr);
25477623264SMatthew G. Knepley     ierr = PetscMalloc1(numNewPoints, &tmpPoints[part]);CHKERRQ(ierr);
25577623264SMatthew G. Knepley     ierr = PetscHashIGetKeys(h, &n, tmpPoints[part]);CHKERRQ(ierr);
25677623264SMatthew G. Knepley     ierr = PetscFree(adj);CHKERRQ(ierr);
25777623264SMatthew G. Knepley     totPoints += numNewPoints;
25877623264SMatthew G. Knepley   }
25977623264SMatthew G. Knepley   ierr = DMPlexSetAdjacencyUseCone(dm, useCone);CHKERRQ(ierr);
26077623264SMatthew G. Knepley   ierr = ISRestoreIndices(origPartition, &points);CHKERRQ(ierr);
26177623264SMatthew G. Knepley   PetscHashIDestroy(h);
26277623264SMatthew G. Knepley   ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
26377623264SMatthew G. Knepley   ierr = PetscMalloc1(totPoints, &newPoints);CHKERRQ(ierr);
26477623264SMatthew G. Knepley   for (part = pStart, q = 0; part < pEnd; ++part) {
26577623264SMatthew G. Knepley     PetscInt numPoints, p;
26677623264SMatthew G. Knepley 
26777623264SMatthew G. Knepley     ierr = PetscSectionGetDof(partSection, part, &numPoints);CHKERRQ(ierr);
26877623264SMatthew G. Knepley     for (p = 0; p < numPoints; ++p, ++q) newPoints[q] = tmpPoints[part][p];
26977623264SMatthew G. Knepley     ierr = PetscFree(tmpPoints[part]);CHKERRQ(ierr);
27077623264SMatthew G. Knepley   }
27177623264SMatthew G. Knepley   ierr = PetscFree(tmpPoints);CHKERRQ(ierr);
27277623264SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), totPoints, newPoints, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
27377623264SMatthew G. Knepley   PetscFunctionReturn(0);
27477623264SMatthew G. Knepley }
27577623264SMatthew G. Knepley 
27677623264SMatthew G. Knepley #undef __FUNCT__
27777623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerRegister"
27877623264SMatthew G. Knepley /*@C
27977623264SMatthew G. Knepley   PetscPartitionerRegister - Adds a new PetscPartitioner implementation
28077623264SMatthew G. Knepley 
28177623264SMatthew G. Knepley   Not Collective
28277623264SMatthew G. Knepley 
28377623264SMatthew G. Knepley   Input Parameters:
28477623264SMatthew G. Knepley + name        - The name of a new user-defined creation routine
28577623264SMatthew G. Knepley - create_func - The creation routine itself
28677623264SMatthew G. Knepley 
28777623264SMatthew G. Knepley   Notes:
28877623264SMatthew G. Knepley   PetscPartitionerRegister() may be called multiple times to add several user-defined PetscPartitioners
28977623264SMatthew G. Knepley 
29077623264SMatthew G. Knepley   Sample usage:
29177623264SMatthew G. Knepley .vb
29277623264SMatthew G. Knepley     PetscPartitionerRegister("my_part", MyPetscPartitionerCreate);
29377623264SMatthew G. Knepley .ve
29477623264SMatthew G. Knepley 
29577623264SMatthew G. Knepley   Then, your PetscPartitioner type can be chosen with the procedural interface via
29677623264SMatthew G. Knepley .vb
29777623264SMatthew G. Knepley     PetscPartitionerCreate(MPI_Comm, PetscPartitioner *);
29877623264SMatthew G. Knepley     PetscPartitionerSetType(PetscPartitioner, "my_part");
29977623264SMatthew G. Knepley .ve
30077623264SMatthew G. Knepley    or at runtime via the option
30177623264SMatthew G. Knepley .vb
30277623264SMatthew G. Knepley     -petscpartitioner_type my_part
30377623264SMatthew G. Knepley .ve
30477623264SMatthew G. Knepley 
30577623264SMatthew G. Knepley   Level: advanced
30677623264SMatthew G. Knepley 
30777623264SMatthew G. Knepley .keywords: PetscPartitioner, register
30877623264SMatthew G. Knepley .seealso: PetscPartitionerRegisterAll(), PetscPartitionerRegisterDestroy()
30977623264SMatthew G. Knepley 
31077623264SMatthew G. Knepley @*/
31177623264SMatthew G. Knepley PetscErrorCode PetscPartitionerRegister(const char sname[], PetscErrorCode (*function)(PetscPartitioner))
31277623264SMatthew G. Knepley {
31377623264SMatthew G. Knepley   PetscErrorCode ierr;
31477623264SMatthew G. Knepley 
31577623264SMatthew G. Knepley   PetscFunctionBegin;
31677623264SMatthew G. Knepley   ierr = PetscFunctionListAdd(&PetscPartitionerList, sname, function);CHKERRQ(ierr);
31777623264SMatthew G. Knepley   PetscFunctionReturn(0);
31877623264SMatthew G. Knepley }
31977623264SMatthew G. Knepley 
32077623264SMatthew G. Knepley #undef __FUNCT__
32177623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerSetType"
32277623264SMatthew G. Knepley /*@C
32377623264SMatthew G. Knepley   PetscPartitionerSetType - Builds a particular PetscPartitioner
32477623264SMatthew G. Knepley 
32577623264SMatthew G. Knepley   Collective on PetscPartitioner
32677623264SMatthew G. Knepley 
32777623264SMatthew G. Knepley   Input Parameters:
32877623264SMatthew G. Knepley + part - The PetscPartitioner object
32977623264SMatthew G. Knepley - name - The kind of partitioner
33077623264SMatthew G. Knepley 
33177623264SMatthew G. Knepley   Options Database Key:
33277623264SMatthew G. Knepley . -petscpartitioner_type <type> - Sets the PetscPartitioner type; use -help for a list of available types
33377623264SMatthew G. Knepley 
33477623264SMatthew G. Knepley   Level: intermediate
33577623264SMatthew G. Knepley 
33677623264SMatthew G. Knepley .keywords: PetscPartitioner, set, type
33777623264SMatthew G. Knepley .seealso: PetscPartitionerGetType(), PetscPartitionerCreate()
33877623264SMatthew G. Knepley @*/
33977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetType(PetscPartitioner part, PetscPartitionerType name)
34077623264SMatthew G. Knepley {
34177623264SMatthew G. Knepley   PetscErrorCode (*r)(PetscPartitioner);
34277623264SMatthew G. Knepley   PetscBool      match;
34377623264SMatthew G. Knepley   PetscErrorCode ierr;
34477623264SMatthew G. Knepley 
34577623264SMatthew G. Knepley   PetscFunctionBegin;
34677623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
34777623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) part, name, &match);CHKERRQ(ierr);
34877623264SMatthew G. Knepley   if (match) PetscFunctionReturn(0);
34977623264SMatthew G. Knepley 
35077623264SMatthew G. Knepley   if (!PetscPartitionerRegisterAllCalled) {ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr);}
35177623264SMatthew G. Knepley   ierr = PetscFunctionListFind(PetscPartitionerList, name, &r);CHKERRQ(ierr);
35277623264SMatthew G. Knepley   if (!r) SETERRQ1(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscPartitioner type: %s", name);
35377623264SMatthew G. Knepley 
35477623264SMatthew G. Knepley   if (part->ops->destroy) {
35577623264SMatthew G. Knepley     ierr              = (*part->ops->destroy)(part);CHKERRQ(ierr);
35677623264SMatthew G. Knepley     part->ops->destroy = NULL;
35777623264SMatthew G. Knepley   }
35877623264SMatthew G. Knepley   ierr = (*r)(part);CHKERRQ(ierr);
35977623264SMatthew G. Knepley   ierr = PetscObjectChangeTypeName((PetscObject) part, name);CHKERRQ(ierr);
36077623264SMatthew G. Knepley   PetscFunctionReturn(0);
36177623264SMatthew G. Knepley }
36277623264SMatthew G. Knepley 
36377623264SMatthew G. Knepley #undef __FUNCT__
36477623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerGetType"
36577623264SMatthew G. Knepley /*@C
36677623264SMatthew G. Knepley   PetscPartitionerGetType - Gets the PetscPartitioner type name (as a string) from the object.
36777623264SMatthew G. Knepley 
36877623264SMatthew G. Knepley   Not Collective
36977623264SMatthew G. Knepley 
37077623264SMatthew G. Knepley   Input Parameter:
37177623264SMatthew G. Knepley . part - The PetscPartitioner
37277623264SMatthew G. Knepley 
37377623264SMatthew G. Knepley   Output Parameter:
37477623264SMatthew G. Knepley . name - The PetscPartitioner type name
37577623264SMatthew G. Knepley 
37677623264SMatthew G. Knepley   Level: intermediate
37777623264SMatthew G. Knepley 
37877623264SMatthew G. Knepley .keywords: PetscPartitioner, get, type, name
37977623264SMatthew G. Knepley .seealso: PetscPartitionerSetType(), PetscPartitionerCreate()
38077623264SMatthew G. Knepley @*/
38177623264SMatthew G. Knepley PetscErrorCode PetscPartitionerGetType(PetscPartitioner part, PetscPartitionerType *name)
38277623264SMatthew G. Knepley {
38377623264SMatthew G. Knepley   PetscErrorCode ierr;
38477623264SMatthew G. Knepley 
38577623264SMatthew G. Knepley   PetscFunctionBegin;
38677623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
38777623264SMatthew G. Knepley   PetscValidPointer(name, 2);
38877623264SMatthew G. Knepley   if (!PetscPartitionerRegisterAllCalled) {ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr);}
38977623264SMatthew G. Knepley   *name = ((PetscObject) part)->type_name;
39077623264SMatthew G. Knepley   PetscFunctionReturn(0);
39177623264SMatthew G. Knepley }
39277623264SMatthew G. Knepley 
39377623264SMatthew G. Knepley #undef __FUNCT__
39477623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView"
39577623264SMatthew G. Knepley /*@C
39677623264SMatthew G. Knepley   PetscPartitionerView - Views a PetscPartitioner
39777623264SMatthew G. Knepley 
39877623264SMatthew G. Knepley   Collective on PetscPartitioner
39977623264SMatthew G. Knepley 
40077623264SMatthew G. Knepley   Input Parameter:
40177623264SMatthew G. Knepley + part - the PetscPartitioner object to view
40277623264SMatthew G. Knepley - v    - the viewer
40377623264SMatthew G. Knepley 
40477623264SMatthew G. Knepley   Level: developer
40577623264SMatthew G. Knepley 
40677623264SMatthew G. Knepley .seealso: PetscPartitionerDestroy()
40777623264SMatthew G. Knepley @*/
40877623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView(PetscPartitioner part, PetscViewer v)
40977623264SMatthew G. Knepley {
41077623264SMatthew G. Knepley   PetscErrorCode ierr;
41177623264SMatthew G. Knepley 
41277623264SMatthew G. Knepley   PetscFunctionBegin;
41377623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
41477623264SMatthew G. Knepley   if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) part), &v);CHKERRQ(ierr);}
41577623264SMatthew G. Knepley   if (part->ops->view) {ierr = (*part->ops->view)(part, v);CHKERRQ(ierr);}
41677623264SMatthew G. Knepley   PetscFunctionReturn(0);
41777623264SMatthew G. Knepley }
41877623264SMatthew G. Knepley 
41977623264SMatthew G. Knepley #undef __FUNCT__
42077623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerViewFromOptions"
42177623264SMatthew G. Knepley /*
42277623264SMatthew G. Knepley   PetscPartitionerViewFromOptions - Processes command line options to determine if/how a PetscPartitioner is to be viewed.
42377623264SMatthew G. Knepley 
42477623264SMatthew G. Knepley   Collective on PetscPartitioner
42577623264SMatthew G. Knepley 
42677623264SMatthew G. Knepley   Input Parameters:
42777623264SMatthew G. Knepley + part   - the PetscPartitioner
42877623264SMatthew G. Knepley . prefix - prefix to use for viewing, or NULL
42977623264SMatthew G. Knepley - optionname - option to activate viewing
43077623264SMatthew G. Knepley 
43177623264SMatthew G. Knepley   Level: intermediate
43277623264SMatthew G. Knepley 
43377623264SMatthew G. Knepley .keywords: PetscPartitioner, view, options, database
43477623264SMatthew G. Knepley .seealso: VecViewFromOptions(), MatViewFromOptions()
43577623264SMatthew G. Knepley */
43677623264SMatthew G. Knepley PetscErrorCode PetscPartitionerViewFromOptions(PetscPartitioner part, const char prefix[], const char optionname[])
43777623264SMatthew G. Knepley {
43877623264SMatthew G. Knepley   PetscViewer       viewer;
43977623264SMatthew G. Knepley   PetscViewerFormat format;
44077623264SMatthew G. Knepley   PetscBool         flg;
44177623264SMatthew G. Knepley   PetscErrorCode    ierr;
44277623264SMatthew G. Knepley 
44377623264SMatthew G. Knepley   PetscFunctionBegin;
44477623264SMatthew G. Knepley   if (prefix) {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) part), prefix,                       optionname, &viewer, &format, &flg);CHKERRQ(ierr);}
44577623264SMatthew G. Knepley   else        {ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject) part), ((PetscObject) part)->prefix, optionname, &viewer, &format, &flg);CHKERRQ(ierr);}
44677623264SMatthew G. Knepley   if (flg) {
44777623264SMatthew G. Knepley     ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr);
44877623264SMatthew G. Knepley     ierr = PetscPartitionerView(part, viewer);CHKERRQ(ierr);
44977623264SMatthew G. Knepley     ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr);
45077623264SMatthew G. Knepley     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
45177623264SMatthew G. Knepley   }
45277623264SMatthew G. Knepley   PetscFunctionReturn(0);
45377623264SMatthew G. Knepley }
45477623264SMatthew G. Knepley #undef __FUNCT__
45577623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerSetTypeFromOptions_Internal"
45677623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetTypeFromOptions_Internal(PetscPartitioner part)
45777623264SMatthew G. Knepley {
45877623264SMatthew G. Knepley   const char    *defaultType;
45977623264SMatthew G. Knepley   char           name[256];
46077623264SMatthew G. Knepley   PetscBool      flg;
46177623264SMatthew G. Knepley   PetscErrorCode ierr;
46277623264SMatthew G. Knepley 
46377623264SMatthew G. Knepley   PetscFunctionBegin;
46477623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
46577623264SMatthew G. Knepley   if (!((PetscObject) part)->type_name) defaultType = PETSCPARTITIONERCHACO;
46677623264SMatthew G. Knepley   else                                  defaultType = ((PetscObject) part)->type_name;
46777623264SMatthew G. Knepley   if (!PetscPartitionerRegisterAllCalled) {ierr = PetscPartitionerRegisterAll();CHKERRQ(ierr);}
46877623264SMatthew G. Knepley 
46977623264SMatthew G. Knepley   ierr = PetscObjectOptionsBegin((PetscObject) part);CHKERRQ(ierr);
47077623264SMatthew G. Knepley   ierr = PetscOptionsFList("-petscpartitioner_type", "Graph partitioner", "PetscPartitionerSetType", PetscPartitionerList, defaultType, name, 256, &flg);CHKERRQ(ierr);
47177623264SMatthew G. Knepley   if (flg) {
47277623264SMatthew G. Knepley     ierr = PetscPartitionerSetType(part, name);CHKERRQ(ierr);
47377623264SMatthew G. Knepley   } else if (!((PetscObject) part)->type_name) {
47477623264SMatthew G. Knepley     ierr = PetscPartitionerSetType(part, defaultType);CHKERRQ(ierr);
47577623264SMatthew G. Knepley   }
47677623264SMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
47777623264SMatthew G. Knepley   PetscFunctionReturn(0);
47877623264SMatthew G. Knepley }
47977623264SMatthew G. Knepley 
48077623264SMatthew G. Knepley #undef __FUNCT__
48177623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerSetFromOptions"
48277623264SMatthew G. Knepley /*@
48377623264SMatthew G. Knepley   PetscPartitionerSetFromOptions - sets parameters in a PetscPartitioner from the options database
48477623264SMatthew G. Knepley 
48577623264SMatthew G. Knepley   Collective on PetscPartitioner
48677623264SMatthew G. Knepley 
48777623264SMatthew G. Knepley   Input Parameter:
48877623264SMatthew G. Knepley . part - the PetscPartitioner object to set options for
48977623264SMatthew G. Knepley 
49077623264SMatthew G. Knepley   Level: developer
49177623264SMatthew G. Knepley 
49277623264SMatthew G. Knepley .seealso: PetscPartitionerView()
49377623264SMatthew G. Knepley @*/
49477623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetFromOptions(PetscPartitioner part)
49577623264SMatthew G. Knepley {
49677623264SMatthew G. Knepley   PetscErrorCode ierr;
49777623264SMatthew G. Knepley 
49877623264SMatthew G. Knepley   PetscFunctionBegin;
49977623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
50077623264SMatthew G. Knepley   ierr = PetscPartitionerSetTypeFromOptions_Internal(part);CHKERRQ(ierr);
50177623264SMatthew G. Knepley 
50277623264SMatthew G. Knepley   ierr = PetscObjectOptionsBegin((PetscObject) part);CHKERRQ(ierr);
50377623264SMatthew G. Knepley   if (part->ops->setfromoptions) {ierr = (*part->ops->setfromoptions)(part);CHKERRQ(ierr);}
50477623264SMatthew G. Knepley   /* process any options handlers added with PetscObjectAddOptionsHandler() */
50577623264SMatthew G. Knepley   ierr = PetscObjectProcessOptionsHandlers((PetscObject) part);CHKERRQ(ierr);
50677623264SMatthew G. Knepley   ierr = PetscOptionsEnd();CHKERRQ(ierr);
50777623264SMatthew G. Knepley   ierr = PetscPartitionerViewFromOptions(part, NULL, "-petscpartitioner_view");CHKERRQ(ierr);
50877623264SMatthew G. Knepley   PetscFunctionReturn(0);
50977623264SMatthew G. Knepley }
51077623264SMatthew G. Knepley 
51177623264SMatthew G. Knepley #undef __FUNCT__
51277623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerSetUp"
51377623264SMatthew G. Knepley /*@C
51477623264SMatthew G. Knepley   PetscPartitionerSetUp - Construct data structures for the PetscPartitioner
51577623264SMatthew G. Knepley 
51677623264SMatthew G. Knepley   Collective on PetscPartitioner
51777623264SMatthew G. Knepley 
51877623264SMatthew G. Knepley   Input Parameter:
51977623264SMatthew G. Knepley . part - the PetscPartitioner object to setup
52077623264SMatthew G. Knepley 
52177623264SMatthew G. Knepley   Level: developer
52277623264SMatthew G. Knepley 
52377623264SMatthew G. Knepley .seealso: PetscPartitionerView(), PetscPartitionerDestroy()
52477623264SMatthew G. Knepley @*/
52577623264SMatthew G. Knepley PetscErrorCode PetscPartitionerSetUp(PetscPartitioner part)
52677623264SMatthew G. Knepley {
52777623264SMatthew G. Knepley   PetscErrorCode ierr;
52877623264SMatthew G. Knepley 
52977623264SMatthew G. Knepley   PetscFunctionBegin;
53077623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
53177623264SMatthew G. Knepley   if (part->ops->setup) {ierr = (*part->ops->setup)(part);CHKERRQ(ierr);}
53277623264SMatthew G. Knepley   PetscFunctionReturn(0);
53377623264SMatthew G. Knepley }
53477623264SMatthew G. Knepley 
53577623264SMatthew G. Knepley #undef __FUNCT__
53677623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerDestroy"
53777623264SMatthew G. Knepley /*@
53877623264SMatthew G. Knepley   PetscPartitionerDestroy - Destroys a PetscPartitioner object
53977623264SMatthew G. Knepley 
54077623264SMatthew G. Knepley   Collective on PetscPartitioner
54177623264SMatthew G. Knepley 
54277623264SMatthew G. Knepley   Input Parameter:
54377623264SMatthew G. Knepley . part - the PetscPartitioner object to destroy
54477623264SMatthew G. Knepley 
54577623264SMatthew G. Knepley   Level: developer
54677623264SMatthew G. Knepley 
54777623264SMatthew G. Knepley .seealso: PetscPartitionerView()
54877623264SMatthew G. Knepley @*/
54977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy(PetscPartitioner *part)
55077623264SMatthew G. Knepley {
55177623264SMatthew G. Knepley   PetscErrorCode ierr;
55277623264SMatthew G. Knepley 
55377623264SMatthew G. Knepley   PetscFunctionBegin;
55477623264SMatthew G. Knepley   if (!*part) PetscFunctionReturn(0);
55577623264SMatthew G. Knepley   PetscValidHeaderSpecific((*part), PETSCPARTITIONER_CLASSID, 1);
55677623264SMatthew G. Knepley 
55777623264SMatthew G. Knepley   if (--((PetscObject)(*part))->refct > 0) {*part = 0; PetscFunctionReturn(0);}
55877623264SMatthew G. Knepley   ((PetscObject) (*part))->refct = 0;
55977623264SMatthew G. Knepley 
56077623264SMatthew G. Knepley   if ((*part)->ops->destroy) {ierr = (*(*part)->ops->destroy)(*part);CHKERRQ(ierr);}
56177623264SMatthew G. Knepley   ierr = PetscHeaderDestroy(part);CHKERRQ(ierr);
56277623264SMatthew G. Knepley   PetscFunctionReturn(0);
56377623264SMatthew G. Knepley }
56477623264SMatthew G. Knepley 
56577623264SMatthew G. Knepley #undef __FUNCT__
56677623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerCreate"
56777623264SMatthew G. Knepley /*@
56877623264SMatthew G. Knepley   PetscPartitionerCreate - Creates an empty PetscPartitioner object. The type can then be set with PetscPartitionerSetType().
56977623264SMatthew G. Knepley 
57077623264SMatthew G. Knepley   Collective on MPI_Comm
57177623264SMatthew G. Knepley 
57277623264SMatthew G. Knepley   Input Parameter:
57377623264SMatthew G. Knepley . comm - The communicator for the PetscPartitioner object
57477623264SMatthew G. Knepley 
57577623264SMatthew G. Knepley   Output Parameter:
57677623264SMatthew G. Knepley . part - The PetscPartitioner object
57777623264SMatthew G. Knepley 
57877623264SMatthew G. Knepley   Level: beginner
57977623264SMatthew G. Knepley 
58077623264SMatthew G. Knepley .seealso: PetscPartitionerSetType(), PETSCPARTITIONERCHACO, PETSCPARTITIONERPARMETIS, PETSCPARTITIONERSHELL
58177623264SMatthew G. Knepley @*/
58277623264SMatthew G. Knepley PetscErrorCode PetscPartitionerCreate(MPI_Comm comm, PetscPartitioner *part)
58377623264SMatthew G. Knepley {
58477623264SMatthew G. Knepley   PetscPartitioner p;
58577623264SMatthew G. Knepley   PetscErrorCode   ierr;
58677623264SMatthew G. Knepley 
58777623264SMatthew G. Knepley   PetscFunctionBegin;
58877623264SMatthew G. Knepley   PetscValidPointer(part, 2);
58977623264SMatthew G. Knepley   *part = NULL;
59077623264SMatthew G. Knepley   ierr = PetscFVInitializePackage();CHKERRQ(ierr);
59177623264SMatthew G. Knepley 
59277623264SMatthew G. Knepley   ierr = PetscHeaderCreate(p, _p_PetscPartitioner, struct _PetscPartitionerOps, PETSCPARTITIONER_CLASSID, "PetscPartitioner", "Graph Partitioner", "PetscPartitioner", comm, PetscPartitionerDestroy, PetscPartitionerView);CHKERRQ(ierr);
59377623264SMatthew G. Knepley   ierr = PetscMemzero(p->ops, sizeof(struct _PetscPartitionerOps));CHKERRQ(ierr);
59477623264SMatthew G. Knepley 
59577623264SMatthew G. Knepley   *part = p;
59677623264SMatthew G. Knepley   PetscFunctionReturn(0);
59777623264SMatthew G. Knepley }
59877623264SMatthew G. Knepley 
59977623264SMatthew G. Knepley #undef __FUNCT__
60077623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerPartition"
60177623264SMatthew G. Knepley /*@
60277623264SMatthew G. Knepley   PetscPartitionerPartition - Create a non-overlapping partition of the cells in the mesh
60377623264SMatthew G. Knepley 
60477623264SMatthew G. Knepley   Collective on DM
60577623264SMatthew G. Knepley 
60677623264SMatthew G. Knepley   Input Parameters:
60777623264SMatthew G. Knepley + part    - The PetscPartitioner
608*f8987ae8SMichael Lange - dm      - The mesh DM
60977623264SMatthew G. Knepley 
61077623264SMatthew G. Knepley   Output Parameters:
61177623264SMatthew G. Knepley + partSection     - The PetscSection giving the division of points by partition
612*f8987ae8SMichael Lange - partition       - The list of points by partition
61377623264SMatthew G. Knepley 
61477623264SMatthew G. Knepley   Note: Instead of cells, points at a given height can be partitioned by calling PetscPartitionerSetPointHeight()
61577623264SMatthew G. Knepley 
61677623264SMatthew G. Knepley   Level: developer
61777623264SMatthew G. Knepley 
61877623264SMatthew G. Knepley .seealso DMPlexDistribute(), PetscPartitionerSetPointHeight(), PetscPartitionerCreate()
6194b15ede2SMatthew G. Knepley @*/
620*f8987ae8SMichael Lange PetscErrorCode PetscPartitionerPartition(PetscPartitioner part, DM dm, PetscSection partSection, IS *partition)
62177623264SMatthew G. Knepley {
62277623264SMatthew G. Knepley   PetscMPIInt    size;
62377623264SMatthew G. Knepley   PetscErrorCode ierr;
62477623264SMatthew G. Knepley 
62577623264SMatthew G. Knepley   PetscFunctionBegin;
62677623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
62777623264SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
62877623264SMatthew G. Knepley   PetscValidHeaderSpecific(partSection, PETSC_SECTION_CLASSID, 4);
62977623264SMatthew G. Knepley   PetscValidPointer(partition, 5);
63077623264SMatthew G. Knepley   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) part), &size);CHKERRQ(ierr);
63177623264SMatthew G. Knepley   if (size == 1) {
63277623264SMatthew G. Knepley     PetscInt *points;
63377623264SMatthew G. Knepley     PetscInt  cStart, cEnd, c;
63477623264SMatthew G. Knepley 
63577623264SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, part->height, &cStart, &cEnd);CHKERRQ(ierr);
63677623264SMatthew G. Knepley     ierr = PetscSectionSetChart(partSection, 0, size);CHKERRQ(ierr);
63777623264SMatthew G. Knepley     ierr = PetscSectionSetDof(partSection, 0, cEnd-cStart);CHKERRQ(ierr);
63877623264SMatthew G. Knepley     ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
63977623264SMatthew G. Knepley     ierr = PetscMalloc1(cEnd-cStart, &points);CHKERRQ(ierr);
64077623264SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) points[c] = c;
64177623264SMatthew G. Knepley     ierr = ISCreateGeneral(PetscObjectComm((PetscObject) part), cEnd-cStart, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
64277623264SMatthew G. Knepley     PetscFunctionReturn(0);
64377623264SMatthew G. Knepley   }
64477623264SMatthew G. Knepley   if (part->height == 0) {
64577623264SMatthew G. Knepley     PetscInt  numVertices;
64677623264SMatthew G. Knepley     PetscInt *start     = NULL;
64777623264SMatthew G. Knepley     PetscInt *adjacency = NULL;
64877623264SMatthew G. Knepley 
64977623264SMatthew G. Knepley     ierr = DMPlexCreateNeighborCSR(dm, 0, &numVertices, &start, &adjacency);CHKERRQ(ierr);
65077623264SMatthew G. Knepley     if (!part->ops->partition) SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_WRONGSTATE, "PetscPartitioner has no type");
65177623264SMatthew G. Knepley     ierr = (*part->ops->partition)(part, dm, size, numVertices, start, adjacency, partSection, partition);CHKERRQ(ierr);
65277623264SMatthew G. Knepley     ierr = PetscFree(start);CHKERRQ(ierr);
65377623264SMatthew G. Knepley     ierr = PetscFree(adjacency);CHKERRQ(ierr);
65477623264SMatthew G. Knepley   } else SETERRQ1(PetscObjectComm((PetscObject) part), PETSC_ERR_ARG_OUTOFRANGE, "Invalid height %D for points to partition", part->height);
65577623264SMatthew G. Knepley   PetscFunctionReturn(0);
65677623264SMatthew G. Knepley }
65777623264SMatthew G. Knepley 
65877623264SMatthew G. Knepley #undef __FUNCT__
65977623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerDestroy_Shell"
66077623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy_Shell(PetscPartitioner part)
66177623264SMatthew G. Knepley {
66277623264SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
66377623264SMatthew G. Knepley   PetscErrorCode          ierr;
66477623264SMatthew G. Knepley 
66577623264SMatthew G. Knepley   PetscFunctionBegin;
66677623264SMatthew G. Knepley   ierr = PetscSectionDestroy(&p->section);CHKERRQ(ierr);
66777623264SMatthew G. Knepley   ierr = ISDestroy(&p->partition);CHKERRQ(ierr);
66877623264SMatthew G. Knepley   ierr = PetscFree(p);CHKERRQ(ierr);
66977623264SMatthew G. Knepley   PetscFunctionReturn(0);
67077623264SMatthew G. Knepley }
67177623264SMatthew G. Knepley 
67277623264SMatthew G. Knepley #undef __FUNCT__
67377623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_Shell_Ascii"
67477623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Shell_Ascii(PetscPartitioner part, PetscViewer viewer)
67577623264SMatthew G. Knepley {
67677623264SMatthew G. Knepley   PetscViewerFormat format;
67777623264SMatthew G. Knepley   PetscErrorCode    ierr;
67877623264SMatthew G. Knepley 
67977623264SMatthew G. Knepley   PetscFunctionBegin;
68077623264SMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
68177623264SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Shell Graph Partitioner:\n");CHKERRQ(ierr);
68277623264SMatthew G. Knepley   PetscFunctionReturn(0);
68377623264SMatthew G. Knepley }
68477623264SMatthew G. Knepley 
68577623264SMatthew G. Knepley #undef __FUNCT__
68677623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_Shell"
68777623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Shell(PetscPartitioner part, PetscViewer viewer)
68877623264SMatthew G. Knepley {
68977623264SMatthew G. Knepley   PetscBool      iascii;
69077623264SMatthew G. Knepley   PetscErrorCode ierr;
69177623264SMatthew G. Knepley 
69277623264SMatthew G. Knepley   PetscFunctionBegin;
69377623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
69477623264SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
69577623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
69677623264SMatthew G. Knepley   if (iascii) {ierr = PetscPartitionerView_Shell_Ascii(part, viewer);CHKERRQ(ierr);}
69777623264SMatthew G. Knepley   PetscFunctionReturn(0);
69877623264SMatthew G. Knepley }
69977623264SMatthew G. Knepley 
70077623264SMatthew G. Knepley #undef __FUNCT__
70177623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerPartition_Shell"
70277623264SMatthew G. Knepley PetscErrorCode PetscPartitionerPartition_Shell(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
70377623264SMatthew G. Knepley {
70477623264SMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
70577623264SMatthew G. Knepley   PetscInt                np;
70677623264SMatthew G. Knepley   PetscErrorCode          ierr;
70777623264SMatthew G. Knepley 
70877623264SMatthew G. Knepley   PetscFunctionBegin;
70977623264SMatthew G. Knepley   ierr = PetscSectionGetChart(p->section, NULL, &np);CHKERRQ(ierr);
71077623264SMatthew G. Knepley   if (nparts != np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of requested partitions %d != configured partitions %d", nparts, np);
71177623264SMatthew G. Knepley   ierr = ISGetLocalSize(p->partition, &np);CHKERRQ(ierr);
71277623264SMatthew G. Knepley   if (numVertices != np) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of input vertices %d != configured vertices %d", numVertices, np);
7135680f57bSMatthew G. Knepley   ierr = PetscSectionCopy(p->section, partSection);CHKERRQ(ierr);
71477623264SMatthew G. Knepley   *partition = p->partition;
71577623264SMatthew G. Knepley   ierr = PetscObjectReference((PetscObject) p->partition);CHKERRQ(ierr);
71677623264SMatthew G. Knepley   PetscFunctionReturn(0);
71777623264SMatthew G. Knepley }
71877623264SMatthew G. Knepley 
71977623264SMatthew G. Knepley #undef __FUNCT__
72077623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerInitialize_Shell"
72177623264SMatthew G. Knepley PetscErrorCode PetscPartitionerInitialize_Shell(PetscPartitioner part)
72277623264SMatthew G. Knepley {
72377623264SMatthew G. Knepley   PetscFunctionBegin;
72477623264SMatthew G. Knepley   part->ops->view      = PetscPartitionerView_Shell;
72577623264SMatthew G. Knepley   part->ops->destroy   = PetscPartitionerDestroy_Shell;
72677623264SMatthew G. Knepley   part->ops->partition = PetscPartitionerPartition_Shell;
72777623264SMatthew G. Knepley   PetscFunctionReturn(0);
72877623264SMatthew G. Knepley }
72977623264SMatthew G. Knepley 
73077623264SMatthew G. Knepley /*MC
73177623264SMatthew G. Knepley   PETSCPARTITIONERSHELL = "shell" - A PetscPartitioner object
73277623264SMatthew G. Knepley 
73377623264SMatthew G. Knepley   Level: intermediate
73477623264SMatthew G. Knepley 
73577623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
73677623264SMatthew G. Knepley M*/
73777623264SMatthew G. Knepley 
73877623264SMatthew G. Knepley #undef __FUNCT__
73977623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerCreate_Shell"
74077623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Shell(PetscPartitioner part)
74177623264SMatthew G. Knepley {
74277623264SMatthew G. Knepley   PetscPartitioner_Shell *p;
74377623264SMatthew G. Knepley   PetscErrorCode          ierr;
74477623264SMatthew G. Knepley 
74577623264SMatthew G. Knepley   PetscFunctionBegin;
74677623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
74777623264SMatthew G. Knepley   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
74877623264SMatthew G. Knepley   part->data = p;
74977623264SMatthew G. Knepley 
75077623264SMatthew G. Knepley   ierr = PetscPartitionerInitialize_Shell(part);CHKERRQ(ierr);
75177623264SMatthew G. Knepley   PetscFunctionReturn(0);
75277623264SMatthew G. Knepley }
75377623264SMatthew G. Knepley 
75477623264SMatthew G. Knepley #undef __FUNCT__
7555680f57bSMatthew G. Knepley #define __FUNCT__ "PetscPartitionerShellSetPartition"
7565680f57bSMatthew G. Knepley /*@C
7575680f57bSMatthew G. Knepley   PetscPartitionerShellSetPartition - Set an artifical partition for a mesh
7585680f57bSMatthew G. Knepley 
7595680f57bSMatthew G. Knepley   Collective on DM
7605680f57bSMatthew G. Knepley 
7615680f57bSMatthew G. Knepley   Input Parameters:
7625680f57bSMatthew G. Knepley + part    - The PetscPartitioner
7635680f57bSMatthew G. Knepley . dm      - The mesh DM
7645680f57bSMatthew G. Knepley - enlarge - Expand each partition with neighbors
7655680f57bSMatthew G. Knepley 
7665680f57bSMatthew G. Knepley   Level: developer
7675680f57bSMatthew G. Knepley 
7685680f57bSMatthew G. Knepley .seealso DMPlexDistribute(), PetscPartitionerCreate()
7695680f57bSMatthew G. Knepley @*/
7705680f57bSMatthew G. Knepley PetscErrorCode PetscPartitionerShellSetPartition(PetscPartitioner part, PetscInt numProcs, const PetscInt sizes[], const PetscInt points[])
7715680f57bSMatthew G. Knepley {
7725680f57bSMatthew G. Knepley   PetscPartitioner_Shell *p = (PetscPartitioner_Shell *) part->data;
7735680f57bSMatthew G. Knepley   PetscInt                proc, numPoints;
7745680f57bSMatthew G. Knepley   PetscErrorCode          ierr;
7755680f57bSMatthew G. Knepley 
7765680f57bSMatthew G. Knepley   PetscFunctionBegin;
7775680f57bSMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
7785680f57bSMatthew G. Knepley   if (sizes) {PetscValidPointer(sizes, 3);}
7795680f57bSMatthew G. Knepley   if (sizes) {PetscValidPointer(points, 4);}
7805680f57bSMatthew G. Knepley   ierr = PetscSectionDestroy(&p->section);CHKERRQ(ierr);
7815680f57bSMatthew G. Knepley   ierr = ISDestroy(&p->partition);CHKERRQ(ierr);
7825680f57bSMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject) part), &p->section);CHKERRQ(ierr);
7835680f57bSMatthew G. Knepley   ierr = PetscSectionSetChart(p->section, 0, numProcs);CHKERRQ(ierr);
7845680f57bSMatthew G. Knepley   if (sizes) {
7855680f57bSMatthew G. Knepley     for (proc = 0; proc < numProcs; ++proc) {
7865680f57bSMatthew G. Knepley       ierr = PetscSectionSetDof(p->section, proc, sizes[proc]);CHKERRQ(ierr);
7875680f57bSMatthew G. Knepley     }
7885680f57bSMatthew G. Knepley   }
7895680f57bSMatthew G. Knepley   ierr = PetscSectionSetUp(p->section);CHKERRQ(ierr);
7905680f57bSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(p->section, &numPoints);CHKERRQ(ierr);
7915680f57bSMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) part), numPoints, points, PETSC_COPY_VALUES, &p->partition);CHKERRQ(ierr);
7925680f57bSMatthew G. Knepley   PetscFunctionReturn(0);
7935680f57bSMatthew G. Knepley }
7945680f57bSMatthew G. Knepley 
7955680f57bSMatthew G. Knepley #undef __FUNCT__
79677623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerDestroy_Chaco"
79777623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy_Chaco(PetscPartitioner part)
79877623264SMatthew G. Knepley {
79977623264SMatthew G. Knepley   PetscPartitioner_Chaco *p = (PetscPartitioner_Chaco *) part->data;
80077623264SMatthew G. Knepley   PetscErrorCode          ierr;
80177623264SMatthew G. Knepley 
80277623264SMatthew G. Knepley   PetscFunctionBegin;
80377623264SMatthew G. Knepley   ierr = PetscFree(p);CHKERRQ(ierr);
80477623264SMatthew G. Knepley   PetscFunctionReturn(0);
80577623264SMatthew G. Knepley }
80677623264SMatthew G. Knepley 
80777623264SMatthew G. Knepley #undef __FUNCT__
80877623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_Chaco_Ascii"
80977623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Chaco_Ascii(PetscPartitioner part, PetscViewer viewer)
81077623264SMatthew G. Knepley {
81177623264SMatthew G. Knepley   PetscViewerFormat format;
81277623264SMatthew G. Knepley   PetscErrorCode    ierr;
81377623264SMatthew G. Knepley 
81477623264SMatthew G. Knepley   PetscFunctionBegin;
81577623264SMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
81677623264SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "Chaco Graph Partitioner:\n");CHKERRQ(ierr);
81777623264SMatthew G. Knepley   PetscFunctionReturn(0);
81877623264SMatthew G. Knepley }
81977623264SMatthew G. Knepley 
82077623264SMatthew G. Knepley #undef __FUNCT__
82177623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_Chaco"
82277623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_Chaco(PetscPartitioner part, PetscViewer viewer)
82377623264SMatthew G. Knepley {
82477623264SMatthew G. Knepley   PetscBool      iascii;
82577623264SMatthew G. Knepley   PetscErrorCode ierr;
82677623264SMatthew G. Knepley 
82777623264SMatthew G. Knepley   PetscFunctionBegin;
82877623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
82977623264SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
83077623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
83177623264SMatthew G. Knepley   if (iascii) {ierr = PetscPartitionerView_Chaco_Ascii(part, viewer);CHKERRQ(ierr);}
83277623264SMatthew G. Knepley   PetscFunctionReturn(0);
83377623264SMatthew G. Knepley }
83477623264SMatthew G. Knepley 
83570034214SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO)
83670034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
83770034214SMatthew G. Knepley #include <unistd.h>
83870034214SMatthew G. Knepley #endif
83970034214SMatthew G. Knepley /* Chaco does not have an include file */
84070034214SMatthew G. Knepley PETSC_EXTERN int interface(int nvtxs, int *start, int *adjacency, int *vwgts,
84170034214SMatthew G. Knepley                        float *ewgts, float *x, float *y, float *z, char *outassignname,
84270034214SMatthew G. Knepley                        char *outfilename, short *assignment, int architecture, int ndims_tot,
84370034214SMatthew G. Knepley                        int mesh_dims[3], double *goal, int global_method, int local_method,
84470034214SMatthew G. Knepley                        int rqi_flag, int vmax, int ndims, double eigtol, long seed);
84570034214SMatthew G. Knepley 
84670034214SMatthew G. Knepley extern int FREE_GRAPH;
84777623264SMatthew G. Knepley #endif
84870034214SMatthew G. Knepley 
84970034214SMatthew G. Knepley #undef __FUNCT__
85077623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerPartition_Chaco"
85177623264SMatthew G. Knepley PetscErrorCode PetscPartitionerPartition_Chaco(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
85270034214SMatthew G. Knepley {
85377623264SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO)
85470034214SMatthew G. Knepley   enum {DEFAULT_METHOD = 1, INERTIAL_METHOD = 3};
85570034214SMatthew G. Knepley   MPI_Comm       comm;
85670034214SMatthew G. Knepley   int            nvtxs          = numVertices; /* number of vertices in full graph */
85770034214SMatthew G. Knepley   int           *vwgts          = NULL;   /* weights for all vertices */
85870034214SMatthew G. Knepley   float         *ewgts          = NULL;   /* weights for all edges */
85970034214SMatthew G. Knepley   float         *x              = NULL, *y = NULL, *z = NULL; /* coordinates for inertial method */
86070034214SMatthew G. Knepley   char          *outassignname  = NULL;   /*  name of assignment output file */
86170034214SMatthew G. Knepley   char          *outfilename    = NULL;   /* output file name */
86270034214SMatthew G. Knepley   int            architecture   = 1;      /* 0 => hypercube, d => d-dimensional mesh */
86370034214SMatthew G. Knepley   int            ndims_tot      = 0;      /* total number of cube dimensions to divide */
86470034214SMatthew G. Knepley   int            mesh_dims[3];            /* dimensions of mesh of processors */
86570034214SMatthew G. Knepley   double        *goal          = NULL;    /* desired set sizes for each set */
86670034214SMatthew G. Knepley   int            global_method = 1;       /* global partitioning algorithm */
86770034214SMatthew G. Knepley   int            local_method  = 1;       /* local partitioning algorithm */
86870034214SMatthew G. Knepley   int            rqi_flag      = 0;       /* should I use RQI/Symmlq eigensolver? */
86970034214SMatthew G. Knepley   int            vmax          = 200;     /* how many vertices to coarsen down to? */
87070034214SMatthew G. Knepley   int            ndims         = 1;       /* number of eigenvectors (2^d sets) */
87170034214SMatthew G. Knepley   double         eigtol        = 0.001;   /* tolerance on eigenvectors */
87270034214SMatthew G. Knepley   long           seed          = 123636512; /* for random graph mutations */
87370034214SMatthew G. Knepley   short int     *assignment;              /* Output partition */
87470034214SMatthew G. Knepley   int            fd_stdout, fd_pipe[2];
87570034214SMatthew G. Knepley   PetscInt      *points;
87670034214SMatthew G. Knepley   int            i, v, p;
87770034214SMatthew G. Knepley   PetscErrorCode ierr;
87870034214SMatthew G. Knepley 
87970034214SMatthew G. Knepley   PetscFunctionBegin;
88070034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
88170034214SMatthew G. Knepley   if (!numVertices) {
88277623264SMatthew G. Knepley     ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
88377623264SMatthew G. Knepley     ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
88470034214SMatthew G. Knepley     ierr = ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
88570034214SMatthew G. Knepley     PetscFunctionReturn(0);
88670034214SMatthew G. Knepley   }
88770034214SMatthew G. Knepley   FREE_GRAPH = 0;                         /* Do not let Chaco free my memory */
88870034214SMatthew G. Knepley   for (i = 0; i < start[numVertices]; ++i) ++adjacency[i];
88970034214SMatthew G. Knepley 
89070034214SMatthew G. Knepley   if (global_method == INERTIAL_METHOD) {
89170034214SMatthew G. Knepley     /* manager.createCellCoordinates(nvtxs, &x, &y, &z); */
89270034214SMatthew G. Knepley     SETERRQ(comm, PETSC_ERR_SUP, "Inertial partitioning not yet supported");
89370034214SMatthew G. Knepley   }
89477623264SMatthew G. Knepley   mesh_dims[0] = nparts;
89570034214SMatthew G. Knepley   mesh_dims[1] = 1;
89670034214SMatthew G. Knepley   mesh_dims[2] = 1;
89770034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &assignment);CHKERRQ(ierr);
89870034214SMatthew G. Knepley   /* Chaco outputs to stdout. We redirect this to a buffer. */
89970034214SMatthew G. Knepley   /* TODO: check error codes for UNIX calls */
90070034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
90170034214SMatthew G. Knepley   {
90270034214SMatthew G. Knepley     int piperet;
90370034214SMatthew G. Knepley     piperet = pipe(fd_pipe);
90470034214SMatthew G. Knepley     if (piperet) SETERRQ(comm,PETSC_ERR_SYS,"Could not create pipe");
90570034214SMatthew G. Knepley     fd_stdout = dup(1);
90670034214SMatthew G. Knepley     close(1);
90770034214SMatthew G. Knepley     dup2(fd_pipe[1], 1);
90870034214SMatthew G. Knepley   }
90970034214SMatthew G. Knepley #endif
91070034214SMatthew G. Knepley   ierr = interface(nvtxs, (int*) start, (int*) adjacency, vwgts, ewgts, x, y, z, outassignname, outfilename,
91170034214SMatthew G. Knepley                    assignment, architecture, ndims_tot, mesh_dims, goal, global_method, local_method, rqi_flag,
91270034214SMatthew G. Knepley                    vmax, ndims, eigtol, seed);
91370034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
91470034214SMatthew G. Knepley   {
91570034214SMatthew G. Knepley     char msgLog[10000];
91670034214SMatthew G. Knepley     int  count;
91770034214SMatthew G. Knepley 
91870034214SMatthew G. Knepley     fflush(stdout);
91970034214SMatthew G. Knepley     count = read(fd_pipe[0], msgLog, (10000-1)*sizeof(char));
92070034214SMatthew G. Knepley     if (count < 0) count = 0;
92170034214SMatthew G. Knepley     msgLog[count] = 0;
92270034214SMatthew G. Knepley     close(1);
92370034214SMatthew G. Knepley     dup2(fd_stdout, 1);
92470034214SMatthew G. Knepley     close(fd_stdout);
92570034214SMatthew G. Knepley     close(fd_pipe[0]);
92670034214SMatthew G. Knepley     close(fd_pipe[1]);
92770034214SMatthew G. Knepley     if (ierr) SETERRQ1(comm, PETSC_ERR_LIB, "Error in Chaco library: %s", msgLog);
92870034214SMatthew G. Knepley   }
92970034214SMatthew G. Knepley #endif
93070034214SMatthew G. Knepley   /* Convert to PetscSection+IS */
93177623264SMatthew G. Knepley   ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
93270034214SMatthew G. Knepley   for (v = 0; v < nvtxs; ++v) {
93377623264SMatthew G. Knepley     ierr = PetscSectionAddDof(partSection, assignment[v], 1);CHKERRQ(ierr);
93470034214SMatthew G. Knepley   }
93577623264SMatthew G. Knepley   ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
93670034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr);
93777623264SMatthew G. Knepley   for (p = 0, i = 0; p < nparts; ++p) {
93870034214SMatthew G. Knepley     for (v = 0; v < nvtxs; ++v) {
93970034214SMatthew G. Knepley       if (assignment[v] == p) points[i++] = v;
94070034214SMatthew G. Knepley     }
94170034214SMatthew G. Knepley   }
94270034214SMatthew G. Knepley   if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs);
94370034214SMatthew G. Knepley   ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
94470034214SMatthew G. Knepley   if (global_method == INERTIAL_METHOD) {
94570034214SMatthew G. Knepley     /* manager.destroyCellCoordinates(nvtxs, &x, &y, &z); */
94670034214SMatthew G. Knepley   }
94770034214SMatthew G. Knepley   ierr = PetscFree(assignment);CHKERRQ(ierr);
94870034214SMatthew G. Knepley   for (i = 0; i < start[numVertices]; ++i) --adjacency[i];
94970034214SMatthew G. Knepley   PetscFunctionReturn(0);
95077623264SMatthew G. Knepley #else
95177623264SMatthew G. Knepley   SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-chaco.");
95270034214SMatthew G. Knepley #endif
95377623264SMatthew G. Knepley }
95477623264SMatthew G. Knepley 
95577623264SMatthew G. Knepley #undef __FUNCT__
95677623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerInitialize_Chaco"
95777623264SMatthew G. Knepley PetscErrorCode PetscPartitionerInitialize_Chaco(PetscPartitioner part)
95877623264SMatthew G. Knepley {
95977623264SMatthew G. Knepley   PetscFunctionBegin;
96077623264SMatthew G. Knepley   part->ops->view      = PetscPartitionerView_Chaco;
96177623264SMatthew G. Knepley   part->ops->destroy   = PetscPartitionerDestroy_Chaco;
96277623264SMatthew G. Knepley   part->ops->partition = PetscPartitionerPartition_Chaco;
96377623264SMatthew G. Knepley   PetscFunctionReturn(0);
96477623264SMatthew G. Knepley }
96577623264SMatthew G. Knepley 
96677623264SMatthew G. Knepley /*MC
96777623264SMatthew G. Knepley   PETSCPARTITIONERCHACO = "chaco" - A PetscPartitioner object using the Chaco library
96877623264SMatthew G. Knepley 
96977623264SMatthew G. Knepley   Level: intermediate
97077623264SMatthew G. Knepley 
97177623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
97277623264SMatthew G. Knepley M*/
97377623264SMatthew G. Knepley 
97477623264SMatthew G. Knepley #undef __FUNCT__
97577623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerCreate_Chaco"
97677623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_Chaco(PetscPartitioner part)
97777623264SMatthew G. Knepley {
97877623264SMatthew G. Knepley   PetscPartitioner_Chaco *p;
97977623264SMatthew G. Knepley   PetscErrorCode          ierr;
98077623264SMatthew G. Knepley 
98177623264SMatthew G. Knepley   PetscFunctionBegin;
98277623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
98377623264SMatthew G. Knepley   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
98477623264SMatthew G. Knepley   part->data = p;
98577623264SMatthew G. Knepley 
98677623264SMatthew G. Knepley   ierr = PetscPartitionerInitialize_Chaco(part);CHKERRQ(ierr);
98777623264SMatthew G. Knepley   ierr = PetscCitationsRegister(ChacoPartitionerCitation, &ChacoPartitionercite);CHKERRQ(ierr);
98877623264SMatthew G. Knepley   PetscFunctionReturn(0);
98977623264SMatthew G. Knepley }
99077623264SMatthew G. Knepley 
99177623264SMatthew G. Knepley #undef __FUNCT__
99277623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerDestroy_ParMetis"
99377623264SMatthew G. Knepley PetscErrorCode PetscPartitionerDestroy_ParMetis(PetscPartitioner part)
99477623264SMatthew G. Knepley {
99577623264SMatthew G. Knepley   PetscPartitioner_ParMetis *p = (PetscPartitioner_ParMetis *) part->data;
99677623264SMatthew G. Knepley   PetscErrorCode             ierr;
99777623264SMatthew G. Knepley 
99877623264SMatthew G. Knepley   PetscFunctionBegin;
99977623264SMatthew G. Knepley   ierr = PetscFree(p);CHKERRQ(ierr);
100077623264SMatthew G. Knepley   PetscFunctionReturn(0);
100177623264SMatthew G. Knepley }
100277623264SMatthew G. Knepley 
100377623264SMatthew G. Knepley #undef __FUNCT__
100477623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_ParMetis_Ascii"
100577623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_ParMetis_Ascii(PetscPartitioner part, PetscViewer viewer)
100677623264SMatthew G. Knepley {
100777623264SMatthew G. Knepley   PetscViewerFormat format;
100877623264SMatthew G. Knepley   PetscErrorCode    ierr;
100977623264SMatthew G. Knepley 
101077623264SMatthew G. Knepley   PetscFunctionBegin;
101177623264SMatthew G. Knepley   ierr = PetscViewerGetFormat(viewer, &format);CHKERRQ(ierr);
101277623264SMatthew G. Knepley   ierr = PetscViewerASCIIPrintf(viewer, "ParMetis Graph Partitioner:\n");CHKERRQ(ierr);
101377623264SMatthew G. Knepley   PetscFunctionReturn(0);
101477623264SMatthew G. Knepley }
101577623264SMatthew G. Knepley 
101677623264SMatthew G. Knepley #undef __FUNCT__
101777623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerView_ParMetis"
101877623264SMatthew G. Knepley PetscErrorCode PetscPartitionerView_ParMetis(PetscPartitioner part, PetscViewer viewer)
101977623264SMatthew G. Knepley {
102077623264SMatthew G. Knepley   PetscBool      iascii;
102177623264SMatthew G. Knepley   PetscErrorCode ierr;
102277623264SMatthew G. Knepley 
102377623264SMatthew G. Knepley   PetscFunctionBegin;
102477623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
102577623264SMatthew G. Knepley   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
102677623264SMatthew G. Knepley   ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
102777623264SMatthew G. Knepley   if (iascii) {ierr = PetscPartitionerView_ParMetis_Ascii(part, viewer);CHKERRQ(ierr);}
102877623264SMatthew G. Knepley   PetscFunctionReturn(0);
102977623264SMatthew G. Knepley }
103070034214SMatthew G. Knepley 
103170034214SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS)
103270034214SMatthew G. Knepley #include <parmetis.h>
103377623264SMatthew G. Knepley #endif
103470034214SMatthew G. Knepley 
103570034214SMatthew G. Knepley #undef __FUNCT__
103677623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerPartition_ParMetis"
103777623264SMatthew G. Knepley PetscErrorCode PetscPartitionerPartition_ParMetis(PetscPartitioner part, DM dm, PetscInt nparts, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection partSection, IS *partition)
103870034214SMatthew G. Knepley {
103977623264SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS)
104070034214SMatthew G. Knepley   MPI_Comm       comm;
104170034214SMatthew G. Knepley   PetscInt       nvtxs      = numVertices; /* The number of vertices in full graph */
104270034214SMatthew G. Knepley   PetscInt      *vtxdist;                  /* Distribution of vertices across processes */
104370034214SMatthew G. Knepley   PetscInt      *xadj       = start;       /* Start of edge list for each vertex */
104470034214SMatthew G. Knepley   PetscInt      *adjncy     = adjacency;   /* Edge lists for all vertices */
104570034214SMatthew G. Knepley   PetscInt      *vwgt       = NULL;        /* Vertex weights */
104670034214SMatthew G. Knepley   PetscInt      *adjwgt     = NULL;        /* Edge weights */
104770034214SMatthew G. Knepley   PetscInt       wgtflag    = 0;           /* Indicates which weights are present */
104870034214SMatthew G. Knepley   PetscInt       numflag    = 0;           /* Indicates initial offset (0 or 1) */
104970034214SMatthew G. Knepley   PetscInt       ncon       = 1;           /* The number of weights per vertex */
105070034214SMatthew G. Knepley   PetscReal     *tpwgts;                   /* The fraction of vertex weights assigned to each partition */
105170034214SMatthew G. Knepley   PetscReal     *ubvec;                    /* The balance intolerance for vertex weights */
105270034214SMatthew G. Knepley   PetscInt       options[5];               /* Options */
105370034214SMatthew G. Knepley   /* Outputs */
105470034214SMatthew G. Knepley   PetscInt       edgeCut;                  /* The number of edges cut by the partition */
105570034214SMatthew G. Knepley   PetscInt      *assignment, *points;
105677623264SMatthew G. Knepley   PetscMPIInt    rank, p, v, i;
105770034214SMatthew G. Knepley   PetscErrorCode ierr;
105870034214SMatthew G. Knepley 
105970034214SMatthew G. Knepley   PetscFunctionBegin;
106077623264SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) part, &comm);CHKERRQ(ierr);
106170034214SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
106270034214SMatthew G. Knepley   options[0] = 0; /* Use all defaults */
106370034214SMatthew G. Knepley   /* Calculate vertex distribution */
106470034214SMatthew G. Knepley   ierr = PetscMalloc4(nparts+1,&vtxdist,nparts*ncon,&tpwgts,ncon,&ubvec,nvtxs,&assignment);CHKERRQ(ierr);
106570034214SMatthew G. Knepley   vtxdist[0] = 0;
106670034214SMatthew G. Knepley   ierr = MPI_Allgather(&nvtxs, 1, MPIU_INT, &vtxdist[1], 1, MPIU_INT, comm);CHKERRQ(ierr);
106770034214SMatthew G. Knepley   for (p = 2; p <= nparts; ++p) {
106870034214SMatthew G. Knepley     vtxdist[p] += vtxdist[p-1];
106970034214SMatthew G. Knepley   }
107070034214SMatthew G. Knepley   /* Calculate weights */
107170034214SMatthew G. Knepley   for (p = 0; p < nparts; ++p) {
107270034214SMatthew G. Knepley     tpwgts[p] = 1.0/nparts;
107370034214SMatthew G. Knepley   }
107470034214SMatthew G. Knepley   ubvec[0] = 1.05;
107570034214SMatthew G. Knepley 
107670034214SMatthew G. Knepley   if (nparts == 1) {
107770034214SMatthew G. Knepley     ierr = PetscMemzero(assignment, nvtxs * sizeof(PetscInt));
107870034214SMatthew G. Knepley   } else {
107970034214SMatthew G. Knepley     if (vtxdist[1] == vtxdist[nparts]) {
108070034214SMatthew G. Knepley       if (!rank) {
108170034214SMatthew G. Knepley         PetscStackPush("METIS_PartGraphKway");
108270034214SMatthew G. Knepley         ierr = METIS_PartGraphKway(&nvtxs, &ncon, xadj, adjncy, vwgt, NULL, adjwgt, &nparts, tpwgts, ubvec, NULL, &edgeCut, assignment);
108370034214SMatthew G. Knepley         PetscStackPop;
108470034214SMatthew G. Knepley         if (ierr != METIS_OK) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in METIS_PartGraphKway()");
108570034214SMatthew G. Knepley       }
108670034214SMatthew G. Knepley     } else {
108770034214SMatthew G. Knepley       PetscStackPush("ParMETIS_V3_PartKway");
108870034214SMatthew G. Knepley       ierr = ParMETIS_V3_PartKway(vtxdist, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgeCut, assignment, &comm);
108970034214SMatthew G. Knepley       PetscStackPop;
109070034214SMatthew G. Knepley       if (ierr != METIS_OK) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in ParMETIS_V3_PartKway()");
109170034214SMatthew G. Knepley     }
109270034214SMatthew G. Knepley   }
109370034214SMatthew G. Knepley   /* Convert to PetscSection+IS */
109477623264SMatthew G. Knepley   ierr = PetscSectionSetChart(partSection, 0, nparts);CHKERRQ(ierr);
109577623264SMatthew G. Knepley   for (v = 0; v < nvtxs; ++v) {ierr = PetscSectionAddDof(partSection, assignment[v], 1);CHKERRQ(ierr);}
109677623264SMatthew G. Knepley   ierr = PetscSectionSetUp(partSection);CHKERRQ(ierr);
109770034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr);
109877623264SMatthew G. Knepley   for (p = 0, i = 0; p < nparts; ++p) {
109970034214SMatthew G. Knepley     for (v = 0; v < nvtxs; ++v) {
110070034214SMatthew G. Knepley       if (assignment[v] == p) points[i++] = v;
110170034214SMatthew G. Knepley     }
110270034214SMatthew G. Knepley   }
110370034214SMatthew G. Knepley   if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs);
110470034214SMatthew G. Knepley   ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
110570034214SMatthew G. Knepley   ierr = PetscFree4(vtxdist,tpwgts,ubvec,assignment);CHKERRQ(ierr);
110670034214SMatthew G. Knepley #else
110777623264SMatthew G. Knepley   SETERRQ(PetscObjectComm((PetscObject) part), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-parmetis.");
110870034214SMatthew G. Knepley #endif
110977623264SMatthew G. Knepley   PetscFunctionReturn(0);
111070034214SMatthew G. Knepley }
111170034214SMatthew G. Knepley 
111277623264SMatthew G. Knepley #undef __FUNCT__
111377623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerInitialize_ParMetis"
111477623264SMatthew G. Knepley PetscErrorCode PetscPartitionerInitialize_ParMetis(PetscPartitioner part)
111577623264SMatthew G. Knepley {
111677623264SMatthew G. Knepley   PetscFunctionBegin;
111777623264SMatthew G. Knepley   part->ops->view      = PetscPartitionerView_ParMetis;
111877623264SMatthew G. Knepley   part->ops->destroy   = PetscPartitionerDestroy_ParMetis;
111977623264SMatthew G. Knepley   part->ops->partition = PetscPartitionerPartition_ParMetis;
112077623264SMatthew G. Knepley   PetscFunctionReturn(0);
112177623264SMatthew G. Knepley }
112277623264SMatthew G. Knepley 
112377623264SMatthew G. Knepley /*MC
112477623264SMatthew G. Knepley   PETSCPARTITIONERPARMETIS = "parmetis" - A PetscPartitioner object using the ParMetis library
112577623264SMatthew G. Knepley 
112677623264SMatthew G. Knepley   Level: intermediate
112777623264SMatthew G. Knepley 
112877623264SMatthew G. Knepley .seealso: PetscPartitionerType, PetscPartitionerCreate(), PetscPartitionerSetType()
112977623264SMatthew G. Knepley M*/
113077623264SMatthew G. Knepley 
113177623264SMatthew G. Knepley #undef __FUNCT__
113277623264SMatthew G. Knepley #define __FUNCT__ "PetscPartitionerCreate_ParMetis"
113377623264SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscPartitionerCreate_ParMetis(PetscPartitioner part)
113477623264SMatthew G. Knepley {
113577623264SMatthew G. Knepley   PetscPartitioner_ParMetis *p;
113677623264SMatthew G. Knepley   PetscErrorCode          ierr;
113777623264SMatthew G. Knepley 
113877623264SMatthew G. Knepley   PetscFunctionBegin;
113977623264SMatthew G. Knepley   PetscValidHeaderSpecific(part, PETSCPARTITIONER_CLASSID, 1);
114077623264SMatthew G. Knepley   ierr       = PetscNewLog(part, &p);CHKERRQ(ierr);
114177623264SMatthew G. Knepley   part->data = p;
114277623264SMatthew G. Knepley 
114377623264SMatthew G. Knepley   ierr = PetscPartitionerInitialize_ParMetis(part);CHKERRQ(ierr);
114477623264SMatthew G. Knepley   ierr = PetscCitationsRegister(ParMetisPartitionerCitation, &ParMetisPartitionercite);CHKERRQ(ierr);
114570034214SMatthew G. Knepley   PetscFunctionReturn(0);
114670034214SMatthew G. Knepley }
114770034214SMatthew G. Knepley 
114870034214SMatthew G. Knepley #undef __FUNCT__
11495680f57bSMatthew G. Knepley #define __FUNCT__ "DMPlexGetPartitioner"
11505680f57bSMatthew G. Knepley /*@
11515680f57bSMatthew G. Knepley   DMPlexGetPartitioner - Get the mesh partitioner
11525680f57bSMatthew G. Knepley 
11535680f57bSMatthew G. Knepley   Not collective
11545680f57bSMatthew G. Knepley 
11555680f57bSMatthew G. Knepley   Input Parameter:
11565680f57bSMatthew G. Knepley . dm - The DM
11575680f57bSMatthew G. Knepley 
11585680f57bSMatthew G. Knepley   Output Parameter:
11595680f57bSMatthew G. Knepley . part - The PetscPartitioner
11605680f57bSMatthew G. Knepley 
11615680f57bSMatthew G. Knepley   Level: developer
11625680f57bSMatthew G. Knepley 
11635680f57bSMatthew G. Knepley .seealso DMPlexDistribute(), PetscPartitionerCreate()
11645680f57bSMatthew G. Knepley @*/
11655680f57bSMatthew G. Knepley PetscErrorCode DMPlexGetPartitioner(DM dm, PetscPartitioner *part)
11665680f57bSMatthew G. Knepley {
11675680f57bSMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
11685680f57bSMatthew G. Knepley 
11695680f57bSMatthew G. Knepley   PetscFunctionBegin;
11705680f57bSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11715680f57bSMatthew G. Knepley   PetscValidPointer(part, 2);
11725680f57bSMatthew G. Knepley   *part = mesh->partitioner;
11735680f57bSMatthew G. Knepley   PetscFunctionReturn(0);
11745680f57bSMatthew G. Knepley }
11755680f57bSMatthew G. Knepley 
11765680f57bSMatthew G. Knepley #undef __FUNCT__
11779ffc88e4SToby Isaac #define __FUNCT__ "DMPlexMarkTreeClosure"
11789ffc88e4SToby Isaac static PetscErrorCode DMPlexMarkTreeClosure(DM dm, PetscSegBuffer segpart, PetscBT bt, PetscInt point, PetscInt *partSize)
11799ffc88e4SToby Isaac {
11809ffc88e4SToby Isaac   PetscInt       parent, closureSize, *closure = NULL, i, pStart, pEnd;
11819ffc88e4SToby Isaac   PetscErrorCode ierr;
11829ffc88e4SToby Isaac 
11839ffc88e4SToby Isaac   PetscFunctionBegin;
11849ffc88e4SToby Isaac   ierr = DMPlexGetTreeParent(dm,point,&parent,NULL);CHKERRQ(ierr);
11859ffc88e4SToby Isaac   if (parent == point) PetscFunctionReturn(0);
11869ffc88e4SToby Isaac   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
11879ffc88e4SToby Isaac   ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
11889ffc88e4SToby Isaac   for (i = 0; i < closureSize; i++) {
11899ffc88e4SToby Isaac     PetscInt cpoint = closure[2*i];
11909ffc88e4SToby Isaac 
11919ffc88e4SToby Isaac     if (!PetscBTLookupSet(bt,cpoint-pStart)) {
11929ffc88e4SToby Isaac       PetscInt *PETSC_RESTRICT pt;
11939ffc88e4SToby Isaac       (*partSize)++;
11949ffc88e4SToby Isaac       ierr = PetscSegBufferGetInts(segpart,1,&pt);CHKERRQ(ierr);
11959ffc88e4SToby Isaac       *pt = cpoint;
11969ffc88e4SToby Isaac     }
11979ffc88e4SToby Isaac     ierr = DMPlexMarkTreeClosure(dm,segpart,bt,cpoint,partSize);CHKERRQ(ierr);
11989ffc88e4SToby Isaac   }
11999ffc88e4SToby Isaac   ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
12009ffc88e4SToby Isaac   PetscFunctionReturn(0);
12019ffc88e4SToby Isaac }
12029ffc88e4SToby Isaac 
12039ffc88e4SToby Isaac #undef __FUNCT__
120470034214SMatthew G. Knepley #define __FUNCT__ "DMPlexCreatePartitionClosure"
120570034214SMatthew G. Knepley PetscErrorCode DMPlexCreatePartitionClosure(DM dm, PetscSection pointSection, IS pointPartition, PetscSection *section, IS *partition)
120670034214SMatthew G. Knepley {
120770034214SMatthew G. Knepley   /* const PetscInt  height = 0; */
120870034214SMatthew G. Knepley   const PetscInt *partArray;
120970034214SMatthew G. Knepley   PetscInt       *allPoints, *packPoints;
121070034214SMatthew G. Knepley   PetscInt        rStart, rEnd, rank, pStart, pEnd, newSize;
121170034214SMatthew G. Knepley   PetscErrorCode  ierr;
121270034214SMatthew G. Knepley   PetscBT         bt;
121370034214SMatthew G. Knepley   PetscSegBuffer  segpack,segpart;
121470034214SMatthew G. Knepley 
121570034214SMatthew G. Knepley   PetscFunctionBegin;
121670034214SMatthew G. Knepley   ierr = PetscSectionGetChart(pointSection, &rStart, &rEnd);CHKERRQ(ierr);
121770034214SMatthew G. Knepley   ierr = ISGetIndices(pointPartition, &partArray);CHKERRQ(ierr);
121870034214SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), section);CHKERRQ(ierr);
121970034214SMatthew G. Knepley   ierr = PetscSectionSetChart(*section, rStart, rEnd);CHKERRQ(ierr);
122070034214SMatthew G. Knepley   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
122170034214SMatthew G. Knepley   ierr = PetscBTCreate(pEnd-pStart,&bt);CHKERRQ(ierr);
122270034214SMatthew G. Knepley   ierr = PetscSegBufferCreate(sizeof(PetscInt),1000,&segpack);CHKERRQ(ierr);
122370034214SMatthew G. Knepley   ierr = PetscSegBufferCreate(sizeof(PetscInt),1000,&segpart);CHKERRQ(ierr);
122470034214SMatthew G. Knepley   for (rank = rStart; rank < rEnd; ++rank) {
122570034214SMatthew G. Knepley     PetscInt partSize = 0, numPoints, offset, p, *PETSC_RESTRICT placePoints;
122670034214SMatthew G. Knepley 
122770034214SMatthew G. Knepley     ierr = PetscSectionGetDof(pointSection, rank, &numPoints);CHKERRQ(ierr);
122870034214SMatthew G. Knepley     ierr = PetscSectionGetOffset(pointSection, rank, &offset);CHKERRQ(ierr);
122970034214SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
123070034214SMatthew G. Knepley       PetscInt  point   = partArray[offset+p], closureSize, c;
123170034214SMatthew G. Knepley       PetscInt *closure = NULL;
123270034214SMatthew G. Knepley 
123370034214SMatthew G. Knepley       /* TODO Include support for height > 0 case */
123470034214SMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
123570034214SMatthew G. Knepley       for (c=0; c<closureSize; c++) {
123670034214SMatthew G. Knepley         PetscInt cpoint = closure[c*2];
123770034214SMatthew G. Knepley         if (!PetscBTLookupSet(bt,cpoint-pStart)) {
123870034214SMatthew G. Knepley           PetscInt *PETSC_RESTRICT pt;
123970034214SMatthew G. Knepley           partSize++;
124070034214SMatthew G. Knepley           ierr = PetscSegBufferGetInts(segpart,1,&pt);CHKERRQ(ierr);
124170034214SMatthew G. Knepley           *pt = cpoint;
124270034214SMatthew G. Knepley         }
12439ffc88e4SToby Isaac         ierr = DMPlexMarkTreeClosure(dm,segpart,bt,cpoint,&partSize);CHKERRQ(ierr);
124470034214SMatthew G. Knepley       }
124570034214SMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
124670034214SMatthew G. Knepley     }
124770034214SMatthew G. Knepley     ierr = PetscSectionSetDof(*section, rank, partSize);CHKERRQ(ierr);
124870034214SMatthew G. Knepley     ierr = PetscSegBufferGetInts(segpack,partSize,&placePoints);CHKERRQ(ierr);
124970034214SMatthew G. Knepley     ierr = PetscSegBufferExtractTo(segpart,placePoints);CHKERRQ(ierr);
125070034214SMatthew G. Knepley     ierr = PetscSortInt(partSize,placePoints);CHKERRQ(ierr);
125170034214SMatthew G. Knepley     for (p=0; p<partSize; p++) {ierr = PetscBTClear(bt,placePoints[p]-pStart);CHKERRQ(ierr);}
125270034214SMatthew G. Knepley   }
125370034214SMatthew G. Knepley   ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
125470034214SMatthew G. Knepley   ierr = PetscSegBufferDestroy(&segpart);CHKERRQ(ierr);
125570034214SMatthew G. Knepley 
125670034214SMatthew G. Knepley   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
125770034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(*section, &newSize);CHKERRQ(ierr);
125870034214SMatthew G. Knepley   ierr = PetscMalloc1(newSize, &allPoints);CHKERRQ(ierr);
125970034214SMatthew G. Knepley 
126070034214SMatthew G. Knepley   ierr = PetscSegBufferExtractInPlace(segpack,&packPoints);CHKERRQ(ierr);
126170034214SMatthew G. Knepley   for (rank = rStart; rank < rEnd; ++rank) {
126270034214SMatthew G. Knepley     PetscInt numPoints, offset;
126370034214SMatthew G. Knepley 
126470034214SMatthew G. Knepley     ierr = PetscSectionGetDof(*section, rank, &numPoints);CHKERRQ(ierr);
126570034214SMatthew G. Knepley     ierr = PetscSectionGetOffset(*section, rank, &offset);CHKERRQ(ierr);
126670034214SMatthew G. Knepley     ierr = PetscMemcpy(&allPoints[offset], packPoints, numPoints * sizeof(PetscInt));CHKERRQ(ierr);
126770034214SMatthew G. Knepley     packPoints += numPoints;
126870034214SMatthew G. Knepley   }
126970034214SMatthew G. Knepley 
127070034214SMatthew G. Knepley   ierr = PetscSegBufferDestroy(&segpack);CHKERRQ(ierr);
127170034214SMatthew G. Knepley   ierr = ISRestoreIndices(pointPartition, &partArray);CHKERRQ(ierr);
127270034214SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), newSize, allPoints, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
127370034214SMatthew G. Knepley   PetscFunctionReturn(0);
127470034214SMatthew G. Knepley }
1275aa3148a8SMichael Lange 
1276aa3148a8SMichael Lange #undef __FUNCT__
12775abbe4feSMichael Lange #define __FUNCT__ "DMPlexPartitionLabelClosure"
12785abbe4feSMichael Lange /*@
12795abbe4feSMichael Lange   DMPlexPartitionLabelClosure - Add the closure of all points to the partition label
12805abbe4feSMichael Lange 
12815abbe4feSMichael Lange   Input Parameters:
12825abbe4feSMichael Lange + dm     - The DM
12835abbe4feSMichael Lange - label  - DMLabel assinging ranks to remote roots
12845abbe4feSMichael Lange 
12855abbe4feSMichael Lange   Level: developer
12865abbe4feSMichael Lange 
12875abbe4feSMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap
12885abbe4feSMichael Lange @*/
12895abbe4feSMichael Lange PetscErrorCode DMPlexPartitionLabelClosure(DM dm, DMLabel label)
12905abbe4feSMichael Lange {
12915abbe4feSMichael Lange   IS              rankIS,   pointIS;
12925abbe4feSMichael Lange   const PetscInt *ranks,   *points;
12935abbe4feSMichael Lange   PetscInt        numRanks, numPoints, r, p, c, closureSize;
12945abbe4feSMichael Lange   PetscInt       *closure = NULL;
12955abbe4feSMichael Lange   PetscErrorCode  ierr;
12965abbe4feSMichael Lange 
12975abbe4feSMichael Lange   PetscFunctionBegin;
12985abbe4feSMichael Lange   ierr = DMLabelGetValueIS(label, &rankIS);CHKERRQ(ierr);
12995abbe4feSMichael Lange   ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr);
13005abbe4feSMichael Lange   ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr);
13015abbe4feSMichael Lange   for (r = 0; r < numRanks; ++r) {
13025abbe4feSMichael Lange     const PetscInt rank = ranks[r];
13035abbe4feSMichael Lange 
13045abbe4feSMichael Lange     ierr = DMLabelGetStratumIS(label, rank, &pointIS);CHKERRQ(ierr);
13055abbe4feSMichael Lange     ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
13065abbe4feSMichael Lange     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
13075abbe4feSMichael Lange     for (p = 0; p < numPoints; ++p) {
13085abbe4feSMichael Lange       ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
13095abbe4feSMichael Lange       for (c = 0; c < closureSize*2; c += 2) {ierr = DMLabelSetValue(label, closure[c], rank);CHKERRQ(ierr);}
13105abbe4feSMichael Lange     }
13115abbe4feSMichael Lange     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
13125abbe4feSMichael Lange     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
13135abbe4feSMichael Lange   }
13145abbe4feSMichael Lange   if (closure) {ierr = DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);}
13155abbe4feSMichael Lange   ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr);
13165abbe4feSMichael Lange   ierr = ISDestroy(&rankIS);CHKERRQ(ierr);
13175abbe4feSMichael Lange   PetscFunctionReturn(0);
13185abbe4feSMichael Lange }
13195abbe4feSMichael Lange 
13205abbe4feSMichael Lange #undef __FUNCT__
13211fd9873aSMichael Lange #define __FUNCT__ "DMPlexPartitionLabelInvert"
13221fd9873aSMichael Lange /*@
13231fd9873aSMichael Lange   DMPlexPartitionLabelInvert - Create a partition label of remote roots from a local root label
13241fd9873aSMichael Lange 
13251fd9873aSMichael Lange   Input Parameters:
13261fd9873aSMichael Lange + dm        - The DM
13271fd9873aSMichael Lange . rootLabel - DMLabel assinging ranks to local roots
13281fd9873aSMichael Lange . processSF - A star forest mapping into the local index on each remote rank
13291fd9873aSMichael Lange 
13301fd9873aSMichael Lange   Output Parameter:
13311fd9873aSMichael Lange - leafLabel - DMLabel assinging ranks to remote roots
13321fd9873aSMichael Lange 
13331fd9873aSMichael Lange   Note: The rootLabel defines a send pattern by mapping local points to remote target ranks. The
13341fd9873aSMichael Lange   resulting leafLabel is a receiver mapping of remote roots to their parent rank.
13351fd9873aSMichael Lange 
13361fd9873aSMichael Lange   Level: developer
13371fd9873aSMichael Lange 
13381fd9873aSMichael Lange .seealso: DMPlexPartitionLabelCreateSF, DMPlexDistribute(), DMPlexCreateOverlap
13391fd9873aSMichael Lange @*/
13401fd9873aSMichael Lange PetscErrorCode DMPlexPartitionLabelInvert(DM dm, DMLabel rootLabel, PetscSF processSF, DMLabel leafLabel)
13411fd9873aSMichael Lange {
13421fd9873aSMichael Lange   MPI_Comm           comm;
13431fd9873aSMichael Lange   PetscMPIInt        rank, numProcs;
13441fd9873aSMichael Lange   PetscInt           p, n, numNeighbors, size, l, nleaves;
13451fd9873aSMichael Lange   PetscSF            sfPoint;
13461fd9873aSMichael Lange   PetscSFNode       *rootPoints, *leafPoints;
13471fd9873aSMichael Lange   PetscSection       rootSection, leafSection;
13481fd9873aSMichael Lange   const PetscSFNode *remote;
13491fd9873aSMichael Lange   const PetscInt    *local, *neighbors;
13501fd9873aSMichael Lange   IS                 valueIS;
13511fd9873aSMichael Lange   PetscErrorCode     ierr;
13521fd9873aSMichael Lange 
13531fd9873aSMichael Lange   PetscFunctionBegin;
13541fd9873aSMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
13551fd9873aSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
13561fd9873aSMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
13571fd9873aSMichael Lange   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
13581fd9873aSMichael Lange 
13591fd9873aSMichael Lange   /* Convert to (point, rank) and use actual owners */
13601fd9873aSMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
13611fd9873aSMichael Lange   ierr = PetscSectionSetChart(rootSection, 0, numProcs);CHKERRQ(ierr);
13621fd9873aSMichael Lange   ierr = DMLabelGetValueIS(rootLabel, &valueIS);CHKERRQ(ierr);
13631fd9873aSMichael Lange   ierr = ISGetLocalSize(valueIS, &numNeighbors);CHKERRQ(ierr);
13641fd9873aSMichael Lange   ierr = ISGetIndices(valueIS, &neighbors);CHKERRQ(ierr);
13651fd9873aSMichael Lange   for (n = 0; n < numNeighbors; ++n) {
13661fd9873aSMichael Lange     PetscInt numPoints;
13671fd9873aSMichael Lange 
13681fd9873aSMichael Lange     ierr = DMLabelGetStratumSize(rootLabel, neighbors[n], &numPoints);CHKERRQ(ierr);
13691fd9873aSMichael Lange     ierr = PetscSectionAddDof(rootSection, neighbors[n], numPoints);CHKERRQ(ierr);
13701fd9873aSMichael Lange   }
13711fd9873aSMichael Lange   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
13721fd9873aSMichael Lange   ierr = PetscSectionGetStorageSize(rootSection, &size);CHKERRQ(ierr);
13731fd9873aSMichael Lange   ierr = PetscMalloc1(size, &rootPoints);CHKERRQ(ierr);
13741fd9873aSMichael Lange   ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr);
13751fd9873aSMichael Lange   for (n = 0; n < numNeighbors; ++n) {
13761fd9873aSMichael Lange     IS              pointIS;
13771fd9873aSMichael Lange     const PetscInt *points;
13781fd9873aSMichael Lange     PetscInt        off, numPoints, p;
13791fd9873aSMichael Lange 
13801fd9873aSMichael Lange     ierr = PetscSectionGetOffset(rootSection, neighbors[n], &off);CHKERRQ(ierr);
13811fd9873aSMichael Lange     ierr = DMLabelGetStratumIS(rootLabel, neighbors[n], &pointIS);CHKERRQ(ierr);
13821fd9873aSMichael Lange     ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
13831fd9873aSMichael Lange     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
13841fd9873aSMichael Lange     for (p = 0; p < numPoints; ++p) {
1385*f8987ae8SMichael Lange       if (local) {ierr = PetscFindInt(points[p], nleaves, local, &l);CHKERRQ(ierr);}
1386*f8987ae8SMichael Lange       else       {l = -1;}
13871fd9873aSMichael Lange       if (l >= 0) {rootPoints[off+p] = remote[l];}
13881fd9873aSMichael Lange       else        {rootPoints[off+p].index = points[p]; rootPoints[off+p].rank = rank;}
13891fd9873aSMichael Lange     }
13901fd9873aSMichael Lange     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
13911fd9873aSMichael Lange     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
13921fd9873aSMichael Lange   }
13931fd9873aSMichael Lange   ierr = ISRestoreIndices(valueIS, &neighbors);CHKERRQ(ierr);
13941fd9873aSMichael Lange   ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
13951fd9873aSMichael Lange   /* Communicate overlap */
13961fd9873aSMichael Lange   ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr);
13971fd9873aSMichael Lange   ierr = DMPlexDistributeData(dm, processSF, rootSection, MPIU_2INT, rootPoints, leafSection, (void**) &leafPoints);CHKERRQ(ierr);
13981fd9873aSMichael Lange   /* Filter remote contributions (ovLeafPoints) into the overlapSF */
13991fd9873aSMichael Lange   ierr = PetscSectionGetStorageSize(leafSection, &size);CHKERRQ(ierr);
14001fd9873aSMichael Lange   for (p = 0; p < size; p++) {
14011fd9873aSMichael Lange     ierr = DMLabelSetValue(leafLabel, leafPoints[p].index, leafPoints[p].rank);CHKERRQ(ierr);
14021fd9873aSMichael Lange   }
14031fd9873aSMichael Lange   ierr = PetscFree(rootPoints);CHKERRQ(ierr);
14041fd9873aSMichael Lange   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
14051fd9873aSMichael Lange   ierr = PetscFree(leafPoints);CHKERRQ(ierr);
14061fd9873aSMichael Lange   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
14071fd9873aSMichael Lange   PetscFunctionReturn(0);
14081fd9873aSMichael Lange }
14091fd9873aSMichael Lange 
14101fd9873aSMichael Lange #undef __FUNCT__
1411aa3148a8SMichael Lange #define __FUNCT__ "DMPlexPartitionLabelCreateSF"
1412aa3148a8SMichael Lange /*@
1413aa3148a8SMichael Lange   DMPlexPartitionLabelCreateSF - Create a star forest from a label that assigns ranks to points
1414aa3148a8SMichael Lange 
1415aa3148a8SMichael Lange   Input Parameters:
1416aa3148a8SMichael Lange + dm    - The DM
1417aa3148a8SMichael Lange . label - DMLabel assinging ranks to remote roots
1418aa3148a8SMichael Lange 
1419aa3148a8SMichael Lange   Output Parameter:
1420aa3148a8SMichael Lange - sf    - The star forest communication context encapsulating the defined mapping
1421aa3148a8SMichael Lange 
1422aa3148a8SMichael Lange   Note: The incoming label is a receiver mapping of remote points to their parent rank.
1423aa3148a8SMichael Lange 
1424aa3148a8SMichael Lange   Level: developer
1425aa3148a8SMichael Lange 
1426aa3148a8SMichael Lange .seealso: DMPlexDistribute(), DMPlexCreateOverlap
1427aa3148a8SMichael Lange @*/
1428aa3148a8SMichael Lange PetscErrorCode DMPlexPartitionLabelCreateSF(DM dm, DMLabel label, PetscSF *sf)
1429aa3148a8SMichael Lange {
1430aa3148a8SMichael Lange   PetscMPIInt    numProcs;
1431aa3148a8SMichael Lange   PetscInt       n, idx, numRemote, p, numPoints, pStart, pEnd;
1432aa3148a8SMichael Lange   PetscSFNode   *remotePoints;
1433aa3148a8SMichael Lange   PetscErrorCode ierr;
1434aa3148a8SMichael Lange 
1435aa3148a8SMichael Lange   PetscFunctionBegin;
1436aa3148a8SMichael Lange   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr);
1437aa3148a8SMichael Lange 
1438aa3148a8SMichael Lange   for (numRemote = 0, n = 0; n < numProcs; ++n) {
1439aa3148a8SMichael Lange     ierr = DMLabelGetStratumSize(label, n, &numPoints);CHKERRQ(ierr);
1440aa3148a8SMichael Lange     numRemote += numPoints;
1441aa3148a8SMichael Lange   }
1442aa3148a8SMichael Lange   ierr = PetscMalloc1(numRemote, &remotePoints);CHKERRQ(ierr);
1443aa3148a8SMichael Lange   for (idx = 0, n = 0; n < numProcs; ++n) {
1444aa3148a8SMichael Lange     IS remoteRootIS;
1445aa3148a8SMichael Lange     const PetscInt *remoteRoots;
1446aa3148a8SMichael Lange     ierr = DMLabelGetStratumSize(label, n, &numPoints);CHKERRQ(ierr);
1447aa3148a8SMichael Lange     if (numPoints <= 0) continue;
1448aa3148a8SMichael Lange     ierr = DMLabelGetStratumIS(label, n, &remoteRootIS);CHKERRQ(ierr);
1449aa3148a8SMichael Lange     ierr = ISGetIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
1450aa3148a8SMichael Lange     for (p = 0; p < numPoints; p++) {
1451aa3148a8SMichael Lange       remotePoints[idx].index = remoteRoots[p];
1452aa3148a8SMichael Lange       remotePoints[idx].rank = n;
1453aa3148a8SMichael Lange       idx++;
1454aa3148a8SMichael Lange     }
1455aa3148a8SMichael Lange     ierr = ISRestoreIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
1456aa3148a8SMichael Lange     ierr = ISDestroy(&remoteRootIS);CHKERRQ(ierr);
1457aa3148a8SMichael Lange   }
1458aa3148a8SMichael Lange   ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), sf);CHKERRQ(ierr);
1459aa3148a8SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
1460aa3148a8SMichael Lange   ierr = PetscSFSetGraph(*sf, pEnd-pStart, numRemote, NULL, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
1461aa3148a8SMichael Lange   PetscFunctionReturn(0);
1462aa3148a8SMichael Lange }
1463