xref: /petsc/src/dm/impls/plex/plexpartition.c (revision 70034214b4f51eb44fb6e1753b0a1c1f20e8240c)
1*70034214SMatthew G. Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2*70034214SMatthew G. Knepley 
3*70034214SMatthew G. Knepley #undef __FUNCT__
4*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexCreateNeighborCSR"
5*70034214SMatthew G. Knepley PetscErrorCode DMPlexCreateNeighborCSR(DM dm, PetscInt cellHeight, PetscInt *numVertices, PetscInt **offsets, PetscInt **adjacency)
6*70034214SMatthew G. Knepley {
7*70034214SMatthew G. Knepley   const PetscInt maxFaceCases = 30;
8*70034214SMatthew G. Knepley   PetscInt       numFaceCases = 0;
9*70034214SMatthew G. Knepley   PetscInt       numFaceVertices[30]; /* maxFaceCases, C89 sucks sucks sucks */
10*70034214SMatthew G. Knepley   PetscInt      *off, *adj;
11*70034214SMatthew G. Knepley   PetscInt      *neighborCells = NULL;
12*70034214SMatthew G. Knepley   PetscInt       dim, cellDim, depth = 0, faceDepth, cStart, cEnd, c, numCells, cell;
13*70034214SMatthew G. Knepley   PetscErrorCode ierr;
14*70034214SMatthew G. Knepley 
15*70034214SMatthew G. Knepley   PetscFunctionBegin;
16*70034214SMatthew G. Knepley   /* For parallel partitioning, I think you have to communicate supports */
17*70034214SMatthew G. Knepley   ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr);
18*70034214SMatthew G. Knepley   cellDim = dim - cellHeight;
19*70034214SMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
20*70034214SMatthew G. Knepley   ierr = DMPlexGetHeightStratum(dm, cellHeight, &cStart, &cEnd);CHKERRQ(ierr);
21*70034214SMatthew G. Knepley   if (cEnd - cStart == 0) {
22*70034214SMatthew G. Knepley     if (numVertices) *numVertices = 0;
23*70034214SMatthew G. Knepley     if (offsets)   *offsets   = NULL;
24*70034214SMatthew G. Knepley     if (adjacency) *adjacency = NULL;
25*70034214SMatthew G. Knepley     PetscFunctionReturn(0);
26*70034214SMatthew G. Knepley   }
27*70034214SMatthew G. Knepley   numCells  = cEnd - cStart;
28*70034214SMatthew G. Knepley   faceDepth = depth - cellHeight;
29*70034214SMatthew G. Knepley   if (dim == depth) {
30*70034214SMatthew G. Knepley     PetscInt f, fStart, fEnd;
31*70034214SMatthew G. Knepley 
32*70034214SMatthew G. Knepley     ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr);
33*70034214SMatthew G. Knepley     /* Count neighboring cells */
34*70034214SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, cellHeight+1, &fStart, &fEnd);CHKERRQ(ierr);
35*70034214SMatthew G. Knepley     for (f = fStart; f < fEnd; ++f) {
36*70034214SMatthew G. Knepley       const PetscInt *support;
37*70034214SMatthew G. Knepley       PetscInt        supportSize;
38*70034214SMatthew G. Knepley       ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr);
39*70034214SMatthew G. Knepley       ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr);
40*70034214SMatthew G. Knepley       if (supportSize == 2) {
41*70034214SMatthew G. Knepley         ++off[support[0]-cStart+1];
42*70034214SMatthew G. Knepley         ++off[support[1]-cStart+1];
43*70034214SMatthew G. Knepley       }
44*70034214SMatthew G. Knepley     }
45*70034214SMatthew G. Knepley     /* Prefix sum */
46*70034214SMatthew G. Knepley     for (c = 1; c <= numCells; ++c) off[c] += off[c-1];
47*70034214SMatthew G. Knepley     if (adjacency) {
48*70034214SMatthew G. Knepley       PetscInt *tmp;
49*70034214SMatthew G. Knepley 
50*70034214SMatthew G. Knepley       ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr);
51*70034214SMatthew G. Knepley       ierr = PetscMalloc1((numCells+1), &tmp);CHKERRQ(ierr);
52*70034214SMatthew G. Knepley       ierr = PetscMemcpy(tmp, off, (numCells+1) * sizeof(PetscInt));CHKERRQ(ierr);
53*70034214SMatthew G. Knepley       /* Get neighboring cells */
54*70034214SMatthew G. Knepley       for (f = fStart; f < fEnd; ++f) {
55*70034214SMatthew G. Knepley         const PetscInt *support;
56*70034214SMatthew G. Knepley         PetscInt        supportSize;
57*70034214SMatthew G. Knepley         ierr = DMPlexGetSupportSize(dm, f, &supportSize);CHKERRQ(ierr);
58*70034214SMatthew G. Knepley         ierr = DMPlexGetSupport(dm, f, &support);CHKERRQ(ierr);
59*70034214SMatthew G. Knepley         if (supportSize == 2) {
60*70034214SMatthew G. Knepley           adj[tmp[support[0]-cStart]++] = support[1];
61*70034214SMatthew G. Knepley           adj[tmp[support[1]-cStart]++] = support[0];
62*70034214SMatthew G. Knepley         }
63*70034214SMatthew G. Knepley       }
64*70034214SMatthew G. Knepley #if defined(PETSC_USE_DEBUG)
65*70034214SMatthew 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);
66*70034214SMatthew G. Knepley #endif
67*70034214SMatthew G. Knepley       ierr = PetscFree(tmp);CHKERRQ(ierr);
68*70034214SMatthew G. Knepley     }
69*70034214SMatthew G. Knepley     if (numVertices) *numVertices = numCells;
70*70034214SMatthew G. Knepley     if (offsets)   *offsets   = off;
71*70034214SMatthew G. Knepley     if (adjacency) *adjacency = adj;
72*70034214SMatthew G. Knepley     PetscFunctionReturn(0);
73*70034214SMatthew G. Knepley   }
74*70034214SMatthew G. Knepley   /* Setup face recognition */
75*70034214SMatthew G. Knepley   if (faceDepth == 1) {
76*70034214SMatthew 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 */
77*70034214SMatthew G. Knepley 
78*70034214SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) {
79*70034214SMatthew G. Knepley       PetscInt corners;
80*70034214SMatthew G. Knepley 
81*70034214SMatthew G. Knepley       ierr = DMPlexGetConeSize(dm, c, &corners);CHKERRQ(ierr);
82*70034214SMatthew G. Knepley       if (!cornersSeen[corners]) {
83*70034214SMatthew G. Knepley         PetscInt nFV;
84*70034214SMatthew G. Knepley 
85*70034214SMatthew G. Knepley         if (numFaceCases >= maxFaceCases) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_PLIB, "Exceeded maximum number of face recognition cases");
86*70034214SMatthew G. Knepley         cornersSeen[corners] = 1;
87*70034214SMatthew G. Knepley 
88*70034214SMatthew G. Knepley         ierr = DMPlexGetNumFaceVertices(dm, cellDim, corners, &nFV);CHKERRQ(ierr);
89*70034214SMatthew G. Knepley 
90*70034214SMatthew G. Knepley         numFaceVertices[numFaceCases++] = nFV;
91*70034214SMatthew G. Knepley       }
92*70034214SMatthew G. Knepley     }
93*70034214SMatthew G. Knepley   }
94*70034214SMatthew G. Knepley   ierr = PetscCalloc1(numCells+1, &off);CHKERRQ(ierr);
95*70034214SMatthew G. Knepley   /* Count neighboring cells */
96*70034214SMatthew G. Knepley   for (cell = cStart; cell < cEnd; ++cell) {
97*70034214SMatthew G. Knepley     PetscInt numNeighbors = PETSC_DETERMINE, n;
98*70034214SMatthew G. Knepley 
99*70034214SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, &numNeighbors, neighborCells);CHKERRQ(ierr);
100*70034214SMatthew G. Knepley     /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */
101*70034214SMatthew G. Knepley     for (n = 0; n < numNeighbors; ++n) {
102*70034214SMatthew G. Knepley       PetscInt        cellPair[2];
103*70034214SMatthew G. Knepley       PetscBool       found    = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE;
104*70034214SMatthew G. Knepley       PetscInt        meetSize = 0;
105*70034214SMatthew G. Knepley       const PetscInt *meet    = NULL;
106*70034214SMatthew G. Knepley 
107*70034214SMatthew G. Knepley       cellPair[0] = cell; cellPair[1] = neighborCells[n];
108*70034214SMatthew G. Knepley       if (cellPair[0] == cellPair[1]) continue;
109*70034214SMatthew G. Knepley       if (!found) {
110*70034214SMatthew G. Knepley         ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
111*70034214SMatthew G. Knepley         if (meetSize) {
112*70034214SMatthew G. Knepley           PetscInt f;
113*70034214SMatthew G. Knepley 
114*70034214SMatthew G. Knepley           for (f = 0; f < numFaceCases; ++f) {
115*70034214SMatthew G. Knepley             if (numFaceVertices[f] == meetSize) {
116*70034214SMatthew G. Knepley               found = PETSC_TRUE;
117*70034214SMatthew G. Knepley               break;
118*70034214SMatthew G. Knepley             }
119*70034214SMatthew G. Knepley           }
120*70034214SMatthew G. Knepley         }
121*70034214SMatthew G. Knepley         ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
122*70034214SMatthew G. Knepley       }
123*70034214SMatthew G. Knepley       if (found) ++off[cell-cStart+1];
124*70034214SMatthew G. Knepley     }
125*70034214SMatthew G. Knepley   }
126*70034214SMatthew G. Knepley   /* Prefix sum */
127*70034214SMatthew G. Knepley   for (cell = 1; cell <= numCells; ++cell) off[cell] += off[cell-1];
128*70034214SMatthew G. Knepley 
129*70034214SMatthew G. Knepley   if (adjacency) {
130*70034214SMatthew G. Knepley     ierr = PetscMalloc1(off[numCells], &adj);CHKERRQ(ierr);
131*70034214SMatthew G. Knepley     /* Get neighboring cells */
132*70034214SMatthew G. Knepley     for (cell = cStart; cell < cEnd; ++cell) {
133*70034214SMatthew G. Knepley       PetscInt numNeighbors = PETSC_DETERMINE, n;
134*70034214SMatthew G. Knepley       PetscInt cellOffset   = 0;
135*70034214SMatthew G. Knepley 
136*70034214SMatthew G. Knepley       ierr = DMPlexGetAdjacency_Internal(dm, cell, PETSC_TRUE, PETSC_FALSE, &numNeighbors, neighborCells);CHKERRQ(ierr);
137*70034214SMatthew G. Knepley       /* Get meet with each cell, and check with recognizer (could optimize to check each pair only once) */
138*70034214SMatthew G. Knepley       for (n = 0; n < numNeighbors; ++n) {
139*70034214SMatthew G. Knepley         PetscInt        cellPair[2];
140*70034214SMatthew G. Knepley         PetscBool       found    = faceDepth > 1 ? PETSC_TRUE : PETSC_FALSE;
141*70034214SMatthew G. Knepley         PetscInt        meetSize = 0;
142*70034214SMatthew G. Knepley         const PetscInt *meet    = NULL;
143*70034214SMatthew G. Knepley 
144*70034214SMatthew G. Knepley         cellPair[0] = cell; cellPair[1] = neighborCells[n];
145*70034214SMatthew G. Knepley         if (cellPair[0] == cellPair[1]) continue;
146*70034214SMatthew G. Knepley         if (!found) {
147*70034214SMatthew G. Knepley           ierr = DMPlexGetMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
148*70034214SMatthew G. Knepley           if (meetSize) {
149*70034214SMatthew G. Knepley             PetscInt f;
150*70034214SMatthew G. Knepley 
151*70034214SMatthew G. Knepley             for (f = 0; f < numFaceCases; ++f) {
152*70034214SMatthew G. Knepley               if (numFaceVertices[f] == meetSize) {
153*70034214SMatthew G. Knepley                 found = PETSC_TRUE;
154*70034214SMatthew G. Knepley                 break;
155*70034214SMatthew G. Knepley               }
156*70034214SMatthew G. Knepley             }
157*70034214SMatthew G. Knepley           }
158*70034214SMatthew G. Knepley           ierr = DMPlexRestoreMeet(dm, 2, cellPair, &meetSize, &meet);CHKERRQ(ierr);
159*70034214SMatthew G. Knepley         }
160*70034214SMatthew G. Knepley         if (found) {
161*70034214SMatthew G. Knepley           adj[off[cell-cStart]+cellOffset] = neighborCells[n];
162*70034214SMatthew G. Knepley           ++cellOffset;
163*70034214SMatthew G. Knepley         }
164*70034214SMatthew G. Knepley       }
165*70034214SMatthew G. Knepley     }
166*70034214SMatthew G. Knepley   }
167*70034214SMatthew G. Knepley   ierr = PetscFree(neighborCells);CHKERRQ(ierr);
168*70034214SMatthew G. Knepley   if (numVertices) *numVertices = numCells;
169*70034214SMatthew G. Knepley   if (offsets)   *offsets   = off;
170*70034214SMatthew G. Knepley   if (adjacency) *adjacency = adj;
171*70034214SMatthew G. Knepley   PetscFunctionReturn(0);
172*70034214SMatthew G. Knepley }
173*70034214SMatthew G. Knepley 
174*70034214SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO)
175*70034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
176*70034214SMatthew G. Knepley #include <unistd.h>
177*70034214SMatthew G. Knepley #endif
178*70034214SMatthew G. Knepley /* Chaco does not have an include file */
179*70034214SMatthew G. Knepley PETSC_EXTERN int interface(int nvtxs, int *start, int *adjacency, int *vwgts,
180*70034214SMatthew G. Knepley                        float *ewgts, float *x, float *y, float *z, char *outassignname,
181*70034214SMatthew G. Knepley                        char *outfilename, short *assignment, int architecture, int ndims_tot,
182*70034214SMatthew G. Knepley                        int mesh_dims[3], double *goal, int global_method, int local_method,
183*70034214SMatthew G. Knepley                        int rqi_flag, int vmax, int ndims, double eigtol, long seed);
184*70034214SMatthew G. Knepley 
185*70034214SMatthew G. Knepley extern int FREE_GRAPH;
186*70034214SMatthew G. Knepley 
187*70034214SMatthew G. Knepley #undef __FUNCT__
188*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexPartition_Chaco"
189*70034214SMatthew G. Knepley PetscErrorCode DMPlexPartition_Chaco(DM dm, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection *partSection, IS *partition)
190*70034214SMatthew G. Knepley {
191*70034214SMatthew G. Knepley   enum {DEFAULT_METHOD = 1, INERTIAL_METHOD = 3};
192*70034214SMatthew G. Knepley   MPI_Comm       comm;
193*70034214SMatthew G. Knepley   int            nvtxs          = numVertices; /* number of vertices in full graph */
194*70034214SMatthew G. Knepley   int           *vwgts          = NULL;   /* weights for all vertices */
195*70034214SMatthew G. Knepley   float         *ewgts          = NULL;   /* weights for all edges */
196*70034214SMatthew G. Knepley   float         *x              = NULL, *y = NULL, *z = NULL; /* coordinates for inertial method */
197*70034214SMatthew G. Knepley   char          *outassignname  = NULL;   /*  name of assignment output file */
198*70034214SMatthew G. Knepley   char          *outfilename    = NULL;   /* output file name */
199*70034214SMatthew G. Knepley   int            architecture   = 1;      /* 0 => hypercube, d => d-dimensional mesh */
200*70034214SMatthew G. Knepley   int            ndims_tot      = 0;      /* total number of cube dimensions to divide */
201*70034214SMatthew G. Knepley   int            mesh_dims[3];            /* dimensions of mesh of processors */
202*70034214SMatthew G. Knepley   double        *goal          = NULL;    /* desired set sizes for each set */
203*70034214SMatthew G. Knepley   int            global_method = 1;       /* global partitioning algorithm */
204*70034214SMatthew G. Knepley   int            local_method  = 1;       /* local partitioning algorithm */
205*70034214SMatthew G. Knepley   int            rqi_flag      = 0;       /* should I use RQI/Symmlq eigensolver? */
206*70034214SMatthew G. Knepley   int            vmax          = 200;     /* how many vertices to coarsen down to? */
207*70034214SMatthew G. Knepley   int            ndims         = 1;       /* number of eigenvectors (2^d sets) */
208*70034214SMatthew G. Knepley   double         eigtol        = 0.001;   /* tolerance on eigenvectors */
209*70034214SMatthew G. Knepley   long           seed          = 123636512; /* for random graph mutations */
210*70034214SMatthew G. Knepley   short int     *assignment;              /* Output partition */
211*70034214SMatthew G. Knepley   int            fd_stdout, fd_pipe[2];
212*70034214SMatthew G. Knepley   PetscInt      *points;
213*70034214SMatthew G. Knepley   PetscMPIInt    commSize;
214*70034214SMatthew G. Knepley   int            i, v, p;
215*70034214SMatthew G. Knepley   PetscErrorCode ierr;
216*70034214SMatthew G. Knepley 
217*70034214SMatthew G. Knepley   PetscFunctionBegin;
218*70034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
219*70034214SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &commSize);CHKERRQ(ierr);
220*70034214SMatthew G. Knepley   if (!numVertices) {
221*70034214SMatthew G. Knepley     ierr = PetscSectionCreate(comm, partSection);CHKERRQ(ierr);
222*70034214SMatthew G. Knepley     ierr = PetscSectionSetChart(*partSection, 0, commSize);CHKERRQ(ierr);
223*70034214SMatthew G. Knepley     ierr = PetscSectionSetUp(*partSection);CHKERRQ(ierr);
224*70034214SMatthew G. Knepley     ierr = ISCreateGeneral(comm, 0, NULL, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
225*70034214SMatthew G. Knepley     PetscFunctionReturn(0);
226*70034214SMatthew G. Knepley   }
227*70034214SMatthew G. Knepley   FREE_GRAPH = 0;                         /* Do not let Chaco free my memory */
228*70034214SMatthew G. Knepley   for (i = 0; i < start[numVertices]; ++i) ++adjacency[i];
229*70034214SMatthew G. Knepley 
230*70034214SMatthew G. Knepley   if (global_method == INERTIAL_METHOD) {
231*70034214SMatthew G. Knepley     /* manager.createCellCoordinates(nvtxs, &x, &y, &z); */
232*70034214SMatthew G. Knepley     SETERRQ(comm, PETSC_ERR_SUP, "Inertial partitioning not yet supported");
233*70034214SMatthew G. Knepley   }
234*70034214SMatthew G. Knepley   mesh_dims[0] = commSize;
235*70034214SMatthew G. Knepley   mesh_dims[1] = 1;
236*70034214SMatthew G. Knepley   mesh_dims[2] = 1;
237*70034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &assignment);CHKERRQ(ierr);
238*70034214SMatthew G. Knepley   /* Chaco outputs to stdout. We redirect this to a buffer. */
239*70034214SMatthew G. Knepley   /* TODO: check error codes for UNIX calls */
240*70034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
241*70034214SMatthew G. Knepley   {
242*70034214SMatthew G. Knepley     int piperet;
243*70034214SMatthew G. Knepley     piperet = pipe(fd_pipe);
244*70034214SMatthew G. Knepley     if (piperet) SETERRQ(comm,PETSC_ERR_SYS,"Could not create pipe");
245*70034214SMatthew G. Knepley     fd_stdout = dup(1);
246*70034214SMatthew G. Knepley     close(1);
247*70034214SMatthew G. Knepley     dup2(fd_pipe[1], 1);
248*70034214SMatthew G. Knepley   }
249*70034214SMatthew G. Knepley #endif
250*70034214SMatthew G. Knepley   ierr = interface(nvtxs, (int*) start, (int*) adjacency, vwgts, ewgts, x, y, z, outassignname, outfilename,
251*70034214SMatthew G. Knepley                    assignment, architecture, ndims_tot, mesh_dims, goal, global_method, local_method, rqi_flag,
252*70034214SMatthew G. Knepley                    vmax, ndims, eigtol, seed);
253*70034214SMatthew G. Knepley #if defined(PETSC_HAVE_UNISTD_H)
254*70034214SMatthew G. Knepley   {
255*70034214SMatthew G. Knepley     char msgLog[10000];
256*70034214SMatthew G. Knepley     int  count;
257*70034214SMatthew G. Knepley 
258*70034214SMatthew G. Knepley     fflush(stdout);
259*70034214SMatthew G. Knepley     count = read(fd_pipe[0], msgLog, (10000-1)*sizeof(char));
260*70034214SMatthew G. Knepley     if (count < 0) count = 0;
261*70034214SMatthew G. Knepley     msgLog[count] = 0;
262*70034214SMatthew G. Knepley     close(1);
263*70034214SMatthew G. Knepley     dup2(fd_stdout, 1);
264*70034214SMatthew G. Knepley     close(fd_stdout);
265*70034214SMatthew G. Knepley     close(fd_pipe[0]);
266*70034214SMatthew G. Knepley     close(fd_pipe[1]);
267*70034214SMatthew G. Knepley     if (ierr) SETERRQ1(comm, PETSC_ERR_LIB, "Error in Chaco library: %s", msgLog);
268*70034214SMatthew G. Knepley   }
269*70034214SMatthew G. Knepley #endif
270*70034214SMatthew G. Knepley   /* Convert to PetscSection+IS */
271*70034214SMatthew G. Knepley   ierr = PetscSectionCreate(comm, partSection);CHKERRQ(ierr);
272*70034214SMatthew G. Knepley   ierr = PetscSectionSetChart(*partSection, 0, commSize);CHKERRQ(ierr);
273*70034214SMatthew G. Knepley   for (v = 0; v < nvtxs; ++v) {
274*70034214SMatthew G. Knepley     ierr = PetscSectionAddDof(*partSection, assignment[v], 1);CHKERRQ(ierr);
275*70034214SMatthew G. Knepley   }
276*70034214SMatthew G. Knepley   ierr = PetscSectionSetUp(*partSection);CHKERRQ(ierr);
277*70034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr);
278*70034214SMatthew G. Knepley   for (p = 0, i = 0; p < commSize; ++p) {
279*70034214SMatthew G. Knepley     for (v = 0; v < nvtxs; ++v) {
280*70034214SMatthew G. Knepley       if (assignment[v] == p) points[i++] = v;
281*70034214SMatthew G. Knepley     }
282*70034214SMatthew G. Knepley   }
283*70034214SMatthew G. Knepley   if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs);
284*70034214SMatthew G. Knepley   ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
285*70034214SMatthew G. Knepley   if (global_method == INERTIAL_METHOD) {
286*70034214SMatthew G. Knepley     /* manager.destroyCellCoordinates(nvtxs, &x, &y, &z); */
287*70034214SMatthew G. Knepley   }
288*70034214SMatthew G. Knepley   ierr = PetscFree(assignment);CHKERRQ(ierr);
289*70034214SMatthew G. Knepley   for (i = 0; i < start[numVertices]; ++i) --adjacency[i];
290*70034214SMatthew G. Knepley   PetscFunctionReturn(0);
291*70034214SMatthew G. Knepley }
292*70034214SMatthew G. Knepley #endif
293*70034214SMatthew G. Knepley 
294*70034214SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS)
295*70034214SMatthew G. Knepley #include <parmetis.h>
296*70034214SMatthew G. Knepley 
297*70034214SMatthew G. Knepley #undef __FUNCT__
298*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexPartition_ParMetis"
299*70034214SMatthew G. Knepley PetscErrorCode DMPlexPartition_ParMetis(DM dm, PetscInt numVertices, PetscInt start[], PetscInt adjacency[], PetscSection *partSection, IS *partition)
300*70034214SMatthew G. Knepley {
301*70034214SMatthew G. Knepley   MPI_Comm       comm;
302*70034214SMatthew G. Knepley   PetscInt       nvtxs      = numVertices; /* The number of vertices in full graph */
303*70034214SMatthew G. Knepley   PetscInt      *vtxdist;                  /* Distribution of vertices across processes */
304*70034214SMatthew G. Knepley   PetscInt      *xadj       = start;       /* Start of edge list for each vertex */
305*70034214SMatthew G. Knepley   PetscInt      *adjncy     = adjacency;   /* Edge lists for all vertices */
306*70034214SMatthew G. Knepley   PetscInt      *vwgt       = NULL;        /* Vertex weights */
307*70034214SMatthew G. Knepley   PetscInt      *adjwgt     = NULL;        /* Edge weights */
308*70034214SMatthew G. Knepley   PetscInt       wgtflag    = 0;           /* Indicates which weights are present */
309*70034214SMatthew G. Knepley   PetscInt       numflag    = 0;           /* Indicates initial offset (0 or 1) */
310*70034214SMatthew G. Knepley   PetscInt       ncon       = 1;           /* The number of weights per vertex */
311*70034214SMatthew G. Knepley   PetscInt       nparts;                   /* The number of partitions */
312*70034214SMatthew G. Knepley   PetscReal     *tpwgts;                   /* The fraction of vertex weights assigned to each partition */
313*70034214SMatthew G. Knepley   PetscReal     *ubvec;                    /* The balance intolerance for vertex weights */
314*70034214SMatthew G. Knepley   PetscInt       options[5];               /* Options */
315*70034214SMatthew G. Knepley   /* Outputs */
316*70034214SMatthew G. Knepley   PetscInt       edgeCut;                  /* The number of edges cut by the partition */
317*70034214SMatthew G. Knepley   PetscInt      *assignment, *points;
318*70034214SMatthew G. Knepley   PetscMPIInt    commSize, rank, p, v, i;
319*70034214SMatthew G. Knepley   PetscErrorCode ierr;
320*70034214SMatthew G. Knepley 
321*70034214SMatthew G. Knepley   PetscFunctionBegin;
322*70034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
323*70034214SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &commSize);CHKERRQ(ierr);
324*70034214SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
325*70034214SMatthew G. Knepley   nparts = commSize;
326*70034214SMatthew G. Knepley   options[0] = 0; /* Use all defaults */
327*70034214SMatthew G. Knepley   /* Calculate vertex distribution */
328*70034214SMatthew G. Knepley   ierr = PetscMalloc4(nparts+1,&vtxdist,nparts*ncon,&tpwgts,ncon,&ubvec,nvtxs,&assignment);CHKERRQ(ierr);
329*70034214SMatthew G. Knepley   vtxdist[0] = 0;
330*70034214SMatthew G. Knepley   ierr = MPI_Allgather(&nvtxs, 1, MPIU_INT, &vtxdist[1], 1, MPIU_INT, comm);CHKERRQ(ierr);
331*70034214SMatthew G. Knepley   for (p = 2; p <= nparts; ++p) {
332*70034214SMatthew G. Knepley     vtxdist[p] += vtxdist[p-1];
333*70034214SMatthew G. Knepley   }
334*70034214SMatthew G. Knepley   /* Calculate weights */
335*70034214SMatthew G. Knepley   for (p = 0; p < nparts; ++p) {
336*70034214SMatthew G. Knepley     tpwgts[p] = 1.0/nparts;
337*70034214SMatthew G. Knepley   }
338*70034214SMatthew G. Knepley   ubvec[0] = 1.05;
339*70034214SMatthew G. Knepley 
340*70034214SMatthew G. Knepley   if (nparts == 1) {
341*70034214SMatthew G. Knepley     ierr = PetscMemzero(assignment, nvtxs * sizeof(PetscInt));
342*70034214SMatthew G. Knepley   } else {
343*70034214SMatthew G. Knepley     if (vtxdist[1] == vtxdist[nparts]) {
344*70034214SMatthew G. Knepley       if (!rank) {
345*70034214SMatthew G. Knepley         PetscStackPush("METIS_PartGraphKway");
346*70034214SMatthew G. Knepley         ierr = METIS_PartGraphKway(&nvtxs, &ncon, xadj, adjncy, vwgt, NULL, adjwgt, &nparts, tpwgts, ubvec, NULL, &edgeCut, assignment);
347*70034214SMatthew G. Knepley         PetscStackPop;
348*70034214SMatthew G. Knepley         if (ierr != METIS_OK) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in METIS_PartGraphKway()");
349*70034214SMatthew G. Knepley       }
350*70034214SMatthew G. Knepley     } else {
351*70034214SMatthew G. Knepley       PetscStackPush("ParMETIS_V3_PartKway");
352*70034214SMatthew G. Knepley       ierr = ParMETIS_V3_PartKway(vtxdist, xadj, adjncy, vwgt, adjwgt, &wgtflag, &numflag, &ncon, &nparts, tpwgts, ubvec, options, &edgeCut, assignment, &comm);
353*70034214SMatthew G. Knepley       PetscStackPop;
354*70034214SMatthew G. Knepley       if (ierr != METIS_OK) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_LIB, "Error in ParMETIS_V3_PartKway()");
355*70034214SMatthew G. Knepley     }
356*70034214SMatthew G. Knepley   }
357*70034214SMatthew G. Knepley   /* Convert to PetscSection+IS */
358*70034214SMatthew G. Knepley   ierr = PetscSectionCreate(comm, partSection);CHKERRQ(ierr);
359*70034214SMatthew G. Knepley   ierr = PetscSectionSetChart(*partSection, 0, commSize);CHKERRQ(ierr);
360*70034214SMatthew G. Knepley   for (v = 0; v < nvtxs; ++v) {
361*70034214SMatthew G. Knepley     ierr = PetscSectionAddDof(*partSection, assignment[v], 1);CHKERRQ(ierr);
362*70034214SMatthew G. Knepley   }
363*70034214SMatthew G. Knepley   ierr = PetscSectionSetUp(*partSection);CHKERRQ(ierr);
364*70034214SMatthew G. Knepley   ierr = PetscMalloc1(nvtxs, &points);CHKERRQ(ierr);
365*70034214SMatthew G. Knepley   for (p = 0, i = 0; p < commSize; ++p) {
366*70034214SMatthew G. Knepley     for (v = 0; v < nvtxs; ++v) {
367*70034214SMatthew G. Knepley       if (assignment[v] == p) points[i++] = v;
368*70034214SMatthew G. Knepley     }
369*70034214SMatthew G. Knepley   }
370*70034214SMatthew G. Knepley   if (i != nvtxs) SETERRQ2(comm, PETSC_ERR_PLIB, "Number of points %D should be %D", i, nvtxs);
371*70034214SMatthew G. Knepley   ierr = ISCreateGeneral(comm, nvtxs, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
372*70034214SMatthew G. Knepley   ierr = PetscFree4(vtxdist,tpwgts,ubvec,assignment);CHKERRQ(ierr);
373*70034214SMatthew G. Knepley   PetscFunctionReturn(0);
374*70034214SMatthew G. Knepley }
375*70034214SMatthew G. Knepley #endif
376*70034214SMatthew G. Knepley 
377*70034214SMatthew G. Knepley #undef __FUNCT__
378*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexEnlargePartition"
379*70034214SMatthew G. Knepley /* Expand the partition by BFS on the adjacency graph */
380*70034214SMatthew G. Knepley PetscErrorCode DMPlexEnlargePartition(DM dm, const PetscInt start[], const PetscInt adjacency[], PetscSection origPartSection, IS origPartition, PetscSection *partSection, IS *partition)
381*70034214SMatthew G. Knepley {
382*70034214SMatthew G. Knepley   PetscHashI      h;
383*70034214SMatthew G. Knepley   const PetscInt *points;
384*70034214SMatthew G. Knepley   PetscInt      **tmpPoints, *newPoints, totPoints = 0;
385*70034214SMatthew G. Knepley   PetscInt        pStart, pEnd, part, q;
386*70034214SMatthew G. Knepley   PetscBool       useCone;
387*70034214SMatthew G. Knepley   PetscErrorCode  ierr;
388*70034214SMatthew G. Knepley 
389*70034214SMatthew G. Knepley   PetscFunctionBegin;
390*70034214SMatthew G. Knepley   PetscHashICreate(h);
391*70034214SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), partSection);CHKERRQ(ierr);
392*70034214SMatthew G. Knepley   ierr = PetscSectionGetChart(origPartSection, &pStart, &pEnd);CHKERRQ(ierr);
393*70034214SMatthew G. Knepley   ierr = PetscSectionSetChart(*partSection, pStart, pEnd);CHKERRQ(ierr);
394*70034214SMatthew G. Knepley   ierr = ISGetIndices(origPartition, &points);CHKERRQ(ierr);
395*70034214SMatthew G. Knepley   ierr = PetscMalloc1((pEnd - pStart), &tmpPoints);CHKERRQ(ierr);
396*70034214SMatthew G. Knepley   ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
397*70034214SMatthew G. Knepley   ierr = DMPlexSetAdjacencyUseCone(dm, PETSC_TRUE);CHKERRQ(ierr);
398*70034214SMatthew G. Knepley   for (part = pStart; part < pEnd; ++part) {
399*70034214SMatthew G. Knepley     PetscInt *adj = NULL;
400*70034214SMatthew G. Knepley     PetscInt  numPoints, nP, numNewPoints, off, p, n = 0;
401*70034214SMatthew G. Knepley 
402*70034214SMatthew G. Knepley     PetscHashIClear(h);
403*70034214SMatthew G. Knepley     ierr = PetscSectionGetDof(origPartSection, part, &numPoints);CHKERRQ(ierr);
404*70034214SMatthew G. Knepley     ierr = PetscSectionGetOffset(origPartSection, part, &off);CHKERRQ(ierr);
405*70034214SMatthew G. Knepley     /* Add all existing points to h */
406*70034214SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
407*70034214SMatthew G. Knepley       const PetscInt point = points[off+p];
408*70034214SMatthew G. Knepley       PetscHashIAdd(h, point, 1);
409*70034214SMatthew G. Knepley     }
410*70034214SMatthew G. Knepley     PetscHashISize(h, nP);
411*70034214SMatthew 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);
412*70034214SMatthew G. Knepley     /* Add all points in next BFS level */
413*70034214SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
414*70034214SMatthew G. Knepley       const PetscInt point   = points[off+p];
415*70034214SMatthew G. Knepley       PetscInt       adjSize = PETSC_DETERMINE, a;
416*70034214SMatthew G. Knepley 
417*70034214SMatthew G. Knepley       ierr = DMPlexGetAdjacency(dm, point, &adjSize, &adj);CHKERRQ(ierr);
418*70034214SMatthew G. Knepley       for (a = 0; a < adjSize; ++a) PetscHashIAdd(h, adj[a], 1);
419*70034214SMatthew G. Knepley     }
420*70034214SMatthew G. Knepley     PetscHashISize(h, numNewPoints);
421*70034214SMatthew G. Knepley     ierr = PetscSectionSetDof(*partSection, part, numNewPoints);CHKERRQ(ierr);
422*70034214SMatthew G. Knepley     ierr = PetscMalloc1(numNewPoints, &tmpPoints[part]);CHKERRQ(ierr);
423*70034214SMatthew G. Knepley     ierr = PetscHashIGetKeys(h, &n, tmpPoints[part]);CHKERRQ(ierr);
424*70034214SMatthew G. Knepley     ierr = PetscFree(adj);CHKERRQ(ierr);
425*70034214SMatthew G. Knepley     totPoints += numNewPoints;
426*70034214SMatthew G. Knepley   }
427*70034214SMatthew G. Knepley   ierr = DMPlexSetAdjacencyUseCone(dm, useCone);CHKERRQ(ierr);
428*70034214SMatthew G. Knepley   ierr = ISRestoreIndices(origPartition, &points);CHKERRQ(ierr);
429*70034214SMatthew G. Knepley   PetscHashIDestroy(h);
430*70034214SMatthew G. Knepley   ierr = PetscSectionSetUp(*partSection);CHKERRQ(ierr);
431*70034214SMatthew G. Knepley   ierr = PetscMalloc1(totPoints, &newPoints);CHKERRQ(ierr);
432*70034214SMatthew G. Knepley   for (part = pStart, q = 0; part < pEnd; ++part) {
433*70034214SMatthew G. Knepley     PetscInt numPoints, p;
434*70034214SMatthew G. Knepley 
435*70034214SMatthew G. Knepley     ierr = PetscSectionGetDof(*partSection, part, &numPoints);CHKERRQ(ierr);
436*70034214SMatthew G. Knepley     for (p = 0; p < numPoints; ++p, ++q) newPoints[q] = tmpPoints[part][p];
437*70034214SMatthew G. Knepley     ierr = PetscFree(tmpPoints[part]);CHKERRQ(ierr);
438*70034214SMatthew G. Knepley   }
439*70034214SMatthew G. Knepley   ierr = PetscFree(tmpPoints);CHKERRQ(ierr);
440*70034214SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), totPoints, newPoints, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
441*70034214SMatthew G. Knepley   PetscFunctionReturn(0);
442*70034214SMatthew G. Knepley }
443*70034214SMatthew G. Knepley 
444*70034214SMatthew G. Knepley #undef __FUNCT__
445*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexCreatePartition"
446*70034214SMatthew G. Knepley /*
447*70034214SMatthew G. Knepley   DMPlexCreatePartition - Create a non-overlapping partition of the points at the given height
448*70034214SMatthew G. Knepley 
449*70034214SMatthew G. Knepley   Collective on DM
450*70034214SMatthew G. Knepley 
451*70034214SMatthew G. Knepley   Input Parameters:
452*70034214SMatthew G. Knepley   + dm - The DM
453*70034214SMatthew G. Knepley   . height - The height for points in the partition
454*70034214SMatthew G. Knepley   - enlarge - Expand each partition with neighbors
455*70034214SMatthew G. Knepley 
456*70034214SMatthew G. Knepley   Output Parameters:
457*70034214SMatthew G. Knepley   + partSection - The PetscSection giving the division of points by partition
458*70034214SMatthew G. Knepley   . partition - The list of points by partition
459*70034214SMatthew G. Knepley   . origPartSection - If enlarge is true, the PetscSection giving the division of points before enlarging by partition, otherwise NULL
460*70034214SMatthew G. Knepley   - origPartition - If enlarge is true, the list of points before enlarging by partition, otherwise NULL
461*70034214SMatthew G. Knepley 
462*70034214SMatthew G. Knepley   Level: developer
463*70034214SMatthew G. Knepley 
464*70034214SMatthew G. Knepley .seealso DMPlexDistribute()
465*70034214SMatthew G. Knepley */
466*70034214SMatthew G. Knepley PetscErrorCode DMPlexCreatePartition(DM dm, const char name[], PetscInt height, PetscBool enlarge, PetscSection *partSection, IS *partition, PetscSection *origPartSection, IS *origPartition)
467*70034214SMatthew G. Knepley {
468*70034214SMatthew G. Knepley   char           partname[1024];
469*70034214SMatthew G. Knepley   PetscBool      isChaco = PETSC_FALSE, isMetis = PETSC_FALSE, flg;
470*70034214SMatthew G. Knepley   PetscMPIInt    size;
471*70034214SMatthew G. Knepley   PetscErrorCode ierr;
472*70034214SMatthew G. Knepley 
473*70034214SMatthew G. Knepley   PetscFunctionBegin;
474*70034214SMatthew G. Knepley   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)dm), &size);CHKERRQ(ierr);
475*70034214SMatthew G. Knepley 
476*70034214SMatthew G. Knepley   *origPartSection = NULL;
477*70034214SMatthew G. Knepley   *origPartition   = NULL;
478*70034214SMatthew G. Knepley   if (size == 1) {
479*70034214SMatthew G. Knepley     PetscInt *points;
480*70034214SMatthew G. Knepley     PetscInt  cStart, cEnd, c;
481*70034214SMatthew G. Knepley 
482*70034214SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
483*70034214SMatthew G. Knepley     ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), partSection);CHKERRQ(ierr);
484*70034214SMatthew G. Knepley     ierr = PetscSectionSetChart(*partSection, 0, size);CHKERRQ(ierr);
485*70034214SMatthew G. Knepley     ierr = PetscSectionSetDof(*partSection, 0, cEnd-cStart);CHKERRQ(ierr);
486*70034214SMatthew G. Knepley     ierr = PetscSectionSetUp(*partSection);CHKERRQ(ierr);
487*70034214SMatthew G. Knepley     ierr = PetscMalloc1((cEnd - cStart), &points);CHKERRQ(ierr);
488*70034214SMatthew G. Knepley     for (c = cStart; c < cEnd; ++c) points[c] = c;
489*70034214SMatthew G. Knepley     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), cEnd-cStart, points, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
490*70034214SMatthew G. Knepley     PetscFunctionReturn(0);
491*70034214SMatthew G. Knepley   }
492*70034214SMatthew G. Knepley   ierr = PetscOptionsGetString(((PetscObject) dm)->prefix, "-dm_plex_partitioner", partname, 1024, &flg);CHKERRQ(ierr);
493*70034214SMatthew G. Knepley   if (flg) name = partname;
494*70034214SMatthew G. Knepley   if (name) {
495*70034214SMatthew G. Knepley     ierr = PetscStrcmp(name, "chaco", &isChaco);CHKERRQ(ierr);
496*70034214SMatthew G. Knepley     ierr = PetscStrcmp(name, "metis", &isMetis);CHKERRQ(ierr);
497*70034214SMatthew G. Knepley   }
498*70034214SMatthew G. Knepley   if (height == 0) {
499*70034214SMatthew G. Knepley     PetscInt  numVertices;
500*70034214SMatthew G. Knepley     PetscInt *start     = NULL;
501*70034214SMatthew G. Knepley     PetscInt *adjacency = NULL;
502*70034214SMatthew G. Knepley 
503*70034214SMatthew G. Knepley     ierr = DMPlexCreateNeighborCSR(dm, 0, &numVertices, &start, &adjacency);CHKERRQ(ierr);
504*70034214SMatthew G. Knepley     if (!name || isChaco) {
505*70034214SMatthew G. Knepley #if defined(PETSC_HAVE_CHACO)
506*70034214SMatthew G. Knepley       ierr = DMPlexPartition_Chaco(dm, numVertices, start, adjacency, partSection, partition);CHKERRQ(ierr);
507*70034214SMatthew G. Knepley #else
508*70034214SMatthew G. Knepley       SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Mesh partitioning needs external package support.\nPlease reconfigure with --download-chaco.");
509*70034214SMatthew G. Knepley #endif
510*70034214SMatthew G. Knepley     } else if (isMetis) {
511*70034214SMatthew G. Knepley #if defined(PETSC_HAVE_PARMETIS)
512*70034214SMatthew G. Knepley       ierr = DMPlexPartition_ParMetis(dm, numVertices, start, adjacency, partSection, partition);CHKERRQ(ierr);
513*70034214SMatthew G. Knepley #endif
514*70034214SMatthew G. Knepley     } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Unknown mesh partitioning package %s", name);
515*70034214SMatthew G. Knepley     if (enlarge) {
516*70034214SMatthew G. Knepley       *origPartSection = *partSection;
517*70034214SMatthew G. Knepley       *origPartition   = *partition;
518*70034214SMatthew G. Knepley 
519*70034214SMatthew G. Knepley       ierr = DMPlexEnlargePartition(dm, start, adjacency, *origPartSection, *origPartition, partSection, partition);CHKERRQ(ierr);
520*70034214SMatthew G. Knepley     }
521*70034214SMatthew G. Knepley     ierr = PetscFree(start);CHKERRQ(ierr);
522*70034214SMatthew G. Knepley     ierr = PetscFree(adjacency);CHKERRQ(ierr);
523*70034214SMatthew G. Knepley # if 0
524*70034214SMatthew G. Knepley   } else if (height == 1) {
525*70034214SMatthew G. Knepley     /* Build the dual graph for faces and partition the hypergraph */
526*70034214SMatthew G. Knepley     PetscInt numEdges;
527*70034214SMatthew G. Knepley 
528*70034214SMatthew G. Knepley     buildFaceCSRV(mesh, mesh->getFactory()->getNumbering(mesh, mesh->depth()-1), &numEdges, &start, &adjacency, GraphPartitioner::zeroBase());
529*70034214SMatthew G. Knepley     GraphPartitioner().partition(numEdges, start, adjacency, partition, manager);
530*70034214SMatthew G. Knepley     destroyCSR(numEdges, start, adjacency);
531*70034214SMatthew G. Knepley #endif
532*70034214SMatthew G. Knepley   } else SETERRQ1(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Invalid partition height %D", height);
533*70034214SMatthew G. Knepley   PetscFunctionReturn(0);
534*70034214SMatthew G. Knepley }
535*70034214SMatthew G. Knepley 
536*70034214SMatthew G. Knepley #undef __FUNCT__
537*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexCreatePartitionClosure"
538*70034214SMatthew G. Knepley PetscErrorCode DMPlexCreatePartitionClosure(DM dm, PetscSection pointSection, IS pointPartition, PetscSection *section, IS *partition)
539*70034214SMatthew G. Knepley {
540*70034214SMatthew G. Knepley   /* const PetscInt  height = 0; */
541*70034214SMatthew G. Knepley   const PetscInt *partArray;
542*70034214SMatthew G. Knepley   PetscInt       *allPoints, *packPoints;
543*70034214SMatthew G. Knepley   PetscInt        rStart, rEnd, rank, pStart, pEnd, newSize;
544*70034214SMatthew G. Knepley   PetscErrorCode  ierr;
545*70034214SMatthew G. Knepley   PetscBT         bt;
546*70034214SMatthew G. Knepley   PetscSegBuffer  segpack,segpart;
547*70034214SMatthew G. Knepley 
548*70034214SMatthew G. Knepley   PetscFunctionBegin;
549*70034214SMatthew G. Knepley   ierr = PetscSectionGetChart(pointSection, &rStart, &rEnd);CHKERRQ(ierr);
550*70034214SMatthew G. Knepley   ierr = ISGetIndices(pointPartition, &partArray);CHKERRQ(ierr);
551*70034214SMatthew G. Knepley   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm), section);CHKERRQ(ierr);
552*70034214SMatthew G. Knepley   ierr = PetscSectionSetChart(*section, rStart, rEnd);CHKERRQ(ierr);
553*70034214SMatthew G. Knepley   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
554*70034214SMatthew G. Knepley   ierr = PetscBTCreate(pEnd-pStart,&bt);CHKERRQ(ierr);
555*70034214SMatthew G. Knepley   ierr = PetscSegBufferCreate(sizeof(PetscInt),1000,&segpack);CHKERRQ(ierr);
556*70034214SMatthew G. Knepley   ierr = PetscSegBufferCreate(sizeof(PetscInt),1000,&segpart);CHKERRQ(ierr);
557*70034214SMatthew G. Knepley   for (rank = rStart; rank < rEnd; ++rank) {
558*70034214SMatthew G. Knepley     PetscInt partSize = 0, numPoints, offset, p, *PETSC_RESTRICT placePoints;
559*70034214SMatthew G. Knepley 
560*70034214SMatthew G. Knepley     ierr = PetscSectionGetDof(pointSection, rank, &numPoints);CHKERRQ(ierr);
561*70034214SMatthew G. Knepley     ierr = PetscSectionGetOffset(pointSection, rank, &offset);CHKERRQ(ierr);
562*70034214SMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
563*70034214SMatthew G. Knepley       PetscInt  point   = partArray[offset+p], closureSize, c;
564*70034214SMatthew G. Knepley       PetscInt *closure = NULL;
565*70034214SMatthew G. Knepley 
566*70034214SMatthew G. Knepley       /* TODO Include support for height > 0 case */
567*70034214SMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
568*70034214SMatthew G. Knepley       for (c=0; c<closureSize; c++) {
569*70034214SMatthew G. Knepley         PetscInt cpoint = closure[c*2];
570*70034214SMatthew G. Knepley         if (!PetscBTLookupSet(bt,cpoint-pStart)) {
571*70034214SMatthew G. Knepley           PetscInt *PETSC_RESTRICT pt;
572*70034214SMatthew G. Knepley           partSize++;
573*70034214SMatthew G. Knepley           ierr = PetscSegBufferGetInts(segpart,1,&pt);CHKERRQ(ierr);
574*70034214SMatthew G. Knepley           *pt = cpoint;
575*70034214SMatthew G. Knepley         }
576*70034214SMatthew G. Knepley       }
577*70034214SMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, point, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
578*70034214SMatthew G. Knepley     }
579*70034214SMatthew G. Knepley     ierr = PetscSectionSetDof(*section, rank, partSize);CHKERRQ(ierr);
580*70034214SMatthew G. Knepley     ierr = PetscSegBufferGetInts(segpack,partSize,&placePoints);CHKERRQ(ierr);
581*70034214SMatthew G. Knepley     ierr = PetscSegBufferExtractTo(segpart,placePoints);CHKERRQ(ierr);
582*70034214SMatthew G. Knepley     ierr = PetscSortInt(partSize,placePoints);CHKERRQ(ierr);
583*70034214SMatthew G. Knepley     for (p=0; p<partSize; p++) {ierr = PetscBTClear(bt,placePoints[p]-pStart);CHKERRQ(ierr);}
584*70034214SMatthew G. Knepley   }
585*70034214SMatthew G. Knepley   ierr = PetscBTDestroy(&bt);CHKERRQ(ierr);
586*70034214SMatthew G. Knepley   ierr = PetscSegBufferDestroy(&segpart);CHKERRQ(ierr);
587*70034214SMatthew G. Knepley 
588*70034214SMatthew G. Knepley   ierr = PetscSectionSetUp(*section);CHKERRQ(ierr);
589*70034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(*section, &newSize);CHKERRQ(ierr);
590*70034214SMatthew G. Knepley   ierr = PetscMalloc1(newSize, &allPoints);CHKERRQ(ierr);
591*70034214SMatthew G. Knepley 
592*70034214SMatthew G. Knepley   ierr = PetscSegBufferExtractInPlace(segpack,&packPoints);CHKERRQ(ierr);
593*70034214SMatthew G. Knepley   for (rank = rStart; rank < rEnd; ++rank) {
594*70034214SMatthew G. Knepley     PetscInt numPoints, offset;
595*70034214SMatthew G. Knepley 
596*70034214SMatthew G. Knepley     ierr = PetscSectionGetDof(*section, rank, &numPoints);CHKERRQ(ierr);
597*70034214SMatthew G. Knepley     ierr = PetscSectionGetOffset(*section, rank, &offset);CHKERRQ(ierr);
598*70034214SMatthew G. Knepley     ierr = PetscMemcpy(&allPoints[offset], packPoints, numPoints * sizeof(PetscInt));CHKERRQ(ierr);
599*70034214SMatthew G. Knepley     packPoints += numPoints;
600*70034214SMatthew G. Knepley   }
601*70034214SMatthew G. Knepley 
602*70034214SMatthew G. Knepley   ierr = PetscSegBufferDestroy(&segpack);CHKERRQ(ierr);
603*70034214SMatthew G. Knepley   ierr = ISRestoreIndices(pointPartition, &partArray);CHKERRQ(ierr);
604*70034214SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), newSize, allPoints, PETSC_OWN_POINTER, partition);CHKERRQ(ierr);
605*70034214SMatthew G. Knepley   PetscFunctionReturn(0);
606*70034214SMatthew G. Knepley }
607