xref: /petsc/src/dm/impls/plex/plexpartition.c (revision 9ffc88e4009c11ed2a29b77417810ef6e4b738fe)
170034214SMatthew G. Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
270034214SMatthew G. Knepley 
370034214SMatthew G. Knepley #undef __FUNCT__
470034214SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateNeighborCSR"
570034214SMatthew G. Knepley PetscErrorCode DMPlexCreateNeighborCSR(DM dm, PetscInt cellHeight, PetscInt *numVertices, PetscInt **offsets, PetscInt **adjacency)
670034214SMatthew G. Knepley {
770034214SMatthew G. Knepley   const PetscInt maxFaceCases = 30;
870034214SMatthew G. Knepley   PetscInt       numFaceCases = 0;
970034214SMatthew G. Knepley   PetscInt       numFaceVertices[30]; /* maxFaceCases, C89 sucks sucks sucks */
1070034214SMatthew G. Knepley   PetscInt      *off, *adj;
1170034214SMatthew G. Knepley   PetscInt      *neighborCells = NULL;
1270034214SMatthew G. Knepley   PetscInt       dim, cellDim, depth = 0, faceDepth, cStart, cEnd, c, numCells, cell;
1370034214SMatthew G. Knepley   PetscErrorCode ierr;
1470034214SMatthew G. Knepley 
1570034214SMatthew G. Knepley   PetscFunctionBegin;
1670034214SMatthew G. Knepley   /* For parallel partitioning, I think you have to communicate supports */
17c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1870034214SMatthew G. Knepley   cellDim = dim - cellHeight;
1970034214SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
2070034214SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
2170034214SMatthew G. Knepley   if (cEnd - cStart == 0) {
2270034214SMatthew G. Knepley     if (numVertices) *numVertices = 0;
2370034214SMatthew G. Knepley     if (offsets)   *offsets   = NULL;
2470034214SMatthew G. Knepley     if (adjacency) *adjacency = NULL;
2570034214SMatthew G. Knepley     PetscFunctionReturn(0);
2670034214SMatthew G. Knepley   }
2770034214SMatthew G. Knepley   numCells  = cEnd - cStart;
2870034214SMatthew G. Knepley   faceDepth = depth - cellHeight;
2970034214SMatthew G. Knepley   if (dim == depth) {
3070034214SMatthew G. Knepley     PetscInt f, fStart, fEnd;
3170034214SMatthew G. Knepley 
3270034214SMatthew G. Knepley     ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr);
3370034214SMatthew G. Knepley     /* Count neighboring cells */
3470034214SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, cellHeight+1, &fStart, &fEnd);CHKERRQ(ierr);
3570034214SMatthew G. Knepley     for (f = fStart; f < fEnd; ++f) {
3670034214SMatthew G. Knepley       const PetscInt *support;
3770034214SMatthew G. Knepley       PetscInt        supportSize;
3870034214SMatthew G. Knepley       ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr);
3970034214SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr);
4070034214SMatthew G. Knepley       if (supportSize == 2) {
41*9ffc88e4SToby Isaac         PetscInt numChildren;
42*9ffc88e4SToby Isaac 
43*9ffc88e4SToby Isaac         ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr);
44*9ffc88e4SToby Isaac         if (!numChildren) {
4570034214SMatthew G. Knepley           ++off[support[0]-cStart+1];
4670034214SMatthew G. Knepley           ++off[support[1]-cStart+1];
4770034214SMatthew G. Knepley         }
4870034214SMatthew G. Knepley       }
49*9ffc88e4SToby Isaac     }
5070034214SMatthew G. Knepley     /* Prefix sum */
5170034214SMatthew G. Knepley     for (c = 1; c <= numCells; ++c) off[c] += off[c-1];
5270034214SMatthew G. Knepley     if (adjacency) {
5370034214SMatthew G. Knepley       PetscInt *tmp;
5470034214SMatthew G. Knepley 
5570034214SMatthew G. Knepley       ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr);
5670034214SMatthew G. Knepley       ierr = PetscMalloc1((numCells+1), &tmp);CHKERRQ(ierr);
5770034214SMatthew G. Knepley       ierr = PetscMemcpy(tmp, off, (numCells+1) * sizeof(PetscInt));CHKERRQ(ierr);
5870034214SMatthew G. Knepley       /* Get neighboring cells */
5970034214SMatthew G. Knepley       for (f = fStart; f < fEnd; ++f) {
6070034214SMatthew G. Knepley         const PetscInt *support;
6170034214SMatthew G. Knepley         PetscInt        supportSize;
6270034214SMatthew G. Knepley         ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr);
6370034214SMatthew G. Knepley         ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr);
6470034214SMatthew G. Knepley         if (supportSize == 2) {
65*9ffc88e4SToby Isaac           PetscInt numChildren;
66*9ffc88e4SToby Isaac 
67*9ffc88e4SToby Isaac           ierr = DMPlexGetTreeChildren(dm,f,&numChildren,NULL);CHKERRQ(ierr);
68*9ffc88e4SToby Isaac           if (!numChildren) {
6970034214SMatthew G. Knepley             adj[tmp[support[0]-cStart]++] = support[1];
7070034214SMatthew G. Knepley             adj[tmp[support[1]-cStart]++] = support[0];
7170034214SMatthew G. Knepley           }
7270034214SMatthew G. Knepley         }
73*9ffc88e4SToby Isaac       }
7470034214SMatthew G. Knepley #if defined(PETSC_USE_DEBUG)
7570034214SMatthew 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);
7670034214SMatthew G. Knepley #endif
7770034214SMatthew G. Knepley       ierr = PetscFree(tmp);CHKERRQ(ierr);
7870034214SMatthew G. Knepley     }
7970034214SMatthew G. Knepley     if (numVertices) *numVertices = numCells;
8070034214SMatthew G. Knepley     if (offsets)   *offsets   = off;
8170034214SMatthew G. Knepley     if (adjacency) *adjacency = adj;
8270034214SMatthew G. Knepley     PetscFunctionReturn(0);
8370034214SMatthew G. Knepley   }
8470034214SMatthew G. Knepley   /* Setup face recognition */
8570034214SMatthew G. Knepley   if (faceDepth == 1) {
8670034214SMatthew 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 */
8770034214SMatthew G. Knepley 
8870034214SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
8970034214SMatthew G. Knepley       PetscInt corners;
9070034214SMatthew G. Knepley 
9170034214SMatthew G. Knepley       ierr = DMPlexGetConeSize(dm, c, &corners);CHKERRQ(ierr);
9270034214SMatthew G. Knepley       if (!cornersSeen[corners]) {
9370034214SMatthew G. Knepley         PetscInt nFV;
9470034214SMatthew G. Knepley 
9570034214SMatthew G. Knepley         if (numFaceCases >= maxFaceCases) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Exceeded maximum number of face recognition cases");
9670034214SMatthew G. Knepley         cornersSeen[corners] = 1;
9770034214SMatthew G. Knepley 
9870034214SMatthew G. Knepley         ierr = DMPlexGetNumFaceVertices(dm, cellDim, corners, &nFV);CHKERRQ(ierr);
9970034214SMatthew G. Knepley 
10070034214SMatthew G. Knepley         numFaceVertices[numFaceCases++] = nFV;
10170034214SMatthew G. Knepley       }
10270034214SMatthew G. Knepley     }
10370034214SMatthew G. Knepley   }
10470034214SMatthew G. Knepley   ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr);
10570034214SMatthew G. Knepley   /* Count neighboring cells */
10670034214SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) {
10770034214SMatthew G. Knepley     PetscInt numNeighbors = PETSC_DETERMINE, n;
10870034214SMatthew G. Knepley 
1098b0b4c70SToby Isaac     ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &numNeighbors, &neighborCells);CHKERRQ(ierr);
11070034214SMatthew G. Knepley     /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */
11170034214SMatthew G. Knepley     for (n = 0; n < numNeighbors; ++n) {
11270034214SMatthew G. Knepley       PetscInt        cellPair[2];
11370034214SMatthew G. Knepley       PetscBool       found    = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE;
11470034214SMatthew G. Knepley       PetscInt        meetSize = 0;
11570034214SMatthew G. Knepley       const PetscInt *meet    = NULL;
11670034214SMatthew G. Knepley 
11770034214SMatthew G. Knepley       cellPair[0] = cell; cellPair[1] = neighborCells[n];
11870034214SMatthew G. Knepley       if (cellPair[0] == cellPair[1]) continue;
11970034214SMatthew G. Knepley       if (!found) {
12070034214SMatthew G. Knepley         ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
12170034214SMatthew G. Knepley         if (meetSize) {
12270034214SMatthew G. Knepley           PetscInt f;
12370034214SMatthew G. Knepley 
12470034214SMatthew G. Knepley           for (f = 0; f < numFaceCases; ++f) {
12570034214SMatthew G. Knepley             if (numFaceVertices[f] == meetSize) {
12670034214SMatthew G. Knepley               found = PETSC_TRUE;
12770034214SMatthew G. Knepley               break;
12870034214SMatthew G. Knepley             }
12970034214SMatthew G. Knepley           }
13070034214SMatthew G. Knepley         }
13170034214SMatthew G. Knepley         ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
13270034214SMatthew G. Knepley       }
13370034214SMatthew G. Knepley       if (found) ++off[cell-cStart+1];
13470034214SMatthew G. Knepley     }
13570034214SMatthew G. Knepley   }
13670034214SMatthew G. Knepley   /* Prefix sum */
13770034214SMatthew G. Knepley   for (cell = 1; cell <= numCells; ++cell) off[cell] += off[cell-1];
13870034214SMatthew G. Knepley 
13970034214SMatthew G. Knepley   if (adjacency) {
14070034214SMatthew G. Knepley     ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr);
14170034214SMatthew G. Knepley     /* Get neighboring cells */
14270034214SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
14370034214SMatthew G. Knepley       PetscInt numNeighbors = PETSC_DETERMINE, n;
14470034214SMatthew G. Knepley       PetscInt cellOffset   = 0;
14570034214SMatthew G. Knepley 
1468b0b4c70SToby Isaac       ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, PETSC_FALSE, &numNeighbors, &neighborCells);CHKERRQ(ierr);
14770034214SMatthew G. Knepley       /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */
14870034214SMatthew G. Knepley       for (n = 0; n < numNeighbors; ++n) {
14970034214SMatthew G. Knepley         PetscInt        cellPair[2];
15070034214SMatthew G. Knepley         PetscBool       found    = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE;
15170034214SMatthew G. Knepley         PetscInt        meetSize = 0;
15270034214SMatthew G. Knepley         const PetscInt *meet    = NULL;
15370034214SMatthew G. Knepley 
15470034214SMatthew G. Knepley         cellPair[0] = cell; cellPair[1] = neighborCells[n];
15570034214SMatthew G. Knepley         if (cellPair[0] == cellPair[1]) continue;
15670034214SMatthew G. Knepley         if (!found) {
15770034214SMatthew G. Knepley           ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
15870034214SMatthew G. Knepley           if (meetSize) {
15970034214SMatthew G. Knepley             PetscInt f;
16070034214SMatthew G. Knepley 
16170034214SMatthew G. Knepley             for (f = 0; f < numFaceCases; ++f) {
16270034214SMatthew G. Knepley               if (numFaceVertices[f] == meetSize) {
16370034214SMatthew G. Knepley                 found = PETSC_TRUE;
16470034214SMatthew G. Knepley                 break;
16570034214SMatthew G. Knepley               }
16670034214SMatthew G. Knepley             }
16770034214SMatthew G. Knepley           }
16870034214SMatthew G. Knepley           ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
16970034214SMatthew G. Knepley         }
17070034214SMatthew G. Knepley         if (found) {
17170034214SMatthew G. Knepley           adj[off[cell-cStart]+cellOffset] = neighborCells[n];
17270034214SMatthew G. Knepley           ++cellOffset;
17370034214SMatthew G. Knepley         }
17470034214SMatthew G. Knepley       }
17570034214SMatthew G. Knepley     }
17670034214SMatthew G. Knepley   }
17770034214SMatthew G. Knepley   ierr = PetscFree(neighborCells);CHKERRQ(ierr);
17870034214SMatthew G. Knepley   if (numVertices) *numVertices = numCells;
17970034214SMatthew G. Knepley   if (offsets)   *offsets   = off;
18070034214SMatthew G. Knepley   if (adjacency) *adjacency = adj;
18170034214SMatthew G. Knepley   PetscFunctionReturn(0);
18270034214SMatthew G. Knepley }
18370034214SMatthew G. Knepley 
18470034214SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO)
18570034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
18670034214SMatthew G. Knepley #include <unistd.h>
18770034214SMatthew G. Knepley #endif
18870034214SMatthew G. Knepley /* Chaco does not have an include file */
18970034214SMatthew G. Knepley PETSC_EXTERN int interface(int nvtxs, int *start, int *adjacency, int *vwgts,
19070034214SMatthew G. Knepley                        float *ewgts, float *x, float *y, float *z, char *outassignname,
19170034214SMatthew G. Knepley                        char *outfilename, short *assignment, int architecture, int ndims_tot,
19270034214SMatthew G. Knepley                        int mesh_dims[3], double *goal, int global_method, int local_method,
19370034214SMatthew G. Knepley                        int rqi_flag, int vmax, int ndims, double eigtol, long seed);
19470034214SMatthew G. Knepley 
19570034214SMatthew G. Knepley extern int FREE_GRAPH;
19670034214SMatthew G. Knepley 
19770034214SMatthew G. Knepley #undef __FUNCT__
19870034214SMatthew G. Knepley #define __FUNCT__ "DMPlexPartition_Chaco"
19970034214SMatthew G. Knepley PetscErrorCode DMPlexPartition_Chaco(DM dm, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection *partSection, IS *partition)
20070034214SMatthew G. Knepley {
20170034214SMatthew G. Knepley   enum {DEFAULT_METHOD = 1, INERTIAL_METHOD = 3};
20270034214SMatthew G. Knepley   MPI_Comm       comm;
20370034214SMatthew G. Knepley   int            nvtxs          = numVertices; /* number of vertices in full graph */
20470034214SMatthew G. Knepley   int           *vwgts          = NULL;   /* weights for all vertices */
20570034214SMatthew G. Knepley   float         *ewgts          = NULL;   /* weights for all edges */
20670034214SMatthew G. Knepley   float         *x              = NULL, *y = NULL, *z = NULL; /* coordinates for inertial method */
20770034214SMatthew G. Knepley   char          *outassignname  = NULL;   /*  name of assignment output file */
20870034214SMatthew G. Knepley   char          *outfilename    = NULL;   /* output file name */
20970034214SMatthew G. Knepley   int            architecture   = 1;      /* 0 => hypercube, d => d-dimensional mesh */
21070034214SMatthew G. Knepley   int            ndims_tot      = 0;      /* total number of cube dimensions to divide */
21170034214SMatthew G. Knepley   int            mesh_dims[3];            /* dimensions of mesh of processors */
21270034214SMatthew G. Knepley   double        *goal          = NULL;    /* desired set sizes for each set */
21370034214SMatthew G. Knepley   int            global_method = 1;       /* global partitioning algorithm */
21470034214SMatthew G. Knepley   int            local_method  = 1;       /* local partitioning algorithm */
21570034214SMatthew G. Knepley   int            rqi_flag      = 0;       /* should I use RQI/Symmlq eigensolver? */
21670034214SMatthew G. Knepley   int            vmax          = 200;     /* how many vertices to coarsen down to? */
21770034214SMatthew G. Knepley   int            ndims         = 1;       /* number of eigenvectors (2^d sets) */
21870034214SMatthew G. Knepley   double         eigtol        = 0.001;   /* tolerance on eigenvectors */
21970034214SMatthew G. Knepley   long           seed          = 123636512; /* for random graph mutations */
22070034214SMatthew G. Knepley   short int     *assignment;              /* Output partition */
22170034214SMatthew G. Knepley   int            fd_stdout, fd_pipe[2];
22270034214SMatthew G. Knepley   PetscInt      *points;
22370034214SMatthew G. Knepley   PetscMPIInt    commSize;
22470034214SMatthew G. Knepley   int            i, v, p;
22570034214SMatthew G. Knepley   PetscErrorCode ierr;
22670034214SMatthew G. Knepley 
22770034214SMatthew G. Knepley   PetscFunctionBegin;
22870034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
22970034214SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &commSize);CHKERRQ(ierr);
23070034214SMatthew G. Knepley   if (!numVertices) {
23170034214SMatthew G. Knepley     ierr = PetscSectionCreate(comm, partSection);CHKERRQ(ierr);
23270034214SMatthew G. Knepley     ierr = PetscSectionSetChart(*partSection, 0, commSize);CHKERRQ(ierr);
23370034214SMatthew G. Knepley     ierr = PetscSectionSetUp(*partSection);CHKERRQ(ierr);
23470034214SMatthew G. Knepley     ierr = ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
23570034214SMatthew G. Knepley     PetscFunctionReturn(0);
23670034214SMatthew G. Knepley   }
23770034214SMatthew G. Knepley   FREE_GRAPH = 0;                         /* Do not let Chaco free my memory */
23870034214SMatthew G. Knepley   for (i = 0; i < start[numVertices]; ++i) ++adjacency[i];
23970034214SMatthew G. Knepley 
24070034214SMatthew G. Knepley   if (global_method == INERTIAL_METHOD) {
24170034214SMatthew G. Knepley     /* manager.createCellCoordinates(nvtxs, &x, &y, &z); */
24270034214SMatthew G. Knepley     SETERRQ(comm, PETSC_ERR_SUP, "Inertial partitioning not yet supported");
24370034214SMatthew G. Knepley   }
24470034214SMatthew G. Knepley   mesh_dims[0] = commSize;
24570034214SMatthew G. Knepley   mesh_dims[1] = 1;
24670034214SMatthew G. Knepley   mesh_dims[2] = 1;
24770034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &assignment);CHKERRQ(ierr);
24870034214SMatthew G. Knepley   /* Chaco outputs to stdout. We redirect this to a buffer. */
24970034214SMatthew G. Knepley   /* TODO: check error codes for UNIX calls */
25070034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
25170034214SMatthew G. Knepley   {
25270034214SMatthew G. Knepley     int piperet;
25370034214SMatthew G. Knepley     piperet = pipe(fd_pipe);
25470034214SMatthew G. Knepley     if (piperet) SETERRQ(comm,PETSC_ERR_SYS,"Could not create pipe");
25570034214SMatthew G. Knepley     fd_stdout = dup(1);
25670034214SMatthew G. Knepley     close(1);
25770034214SMatthew G. Knepley     dup2(fd_pipe[1], 1);
25870034214SMatthew G. Knepley   }
25970034214SMatthew G. Knepley #endif
26070034214SMatthew G. Knepley   ierr = interface(nvtxs, (int*) start, (int*) adjacency, vwgts, ewgts, x, y, z, outassignname, outfilename,
26170034214SMatthew G. Knepley                    assignment, architecture, ndims_tot, mesh_dims, goal, global_method, local_method, rqi_flag,
26270034214SMatthew G. Knepley                    vmax, ndims, eigtol, seed);
26370034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
26470034214SMatthew G. Knepley   {
26570034214SMatthew G. Knepley     char msgLog[10000];
26670034214SMatthew G. Knepley     int  count;
26770034214SMatthew G. Knepley 
26870034214SMatthew G. Knepley     fflush(stdout);
26970034214SMatthew G. Knepley     count = read(fd_pipe[0], msgLog, (10000-1)*sizeof(char));
27070034214SMatthew G. Knepley     if (count < 0) count = 0;
27170034214SMatthew G. Knepley     msgLog[count] = 0;
27270034214SMatthew G. Knepley     close(1);
27370034214SMatthew G. Knepley     dup2(fd_stdout, 1);
27470034214SMatthew G. Knepley     close(fd_stdout);
27570034214SMatthew G. Knepley     close(fd_pipe[0]);
27670034214SMatthew G. Knepley     close(fd_pipe[1]);
27770034214SMatthew G. Knepley     if (ierr) SETERRQ1(comm, PETSC_ERR_LIB, "Error in Chaco library: %s", msgLog);
27870034214SMatthew G. Knepley   }
27970034214SMatthew G. Knepley #endif
28070034214SMatthew G. Knepley   /* Convert to PetscSection+IS */
28170034214SMatthew G. Knepley   ierr = PetscSectionCreate(comm, partSection);CHKERRQ(ierr);
28270034214SMatthew G. Knepley   ierr = PetscSectionSetChart(*partSection, 0, commSize);CHKERRQ(ierr);
28370034214SMatthew G. Knepley   for (v = 0; v < nvtxs; ++v) {
28470034214SMatthew G. Knepley     ierr = PetscSectionAddDof(*partSection, assignment[v], 1);CHKERRQ(ierr);
28570034214SMatthew G. Knepley   }
28670034214SMatthew G. Knepley   ierr = PetscSectionSetUp(*partSection);CHKERRQ(ierr);
28770034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr);
28870034214SMatthew G. Knepley   for (p = 0, i = 0; p < commSize; ++p) {
28970034214SMatthew G. Knepley     for (v = 0; v < nvtxs; ++v) {
29070034214SMatthew G. Knepley       if (assignment[v] == p) points[i++] = v;
29170034214SMatthew G. Knepley     }
29270034214SMatthew G. Knepley   }
29370034214SMatthew G. Knepley   if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs);
29470034214SMatthew G. Knepley   ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
29570034214SMatthew G. Knepley   if (global_method == INERTIAL_METHOD) {
29670034214SMatthew G. Knepley     /* manager.destroyCellCoordinates(nvtxs, &x, &y, &z); */
29770034214SMatthew G. Knepley   }
29870034214SMatthew G. Knepley   ierr = PetscFree(assignment);CHKERRQ(ierr);
29970034214SMatthew G. Knepley   for (i = 0; i < start[numVertices]; ++i) --adjacency[i];
30070034214SMatthew G. Knepley   PetscFunctionReturn(0);
30170034214SMatthew G. Knepley }
30270034214SMatthew G. Knepley #endif
30370034214SMatthew G. Knepley 
30470034214SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS)
30570034214SMatthew G. Knepley #include <parmetis.h>
30670034214SMatthew G. Knepley 
30770034214SMatthew G. Knepley #undef __FUNCT__
30870034214SMatthew G. Knepley #define __FUNCT__ "DMPlexPartition_ParMetis"
30970034214SMatthew G. Knepley PetscErrorCode DMPlexPartition_ParMetis(DM dm, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection *partSection, IS *partition)
31070034214SMatthew G. Knepley {
31170034214SMatthew G. Knepley   MPI_Comm       comm;
31270034214SMatthew G. Knepley   PetscInt       nvtxs      = numVertices; /* The number of vertices in full graph */
31370034214SMatthew G. Knepley   PetscInt      *vtxdist;                  /* Distribution of vertices across processes */
31470034214SMatthew G. Knepley   PetscInt      *xadj       = start;       /* Start of edge list for each vertex */
31570034214SMatthew G. Knepley   PetscInt      *adjncy     = adjacency;   /* Edge lists for all vertices */
31670034214SMatthew G. Knepley   PetscInt      *vwgt       = NULL;        /* Vertex weights */
31770034214SMatthew G. Knepley   PetscInt      *adjwgt     = NULL;        /* Edge weights */
31870034214SMatthew G. Knepley   PetscInt       wgtflag    = 0;           /* Indicates which weights are present */
31970034214SMatthew G. Knepley   PetscInt       numflag    = 0;           /* Indicates initial offset (0 or 1) */
32070034214SMatthew G. Knepley   PetscInt       ncon       = 1;           /* The number of weights per vertex */
32170034214SMatthew G. Knepley   PetscInt       nparts;                   /* The number of partitions */
32270034214SMatthew G. Knepley   PetscReal     *tpwgts;                   /* The fraction of vertex weights assigned to each partition */
32370034214SMatthew G. Knepley   PetscReal     *ubvec;                    /* The balance intolerance for vertex weights */
32470034214SMatthew G. Knepley   PetscInt       options[5];               /* Options */
32570034214SMatthew G. Knepley   /* Outputs */
32670034214SMatthew G. Knepley   PetscInt       edgeCut;                  /* The number of edges cut by the partition */
32770034214SMatthew G. Knepley   PetscInt      *assignment, *points;
32870034214SMatthew G. Knepley   PetscMPIInt    commSize, rank, p, v, i;
32970034214SMatthew G. Knepley   PetscErrorCode ierr;
33070034214SMatthew G. Knepley 
33170034214SMatthew G. Knepley   PetscFunctionBegin;
33270034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
33370034214SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &commSize);CHKERRQ(ierr);
33470034214SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
33570034214SMatthew G. Knepley   nparts = commSize;
33670034214SMatthew G. Knepley   options[0] = 0; /* Use all defaults */
33770034214SMatthew G. Knepley   /* Calculate vertex distribution */
33870034214SMatthew G. Knepley   ierr = PetscMalloc4(nparts+1,&vtxdist,nparts*ncon,&tpwgts,ncon,&ubvec,nvtxs,&assignment);CHKERRQ(ierr);
33970034214SMatthew G. Knepley   vtxdist[0] = 0;
34070034214SMatthew G. Knepley   ierr = MPI_Allgather(&nvtxs, 1, MPIU_INT, &vtxdist[1], 1, MPIU_INT, comm);CHKERRQ(ierr);
34170034214SMatthew G. Knepley   for (p = 2; p <= nparts; ++p) {
34270034214SMatthew G. Knepley     vtxdist[p] += vtxdist[p-1];
34370034214SMatthew G. Knepley   }
34470034214SMatthew G. Knepley   /* Calculate weights */
34570034214SMatthew G. Knepley   for (p = 0; p < nparts; ++p) {
34670034214SMatthew G. Knepley     tpwgts[p] = 1.0/nparts;
34770034214SMatthew G. Knepley   }
34870034214SMatthew G. Knepley   ubvec[0] = 1.05;
34970034214SMatthew G. Knepley 
35070034214SMatthew G. Knepley   if (nparts == 1) {
35170034214SMatthew G. Knepley     ierr = PetscMemzero(assignment, nvtxs * sizeof(PetscInt));
35270034214SMatthew G. Knepley   } else {
35370034214SMatthew G. Knepley     if (vtxdist[1] == vtxdist[nparts]) {
35470034214SMatthew G. Knepley       if (!rank) {
35570034214SMatthew G. Knepley         PetscStackPush("METIS_PartGraphKway");
35670034214SMatthew G. Knepley         ierr = METIS_PartGraphKway(&nvtxs, &ncon, xadj, adjncy, vwgt, NULL, adjwgt, &nparts, tpwgts, ubvec, NULL, &edgeCut, assignment);
35770034214SMatthew G. Knepley         PetscStackPop;
35870034214SMatthew G. Knepley         if (ierr != METIS_OK) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in METIS_PartGraphKway()");
35970034214SMatthew G. Knepley       }
36070034214SMatthew G. Knepley     } else {
36170034214SMatthew G. Knepley       PetscStackPush("ParMETIS_V3_PartKway");
36270034214SMatthew G. Knepley       ierr = ParMETIS_V3_PartKway(vtxdist, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgeCut, assignment, &comm);
36370034214SMatthew G. Knepley       PetscStackPop;
36470034214SMatthew G. Knepley       if (ierr != METIS_OK) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in ParMETIS_V3_PartKway()");
36570034214SMatthew G. Knepley     }
36670034214SMatthew G. Knepley   }
36770034214SMatthew G. Knepley   /* Convert to PetscSection+IS */
36870034214SMatthew G. Knepley   ierr = PetscSectionCreate(comm, partSection);CHKERRQ(ierr);
36970034214SMatthew G. Knepley   ierr = PetscSectionSetChart(*partSection, 0, commSize);CHKERRQ(ierr);
37070034214SMatthew G. Knepley   for (v = 0; v < nvtxs; ++v) {
37170034214SMatthew G. Knepley     ierr = PetscSectionAddDof(*partSection, assignment[v], 1);CHKERRQ(ierr);
37270034214SMatthew G. Knepley   }
37370034214SMatthew G. Knepley   ierr = PetscSectionSetUp(*partSection);CHKERRQ(ierr);
37470034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr);
37570034214SMatthew G. Knepley   for (p = 0, i = 0; p < commSize; ++p) {
37670034214SMatthew G. Knepley     for (v = 0; v < nvtxs; ++v) {
37770034214SMatthew G. Knepley       if (assignment[v] == p) points[i++] = v;
37870034214SMatthew G. Knepley     }
37970034214SMatthew G. Knepley   }
38070034214SMatthew G. Knepley   if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs);
38170034214SMatthew G. Knepley   ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
38270034214SMatthew G. Knepley   ierr = PetscFree4(vtxdist,tpwgts,ubvec,assignment);CHKERRQ(ierr);
38370034214SMatthew G. Knepley   PetscFunctionReturn(0);
38470034214SMatthew G. Knepley }
38570034214SMatthew G. Knepley #endif
38670034214SMatthew G. Knepley 
38770034214SMatthew G. Knepley #undef __FUNCT__
38870034214SMatthew G. Knepley #define __FUNCT__ "DMPlexEnlargePartition"
38970034214SMatthew G. Knepley /* Expand the partition by BFS on the adjacency graph */
39070034214SMatthew G. Knepley PetscErrorCode DMPlexEnlargePartition(DM dm, const PetscInt start[], const PetscInt adjacency[], PetscSection origPartSection, IS origPartition, PetscSection *partSection, IS *partition)
39170034214SMatthew G. Knepley {
39270034214SMatthew G. Knepley   PetscHashI      h;
39370034214SMatthew G. Knepley   const PetscInt *points;
39470034214SMatthew G. Knepley   PetscInt      **tmpPoints, *newPoints, totPoints = 0;
39570034214SMatthew G. Knepley   PetscInt        pStart, pEnd, part, q;
39670034214SMatthew G. Knepley   PetscBool       useCone;
39770034214SMatthew G. Knepley   PetscErrorCode  ierr;
39870034214SMatthew G. Knepley 
39970034214SMatthew G. Knepley   PetscFunctionBegin;
40070034214SMatthew G. Knepley   PetscHashICreate(h);
40170034214SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), partSection);CHKERRQ(ierr);
40270034214SMatthew G. Knepley   ierr = PetscSectionGetChart(origPartSection, &pStart, &pEnd);CHKERRQ(ierr);
40370034214SMatthew G. Knepley   ierr = PetscSectionSetChart(*partSection, pStart, pEnd);CHKERRQ(ierr);
40470034214SMatthew G. Knepley   ierr = ISGetIndices(origPartition, &points);CHKERRQ(ierr);
40570034214SMatthew G. Knepley   ierr = PetscMalloc1((pEnd - pStart), &tmpPoints);CHKERRQ(ierr);
40670034214SMatthew G. Knepley   ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
40770034214SMatthew G. Knepley   ierr = DMPlexSetAdjacencyUseCone(dm, PETSC_TRUE);CHKERRQ(ierr);
40870034214SMatthew G. Knepley   for (part = pStart; part < pEnd; ++part) {
40970034214SMatthew G. Knepley     PetscInt *adj = NULL;
41070034214SMatthew G. Knepley     PetscInt  numPoints, nP, numNewPoints, off, p, n = 0;
41170034214SMatthew G. Knepley 
41270034214SMatthew G. Knepley     PetscHashIClear(h);
41370034214SMatthew G. Knepley     ierr = PetscSectionGetDof(origPartSection, part, &numPoints);CHKERRQ(ierr);
41470034214SMatthew G. Knepley     ierr = PetscSectionGetOffset(origPartSection, part, &off);CHKERRQ(ierr);
41570034214SMatthew G. Knepley     /* Add all existing points to h */
41670034214SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
41770034214SMatthew G. Knepley       const PetscInt point = points[off+p];
41870034214SMatthew G. Knepley       PetscHashIAdd(h, point, 1);
41970034214SMatthew G. Knepley     }
42070034214SMatthew G. Knepley     PetscHashISize(h, nP);
42170034214SMatthew 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);
42270034214SMatthew G. Knepley     /* Add all points in next BFS level */
42370034214SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
42470034214SMatthew G. Knepley       const PetscInt point   = points[off+p];
42570034214SMatthew G. Knepley       PetscInt       adjSize = PETSC_DETERMINE, a;
42670034214SMatthew G. Knepley 
42770034214SMatthew G. Knepley       ierr = DMPlexGetAdjacency(dm, point, &adjSize, &adj);CHKERRQ(ierr);
42870034214SMatthew G. Knepley       for (a = 0; a < adjSize; ++a) PetscHashIAdd(h, adj[a], 1);
42970034214SMatthew G. Knepley     }
43070034214SMatthew G. Knepley     PetscHashISize(h, numNewPoints);
43170034214SMatthew G. Knepley     ierr = PetscSectionSetDof(*partSection, part, numNewPoints);CHKERRQ(ierr);
43270034214SMatthew G. Knepley     ierr = PetscMalloc1(numNewPoints, &tmpPoints[part]);CHKERRQ(ierr);
43370034214SMatthew G. Knepley     ierr = PetscHashIGetKeys(h, &n, tmpPoints[part]);CHKERRQ(ierr);
43470034214SMatthew G. Knepley     ierr = PetscFree(adj);CHKERRQ(ierr);
43570034214SMatthew G. Knepley     totPoints += numNewPoints;
43670034214SMatthew G. Knepley   }
43770034214SMatthew G. Knepley   ierr = DMPlexSetAdjacencyUseCone(dm, useCone);CHKERRQ(ierr);
43870034214SMatthew G. Knepley   ierr = ISRestoreIndices(origPartition, &points);CHKERRQ(ierr);
43970034214SMatthew G. Knepley   PetscHashIDestroy(h);
44070034214SMatthew G. Knepley   ierr = PetscSectionSetUp(*partSection);CHKERRQ(ierr);
44170034214SMatthew G. Knepley   ierr = PetscMalloc1(totPoints, &newPoints);CHKERRQ(ierr);
44270034214SMatthew G. Knepley   for (part = pStart, q = 0; part < pEnd; ++part) {
44370034214SMatthew G. Knepley     PetscInt numPoints, p;
44470034214SMatthew G. Knepley 
44570034214SMatthew G. Knepley     ierr = PetscSectionGetDof(*partSection, part, &numPoints);CHKERRQ(ierr);
44670034214SMatthew G. Knepley     for (p = 0; p < numPoints; ++p, ++q) newPoints[q] = tmpPoints[part][p];
44770034214SMatthew G. Knepley     ierr = PetscFree(tmpPoints[part]);CHKERRQ(ierr);
44870034214SMatthew G. Knepley   }
44970034214SMatthew G. Knepley   ierr = PetscFree(tmpPoints);CHKERRQ(ierr);
45070034214SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), totPoints, newPoints, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
45170034214SMatthew G. Knepley   PetscFunctionReturn(0);
45270034214SMatthew G. Knepley }
45370034214SMatthew G. Knepley 
45470034214SMatthew G. Knepley #undef __FUNCT__
45570034214SMatthew G. Knepley #define __FUNCT__ "DMPlexCreatePartition"
45670034214SMatthew G. Knepley /*
45770034214SMatthew G. Knepley   DMPlexCreatePartition - Create a non-overlapping partition of the points at the given height
45870034214SMatthew G. Knepley 
45970034214SMatthew G. Knepley   Collective on DM
46070034214SMatthew G. Knepley 
46170034214SMatthew G. Knepley   Input Parameters:
46270034214SMatthew G. Knepley   + dm - The DM
46370034214SMatthew G. Knepley   . height - The height for points in the partition
46470034214SMatthew G. Knepley   - enlarge - Expand each partition with neighbors
46570034214SMatthew G. Knepley 
46670034214SMatthew G. Knepley   Output Parameters:
46770034214SMatthew G. Knepley   + partSection - The PetscSection giving the division of points by partition
46870034214SMatthew G. Knepley   . partition - The list of points by partition
46970034214SMatthew G. Knepley   . origPartSection - If enlarge is true, the PetscSection giving the division of points before enlarging by partition, otherwise NULL
47070034214SMatthew G. Knepley   - origPartition - If enlarge is true, the list of points before enlarging by partition, otherwise NULL
47170034214SMatthew G. Knepley 
47270034214SMatthew G. Knepley   Level: developer
47370034214SMatthew G. Knepley 
47470034214SMatthew G. Knepley .seealso DMPlexDistribute()
47570034214SMatthew G. Knepley */
47670034214SMatthew G. Knepley PetscErrorCode DMPlexCreatePartition(DM dm, const char name[], PetscInt height, PetscBool enlarge, PetscSection *partSection, IS *partition, PetscSection *origPartSection, IS *origPartition)
47770034214SMatthew G. Knepley {
47870034214SMatthew G. Knepley   char           partname[1024];
47970034214SMatthew G. Knepley   PetscBool      isChaco = PETSC_FALSE, isMetis = PETSC_FALSE, flg;
48070034214SMatthew G. Knepley   PetscMPIInt    size;
48170034214SMatthew G. Knepley   PetscErrorCode ierr;
48270034214SMatthew G. Knepley 
48370034214SMatthew G. Knepley   PetscFunctionBegin;
48470034214SMatthew G. Knepley   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);CHKERRQ(ierr);
48570034214SMatthew G. Knepley 
48670034214SMatthew G. Knepley   *origPartSection = NULL;
48770034214SMatthew G. Knepley   *origPartition   = NULL;
48870034214SMatthew G. Knepley   if (size == 1) {
48970034214SMatthew G. Knepley     PetscInt *points;
49070034214SMatthew G. Knepley     PetscInt  cStart, cEnd, c;
49170034214SMatthew G. Knepley 
49270034214SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
49370034214SMatthew G. Knepley     ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), partSection);CHKERRQ(ierr);
49470034214SMatthew G. Knepley     ierr = PetscSectionSetChart(*partSection, 0, size);CHKERRQ(ierr);
49570034214SMatthew G. Knepley     ierr = PetscSectionSetDof(*partSection, 0, cEnd-cStart);CHKERRQ(ierr);
49670034214SMatthew G. Knepley     ierr = PetscSectionSetUp(*partSection);CHKERRQ(ierr);
49770034214SMatthew G. Knepley     ierr = PetscMalloc1((cEnd - cStart), &points);CHKERRQ(ierr);
49870034214SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) points[c] = c;
49970034214SMatthew G. Knepley     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), cEnd-cStart, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
50070034214SMatthew G. Knepley     PetscFunctionReturn(0);
50170034214SMatthew G. Knepley   }
50270034214SMatthew G. Knepley   ierr = PetscOptionsGetString(((PetscObject) dm)->prefix, "-dm_plex_partitioner", partname, 1024, &flg);CHKERRQ(ierr);
50370034214SMatthew G. Knepley   if (flg) name = partname;
50470034214SMatthew G. Knepley   if (name) {
50570034214SMatthew G. Knepley     ierr = PetscStrcmp(name, "chaco", &isChaco);CHKERRQ(ierr);
50670034214SMatthew G. Knepley     ierr = PetscStrcmp(name, "metis", &isMetis);CHKERRQ(ierr);
50770034214SMatthew G. Knepley   }
50870034214SMatthew G. Knepley   if (height == 0) {
50970034214SMatthew G. Knepley     PetscInt  numVertices;
51070034214SMatthew G. Knepley     PetscInt *start     = NULL;
51170034214SMatthew G. Knepley     PetscInt *adjacency = NULL;
51270034214SMatthew G. Knepley 
51370034214SMatthew G. Knepley     ierr = DMPlexCreateNeighborCSR(dm, 0, &numVertices, &start, &adjacency);CHKERRQ(ierr);
51470034214SMatthew G. Knepley     if (!name || isChaco) {
51570034214SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO)
51670034214SMatthew G. Knepley       ierr = DMPlexPartition_Chaco(dm, numVertices, start, adjacency, partSection, partition);CHKERRQ(ierr);
51770034214SMatthew G. Knepley #else
51870034214SMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-chaco.");
51970034214SMatthew G. Knepley #endif
52070034214SMatthew G. Knepley     } else if (isMetis) {
52170034214SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS)
52270034214SMatthew G. Knepley       ierr = DMPlexPartition_ParMetis(dm, numVertices, start, adjacency, partSection, partition);CHKERRQ(ierr);
52370034214SMatthew G. Knepley #endif
52470034214SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Unknown mesh partitioning package %s", name);
52570034214SMatthew G. Knepley     if (enlarge) {
52670034214SMatthew G. Knepley       *origPartSection = *partSection;
52770034214SMatthew G. Knepley       *origPartition   = *partition;
52870034214SMatthew G. Knepley 
52970034214SMatthew G. Knepley       ierr = DMPlexEnlargePartition(dm, start, adjacency, *origPartSection, *origPartition, partSection, partition);CHKERRQ(ierr);
53070034214SMatthew G. Knepley     }
53170034214SMatthew G. Knepley     ierr = PetscFree(start);CHKERRQ(ierr);
53270034214SMatthew G. Knepley     ierr = PetscFree(adjacency);CHKERRQ(ierr);
53370034214SMatthew G. Knepley # if 0
53470034214SMatthew G. Knepley   } else if (height == 1) {
53570034214SMatthew G. Knepley     /* Build the dual graph for faces and partition the hypergraph */
53670034214SMatthew G. Knepley     PetscInt numEdges;
53770034214SMatthew G. Knepley 
53870034214SMatthew G. Knepley     buildFaceCSRV(mesh, mesh->getFactory()->getNumbering(mesh, mesh->depth()-1), &numEdges, &start, &adjacency, GraphPartitioner::zeroBase());
53970034214SMatthew G. Knepley     GraphPartitioner().partition(numEdges, start, adjacency, partition, manager);
54070034214SMatthew G. Knepley     destroyCSR(numEdges, start, adjacency);
54170034214SMatthew G. Knepley #endif
54270034214SMatthew G. Knepley   } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid partition height %D", height);
54370034214SMatthew G. Knepley   PetscFunctionReturn(0);
54470034214SMatthew G. Knepley }
54570034214SMatthew G. Knepley 
54670034214SMatthew G. Knepley #undef __FUNCT__
547*9ffc88e4SToby Isaac #define __FUNCT__ "DMPlexMarkTreeClosure"
548*9ffc88e4SToby Isaac static PetscErrorCode DMPlexMarkTreeClosure(DM dm, PetscSegBuffer segpart, PetscBT bt, PetscInt point, PetscInt *partSize)
549*9ffc88e4SToby Isaac {
550*9ffc88e4SToby Isaac   PetscInt       parent, closureSize, *closure = NULL, i, pStart, pEnd;
551*9ffc88e4SToby Isaac   PetscErrorCode ierr;
552*9ffc88e4SToby Isaac 
553*9ffc88e4SToby Isaac   PetscFunctionBegin;
554*9ffc88e4SToby Isaac   ierr = DMPlexGetTreeParent(dm,point,&parent,NULL);CHKERRQ(ierr);
555*9ffc88e4SToby Isaac   if (parent == point) PetscFunctionReturn(0);
556*9ffc88e4SToby Isaac   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
557*9ffc88e4SToby Isaac   ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
558*9ffc88e4SToby Isaac   for (i = 0; i < closureSize; i++) {
559*9ffc88e4SToby Isaac     PetscInt cpoint = closure[2*i];
560*9ffc88e4SToby Isaac 
561*9ffc88e4SToby Isaac     if (!PetscBTLookupSet(bt,cpoint-pStart)) {
562*9ffc88e4SToby Isaac       PetscInt *PETSC_RESTRICT pt;
563*9ffc88e4SToby Isaac       (*partSize)++;
564*9ffc88e4SToby Isaac       ierr = PetscSegBufferGetInts(segpart,1,&pt);CHKERRQ(ierr);
565*9ffc88e4SToby Isaac       *pt = cpoint;
566*9ffc88e4SToby Isaac     }
567*9ffc88e4SToby Isaac     ierr = DMPlexMarkTreeClosure(dm,segpart,bt,cpoint,partSize);CHKERRQ(ierr);
568*9ffc88e4SToby Isaac   }
569*9ffc88e4SToby Isaac   ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
570*9ffc88e4SToby Isaac   PetscFunctionReturn(0);
571*9ffc88e4SToby Isaac }
572*9ffc88e4SToby Isaac 
573*9ffc88e4SToby Isaac #undef __FUNCT__
57470034214SMatthew G. Knepley #define __FUNCT__ "DMPlexCreatePartitionClosure"
57570034214SMatthew G. Knepley PetscErrorCode DMPlexCreatePartitionClosure(DM dm, PetscSection pointSection, IS pointPartition, PetscSection *section, IS *partition)
57670034214SMatthew G. Knepley {
57770034214SMatthew G. Knepley   /* const PetscInt  height = 0; */
57870034214SMatthew G. Knepley   const PetscInt *partArray;
57970034214SMatthew G. Knepley   PetscInt       *allPoints, *packPoints;
58070034214SMatthew G. Knepley   PetscInt        rStart, rEnd, rank, pStart, pEnd, newSize;
58170034214SMatthew G. Knepley   PetscErrorCode  ierr;
58270034214SMatthew G. Knepley   PetscBT         bt;
58370034214SMatthew G. Knepley   PetscSegBuffer  segpack,segpart;
58470034214SMatthew G. Knepley 
58570034214SMatthew G. Knepley   PetscFunctionBegin;
58670034214SMatthew G. Knepley   ierr = PetscSectionGetChart(pointSection, &rStart, &rEnd);CHKERRQ(ierr);
58770034214SMatthew G. Knepley   ierr = ISGetIndices(pointPartition, &partArray);CHKERRQ(ierr);
58870034214SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), section);CHKERRQ(ierr);
58970034214SMatthew G. Knepley   ierr = PetscSectionSetChart(*section, rStart, rEnd);CHKERRQ(ierr);
59070034214SMatthew G. Knepley   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
59170034214SMatthew G. Knepley   ierr = PetscBTCreate(pEnd-pStart,&bt);CHKERRQ(ierr);
59270034214SMatthew G. Knepley   ierr = PetscSegBufferCreate(sizeof(PetscInt),1000,&segpack);CHKERRQ(ierr);
59370034214SMatthew G. Knepley   ierr = PetscSegBufferCreate(sizeof(PetscInt),1000,&segpart);CHKERRQ(ierr);
59470034214SMatthew G. Knepley   for (rank = rStart; rank < rEnd; ++rank) {
59570034214SMatthew G. Knepley     PetscInt partSize = 0, numPoints, offset, p, *PETSC_RESTRICT placePoints;
59670034214SMatthew G. Knepley 
59770034214SMatthew G. Knepley     ierr = PetscSectionGetDof(pointSection, rank, &numPoints);CHKERRQ(ierr);
59870034214SMatthew G. Knepley     ierr = PetscSectionGetOffset(pointSection, rank, &offset);CHKERRQ(ierr);
59970034214SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
60070034214SMatthew G. Knepley       PetscInt  point   = partArray[offset+p], closureSize, c;
60170034214SMatthew G. Knepley       PetscInt *closure = NULL;
60270034214SMatthew G. Knepley 
60370034214SMatthew G. Knepley       /* TODO Include support for height > 0 case */
60470034214SMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
60570034214SMatthew G. Knepley       for (c=0; c<closureSize; c++) {
60670034214SMatthew G. Knepley         PetscInt cpoint = closure[c*2];
60770034214SMatthew G. Knepley         if (!PetscBTLookupSet(bt,cpoint-pStart)) {
60870034214SMatthew G. Knepley           PetscInt *PETSC_RESTRICT pt;
60970034214SMatthew G. Knepley           partSize++;
61070034214SMatthew G. Knepley           ierr = PetscSegBufferGetInts(segpart,1,&pt);CHKERRQ(ierr);
61170034214SMatthew G. Knepley           *pt = cpoint;
61270034214SMatthew G. Knepley         }
613*9ffc88e4SToby Isaac         ierr = DMPlexMarkTreeClosure(dm,segpart,bt,cpoint,&partSize);CHKERRQ(ierr);
61470034214SMatthew G. Knepley       }
61570034214SMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
61670034214SMatthew G. Knepley     }
61770034214SMatthew G. Knepley     ierr = PetscSectionSetDof(*section, rank, partSize);CHKERRQ(ierr);
61870034214SMatthew G. Knepley     ierr = PetscSegBufferGetInts(segpack,partSize,&placePoints);CHKERRQ(ierr);
61970034214SMatthew G. Knepley     ierr = PetscSegBufferExtractTo(segpart,placePoints);CHKERRQ(ierr);
62070034214SMatthew G. Knepley     ierr = PetscSortInt(partSize,placePoints);CHKERRQ(ierr);
62170034214SMatthew G. Knepley     for (p=0; p<partSize; p++) {ierr = PetscBTClear(bt,placePoints[p]-pStart);CHKERRQ(ierr);}
62270034214SMatthew G. Knepley   }
62370034214SMatthew G. Knepley   ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
62470034214SMatthew G. Knepley   ierr = PetscSegBufferDestroy(&segpart);CHKERRQ(ierr);
62570034214SMatthew G. Knepley 
62670034214SMatthew G. Knepley   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
62770034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(*section, &newSize);CHKERRQ(ierr);
62870034214SMatthew G. Knepley   ierr = PetscMalloc1(newSize, &allPoints);CHKERRQ(ierr);
62970034214SMatthew G. Knepley 
63070034214SMatthew G. Knepley   ierr = PetscSegBufferExtractInPlace(segpack,&packPoints);CHKERRQ(ierr);
63170034214SMatthew G. Knepley   for (rank = rStart; rank < rEnd; ++rank) {
63270034214SMatthew G. Knepley     PetscInt numPoints, offset;
63370034214SMatthew G. Knepley 
63470034214SMatthew G. Knepley     ierr = PetscSectionGetDof(*section, rank, &numPoints);CHKERRQ(ierr);
63570034214SMatthew G. Knepley     ierr = PetscSectionGetOffset(*section, rank, &offset);CHKERRQ(ierr);
63670034214SMatthew G. Knepley     ierr = PetscMemcpy(&allPoints[offset], packPoints, numPoints * sizeof(PetscInt));CHKERRQ(ierr);
63770034214SMatthew G. Knepley     packPoints += numPoints;
63870034214SMatthew G. Knepley   }
63970034214SMatthew G. Knepley 
64070034214SMatthew G. Knepley   ierr = PetscSegBufferDestroy(&segpack);CHKERRQ(ierr);
64170034214SMatthew G. Knepley   ierr = ISRestoreIndices(pointPartition, &partArray);CHKERRQ(ierr);
64270034214SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), newSize, allPoints, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
64370034214SMatthew G. Knepley   PetscFunctionReturn(0);
64470034214SMatthew G. Knepley }
645