xref: /petsc/src/dm/impls/plex/plexdistribute.c (revision 8276244a08780f88cbecdeda759dff4b2f999d0f)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 
3 #undef __FUNCT__
4 #define __FUNCT__ "DMPlexSetAdjacencyUseCone"
5 /*@
6   DMPlexSetAdjacencyUseCone - Define adjacency in the mesh using either the cone or the support first
7 
8   Input Parameters:
9 + dm      - The DM object
10 - useCone - Flag to use the cone first
11 
12   Level: intermediate
13 
14   Notes:
15 $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
16 $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
17 $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
18 
19 .seealso: DMPlexGetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
20 @*/
21 PetscErrorCode DMPlexSetAdjacencyUseCone(DM dm, PetscBool useCone)
22 {
23   DM_Plex *mesh = (DM_Plex *) dm->data;
24 
25   PetscFunctionBegin;
26   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
27   mesh->useCone = useCone;
28   PetscFunctionReturn(0);
29 }
30 
31 #undef __FUNCT__
32 #define __FUNCT__ "DMPlexGetAdjacencyUseCone"
33 /*@
34   DMPlexGetAdjacencyUseCone - Query whether adjacency in the mesh uses the cone or the support first
35 
36   Input Parameter:
37 . dm      - The DM object
38 
39   Output Parameter:
40 . useCone - Flag to use the cone first
41 
42   Level: intermediate
43 
44   Notes:
45 $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
46 $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
47 $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
48 
49 .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
50 @*/
51 PetscErrorCode DMPlexGetAdjacencyUseCone(DM dm, PetscBool *useCone)
52 {
53   DM_Plex *mesh = (DM_Plex *) dm->data;
54 
55   PetscFunctionBegin;
56   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
57   PetscValidIntPointer(useCone, 2);
58   *useCone = mesh->useCone;
59   PetscFunctionReturn(0);
60 }
61 
62 #undef __FUNCT__
63 #define __FUNCT__ "DMPlexSetAdjacencyUseClosure"
64 /*@
65   DMPlexSetAdjacencyUseClosure - Define adjacency in the mesh using the transitive closure
66 
67   Input Parameters:
68 + dm      - The DM object
69 - useClosure - Flag to use the closure
70 
71   Level: intermediate
72 
73   Notes:
74 $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
75 $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
76 $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
77 
78 .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
79 @*/
80 PetscErrorCode DMPlexSetAdjacencyUseClosure(DM dm, PetscBool useClosure)
81 {
82   DM_Plex *mesh = (DM_Plex *) dm->data;
83 
84   PetscFunctionBegin;
85   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
86   mesh->useClosure = useClosure;
87   PetscFunctionReturn(0);
88 }
89 
90 #undef __FUNCT__
91 #define __FUNCT__ "DMPlexGetAdjacencyUseClosure"
92 /*@
93   DMPlexGetAdjacencyUseClosure - Query whether adjacency in the mesh uses the transitive closure
94 
95   Input Parameter:
96 . dm      - The DM object
97 
98   Output Parameter:
99 . useClosure - Flag to use the closure
100 
101   Level: intermediate
102 
103   Notes:
104 $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
105 $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
106 $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
107 
108 .seealso: DMPlexSetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
109 @*/
110 PetscErrorCode DMPlexGetAdjacencyUseClosure(DM dm, PetscBool *useClosure)
111 {
112   DM_Plex *mesh = (DM_Plex *) dm->data;
113 
114   PetscFunctionBegin;
115   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
116   PetscValidIntPointer(useClosure, 2);
117   *useClosure = mesh->useClosure;
118   PetscFunctionReturn(0);
119 }
120 
121 #undef __FUNCT__
122 #define __FUNCT__ "DMPlexSetAdjacencyUseAnchors"
123 /*@
124   DMPlexSetAdjacencyUseAnchors - Define adjacency in the mesh using the point-to-point constraints.
125 
126   Input Parameters:
127 + dm      - The DM object
128 - useAnchors - Flag to use the constraints.  If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place.
129 
130   Level: intermediate
131 
132 .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
133 @*/
134 PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM dm, PetscBool useAnchors)
135 {
136   DM_Plex *mesh = (DM_Plex *) dm->data;
137 
138   PetscFunctionBegin;
139   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
140   mesh->useAnchors = useAnchors;
141   PetscFunctionReturn(0);
142 }
143 
144 #undef __FUNCT__
145 #define __FUNCT__ "DMPlexGetAdjacencyUseAnchors"
146 /*@
147   DMPlexGetAdjacencyUseAnchors - Query whether adjacency in the mesh uses the point-to-point constraints.
148 
149   Input Parameter:
150 . dm      - The DM object
151 
152   Output Parameter:
153 . useAnchors - Flag to use the closure.  If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place.
154 
155   Level: intermediate
156 
157 .seealso: DMPlexSetAdjacencyUseAnchors(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
158 @*/
159 PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM dm, PetscBool *useAnchors)
160 {
161   DM_Plex *mesh = (DM_Plex *) dm->data;
162 
163   PetscFunctionBegin;
164   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
165   PetscValidIntPointer(useAnchors, 2);
166   *useAnchors = mesh->useAnchors;
167   PetscFunctionReturn(0);
168 }
169 
170 #undef __FUNCT__
171 #define __FUNCT__ "DMPlexGetAdjacency_Cone_Internal"
172 static PetscErrorCode DMPlexGetAdjacency_Cone_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
173 {
174   const PetscInt *cone = NULL;
175   PetscInt        numAdj = 0, maxAdjSize = *adjSize, coneSize, c;
176   PetscErrorCode  ierr;
177 
178   PetscFunctionBeginHot;
179   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
180   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
181   for (c = 0; c <= coneSize; ++c) {
182     const PetscInt  point   = !c ? p : cone[c-1];
183     const PetscInt *support = NULL;
184     PetscInt        supportSize, s, q;
185 
186     ierr = DMPlexGetSupportSize(dm, point, &supportSize);CHKERRQ(ierr);
187     ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr);
188     for (s = 0; s < supportSize; ++s) {
189       for (q = 0; q < numAdj || (adj[numAdj++] = support[s],0); ++q) {
190         if (support[s] == adj[q]) break;
191       }
192       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
193     }
194   }
195   *adjSize = numAdj;
196   PetscFunctionReturn(0);
197 }
198 
199 #undef __FUNCT__
200 #define __FUNCT__ "DMPlexGetAdjacency_Support_Internal"
201 static PetscErrorCode DMPlexGetAdjacency_Support_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
202 {
203   const PetscInt *support = NULL;
204   PetscInt        numAdj   = 0, maxAdjSize = *adjSize, supportSize, s;
205   PetscErrorCode  ierr;
206 
207   PetscFunctionBeginHot;
208   ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr);
209   ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr);
210   for (s = 0; s <= supportSize; ++s) {
211     const PetscInt  point = !s ? p : support[s-1];
212     const PetscInt *cone  = NULL;
213     PetscInt        coneSize, c, q;
214 
215     ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
216     ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
217     for (c = 0; c < coneSize; ++c) {
218       for (q = 0; q < numAdj || (adj[numAdj++] = cone[c],0); ++q) {
219         if (cone[c] == adj[q]) break;
220       }
221       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
222     }
223   }
224   *adjSize = numAdj;
225   PetscFunctionReturn(0);
226 }
227 
228 #undef __FUNCT__
229 #define __FUNCT__ "DMPlexGetAdjacency_Transitive_Internal"
230 static PetscErrorCode DMPlexGetAdjacency_Transitive_Internal(DM dm, PetscInt p, PetscBool useClosure, PetscInt *adjSize, PetscInt adj[])
231 {
232   PetscInt      *star = NULL;
233   PetscInt       numAdj = 0, maxAdjSize = *adjSize, starSize, s;
234   PetscErrorCode ierr;
235 
236   PetscFunctionBeginHot;
237   ierr = DMPlexGetTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr);
238   for (s = 0; s < starSize*2; s += 2) {
239     const PetscInt *closure = NULL;
240     PetscInt        closureSize, c, q;
241 
242     ierr = DMPlexGetTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
243     for (c = 0; c < closureSize*2; c += 2) {
244       for (q = 0; q < numAdj || (adj[numAdj++] = closure[c],0); ++q) {
245         if (closure[c] == adj[q]) break;
246       }
247       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
248     }
249     ierr = DMPlexRestoreTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
250   }
251   ierr = DMPlexRestoreTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr);
252   *adjSize = numAdj;
253   PetscFunctionReturn(0);
254 }
255 
256 #undef __FUNCT__
257 #define __FUNCT__ "DMPlexGetAdjacency_Internal"
258 PetscErrorCode DMPlexGetAdjacency_Internal(DM dm, PetscInt p, PetscBool useCone, PetscBool useTransitiveClosure, PetscBool useAnchors, PetscInt *adjSize, PetscInt *adj[])
259 {
260   static PetscInt asiz = 0;
261   PetscInt maxAnchors = 1;
262   PetscInt aStart = -1, aEnd = -1;
263   PetscInt maxAdjSize;
264   PetscSection aSec = NULL;
265   IS aIS = NULL;
266   const PetscInt *anchors;
267   PetscErrorCode  ierr;
268 
269   PetscFunctionBeginHot;
270   if (useAnchors) {
271     ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr);
272     if (aSec) {
273       ierr = PetscSectionGetMaxDof(aSec,&maxAnchors);CHKERRQ(ierr);
274       maxAnchors = PetscMax(1,maxAnchors);
275       ierr = PetscSectionGetChart(aSec,&aStart,&aEnd);CHKERRQ(ierr);
276       ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr);
277     }
278   }
279   if (!*adj) {
280     PetscInt depth, coneSeries, supportSeries, maxC, maxS, pStart, pEnd;
281 
282     ierr  = DMPlexGetChart(dm, &pStart,&pEnd);CHKERRQ(ierr);
283     ierr  = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
284     ierr  = DMPlexGetMaxSizes(dm, &maxC, &maxS);CHKERRQ(ierr);
285     coneSeries    = (maxC > 1) ? ((PetscPowInt(maxC,depth+1)-1)/(maxC-1)) : depth+1;
286     supportSeries = (maxS > 1) ? ((PetscPowInt(maxS,depth+1)-1)/(maxS-1)) : depth+1;
287     asiz  = PetscMax(PetscPowInt(maxS,depth)*coneSeries,PetscPowInt(maxC,depth)*supportSeries);
288     asiz *= maxAnchors;
289     asiz  = PetscMin(asiz,pEnd-pStart);
290     ierr  = PetscMalloc1(asiz,adj);CHKERRQ(ierr);
291   }
292   if (*adjSize < 0) *adjSize = asiz;
293   maxAdjSize = *adjSize;
294   if (useTransitiveClosure) {
295     ierr = DMPlexGetAdjacency_Transitive_Internal(dm, p, useCone, adjSize, *adj);CHKERRQ(ierr);
296   } else if (useCone) {
297     ierr = DMPlexGetAdjacency_Cone_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
298   } else {
299     ierr = DMPlexGetAdjacency_Support_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
300   }
301   if (useAnchors && aSec) {
302     PetscInt origSize = *adjSize;
303     PetscInt numAdj = origSize;
304     PetscInt i = 0, j;
305     PetscInt *orig = *adj;
306 
307     while (i < origSize) {
308       PetscInt p = orig[i];
309       PetscInt aDof = 0;
310 
311       if (p >= aStart && p < aEnd) {
312         ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr);
313       }
314       if (aDof) {
315         PetscInt aOff;
316         PetscInt s, q;
317 
318         for (j = i + 1; j < numAdj; j++) {
319           orig[j - 1] = orig[j];
320         }
321         origSize--;
322         numAdj--;
323         ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr);
324         for (s = 0; s < aDof; ++s) {
325           for (q = 0; q < numAdj || (orig[numAdj++] = anchors[aOff+s],0); ++q) {
326             if (anchors[aOff+s] == orig[q]) break;
327           }
328           if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
329         }
330       }
331       else {
332         i++;
333       }
334     }
335     *adjSize = numAdj;
336     ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr);
337   }
338   PetscFunctionReturn(0);
339 }
340 
341 #undef __FUNCT__
342 #define __FUNCT__ "DMPlexGetAdjacency"
343 /*@
344   DMPlexGetAdjacency - Return all points adjacent to the given point
345 
346   Input Parameters:
347 + dm - The DM object
348 . p  - The point
349 . adjSize - The maximum size of adj if it is non-NULL, or PETSC_DETERMINE
350 - adj - Either NULL so that the array is allocated, or an existing array with size adjSize
351 
352   Output Parameters:
353 + adjSize - The number of adjacent points
354 - adj - The adjacent points
355 
356   Level: advanced
357 
358   Notes: The user must PetscFree the adj array if it was not passed in.
359 
360 .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexDistribute(), DMCreateMatrix(), DMPlexPreallocateOperator()
361 @*/
362 PetscErrorCode DMPlexGetAdjacency(DM dm, PetscInt p, PetscInt *adjSize, PetscInt *adj[])
363 {
364   DM_Plex       *mesh = (DM_Plex *) dm->data;
365   PetscErrorCode ierr;
366 
367   PetscFunctionBeginHot;
368   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
369   PetscValidPointer(adjSize,3);
370   PetscValidPointer(adj,4);
371   ierr = DMPlexGetAdjacency_Internal(dm, p, mesh->useCone, mesh->useClosure, mesh->useAnchors, adjSize, adj);CHKERRQ(ierr);
372   PetscFunctionReturn(0);
373 }
374 
375 #undef __FUNCT__
376 #define __FUNCT__ "DMPlexCreateTwoSidedProcessSF"
377 /*@
378   DMPlexCreateTwoSidedProcessSF - Create an SF which just has process connectivity
379 
380   Collective on DM
381 
382   Input Parameters:
383 + dm      - The DM
384 - sfPoint - The PetscSF which encodes point connectivity
385 
386   Output Parameters:
387 + processRanks - A list of process neighbors, or NULL
388 - sfProcess    - An SF encoding the two-sided process connectivity, or NULL
389 
390   Level: developer
391 
392 .seealso: PetscSFCreate(), DMPlexCreateProcessSF()
393 @*/
394 PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM dm, PetscSF sfPoint, PetscSection rootRankSection, IS rootRanks, PetscSection leafRankSection, IS leafRanks, IS *processRanks, PetscSF *sfProcess)
395 {
396   const PetscSFNode *remotePoints;
397   PetscInt          *localPointsNew;
398   PetscSFNode       *remotePointsNew;
399   const PetscInt    *nranks;
400   PetscInt          *ranksNew;
401   PetscBT            neighbors;
402   PetscInt           pStart, pEnd, p, numLeaves, l, numNeighbors, n;
403   PetscMPIInt        numProcs, proc, rank;
404   PetscErrorCode     ierr;
405 
406   PetscFunctionBegin;
407   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
408   PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2);
409   if (processRanks) {PetscValidPointer(processRanks, 3);}
410   if (sfProcess)    {PetscValidPointer(sfProcess, 4);}
411   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr);
412   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
413   ierr = PetscSFGetGraph(sfPoint, NULL, &numLeaves, NULL, &remotePoints);CHKERRQ(ierr);
414   ierr = PetscBTCreate(numProcs, &neighbors);CHKERRQ(ierr);
415   ierr = PetscBTMemzero(numProcs, neighbors);CHKERRQ(ierr);
416   /* Compute root-to-leaf process connectivity */
417   ierr = PetscSectionGetChart(rootRankSection, &pStart, &pEnd);CHKERRQ(ierr);
418   ierr = ISGetIndices(rootRanks, &nranks);CHKERRQ(ierr);
419   for (p = pStart; p < pEnd; ++p) {
420     PetscInt ndof, noff, n;
421 
422     ierr = PetscSectionGetDof(rootRankSection, p, &ndof);CHKERRQ(ierr);
423     ierr = PetscSectionGetOffset(rootRankSection, p, &noff);CHKERRQ(ierr);
424     for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);CHKERRQ(ierr);}
425   }
426   ierr = ISRestoreIndices(rootRanks, &nranks);CHKERRQ(ierr);
427   /* Compute leaf-to-neighbor process connectivity */
428   ierr = PetscSectionGetChart(leafRankSection, &pStart, &pEnd);CHKERRQ(ierr);
429   ierr = ISGetIndices(leafRanks, &nranks);CHKERRQ(ierr);
430   for (p = pStart; p < pEnd; ++p) {
431     PetscInt ndof, noff, n;
432 
433     ierr = PetscSectionGetDof(leafRankSection, p, &ndof);CHKERRQ(ierr);
434     ierr = PetscSectionGetOffset(leafRankSection, p, &noff);CHKERRQ(ierr);
435     for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);CHKERRQ(ierr);}
436   }
437   ierr = ISRestoreIndices(leafRanks, &nranks);CHKERRQ(ierr);
438   /* Compute leaf-to-root process connectivity */
439   for (l = 0; l < numLeaves; ++l) {PetscBTSet(neighbors, remotePoints[l].rank);}
440   /* Calculate edges */
441   PetscBTClear(neighbors, rank);
442   for(proc = 0, numNeighbors = 0; proc < numProcs; ++proc) {if (PetscBTLookup(neighbors, proc)) ++numNeighbors;}
443   ierr = PetscMalloc1(numNeighbors, &ranksNew);CHKERRQ(ierr);
444   ierr = PetscMalloc1(numNeighbors, &localPointsNew);CHKERRQ(ierr);
445   ierr = PetscMalloc1(numNeighbors, &remotePointsNew);CHKERRQ(ierr);
446   for(proc = 0, n = 0; proc < numProcs; ++proc) {
447     if (PetscBTLookup(neighbors, proc)) {
448       ranksNew[n]              = proc;
449       localPointsNew[n]        = proc;
450       remotePointsNew[n].index = rank;
451       remotePointsNew[n].rank  = proc;
452       ++n;
453     }
454   }
455   ierr = PetscBTDestroy(&neighbors);CHKERRQ(ierr);
456   if (processRanks) {ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), numNeighbors, ranksNew, PETSC_OWN_POINTER, processRanks);CHKERRQ(ierr);}
457   else              {ierr = PetscFree(ranksNew);CHKERRQ(ierr);}
458   if (sfProcess) {
459     ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess);CHKERRQ(ierr);
460     ierr = PetscObjectSetName((PetscObject) *sfProcess, "Two-Sided Process SF");CHKERRQ(ierr);
461     ierr = PetscSFSetFromOptions(*sfProcess);CHKERRQ(ierr);
462     ierr = PetscSFSetGraph(*sfProcess, numProcs, numNeighbors, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr);
463   }
464   PetscFunctionReturn(0);
465 }
466 
467 #undef __FUNCT__
468 #define __FUNCT__ "DMPlexDistributeOwnership"
469 /*@
470   DMPlexDistributeOwnership - Compute owner information for shared points. This basically gets two-sided for an SF.
471 
472   Collective on DM
473 
474   Input Parameter:
475 . dm - The DM
476 
477   Output Parameters:
478 + rootSection - The number of leaves for a given root point
479 . rootrank    - The rank of each edge into the root point
480 . leafSection - The number of processes sharing a given leaf point
481 - leafrank    - The rank of each process sharing a leaf point
482 
483   Level: developer
484 
485 .seealso: DMPlexCreateOverlap()
486 @*/
487 PetscErrorCode DMPlexDistributeOwnership(DM dm, PetscSection rootSection, IS *rootrank, PetscSection leafSection, IS *leafrank)
488 {
489   MPI_Comm        comm;
490   PetscSF         sfPoint;
491   const PetscInt *rootdegree;
492   PetscInt       *myrank, *remoterank;
493   PetscInt        pStart, pEnd, p, nedges;
494   PetscMPIInt     rank;
495   PetscErrorCode  ierr;
496 
497   PetscFunctionBegin;
498   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
499   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
500   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
501   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
502   /* Compute number of leaves for each root */
503   ierr = PetscObjectSetName((PetscObject) rootSection, "Root Section");CHKERRQ(ierr);
504   ierr = PetscSectionSetChart(rootSection, pStart, pEnd);CHKERRQ(ierr);
505   ierr = PetscSFComputeDegreeBegin(sfPoint, &rootdegree);CHKERRQ(ierr);
506   ierr = PetscSFComputeDegreeEnd(sfPoint, &rootdegree);CHKERRQ(ierr);
507   for (p = pStart; p < pEnd; ++p) {ierr = PetscSectionSetDof(rootSection, p, rootdegree[p-pStart]);CHKERRQ(ierr);}
508   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
509   /* Gather rank of each leaf to root */
510   ierr = PetscSectionGetStorageSize(rootSection, &nedges);CHKERRQ(ierr);
511   ierr = PetscMalloc1(pEnd-pStart, &myrank);CHKERRQ(ierr);
512   ierr = PetscMalloc1(nedges,  &remoterank);CHKERRQ(ierr);
513   for (p = 0; p < pEnd-pStart; ++p) myrank[p] = rank;
514   ierr = PetscSFGatherBegin(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr);
515   ierr = PetscSFGatherEnd(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr);
516   ierr = PetscFree(myrank);CHKERRQ(ierr);
517   ierr = ISCreateGeneral(comm, nedges, remoterank, PETSC_OWN_POINTER, rootrank);CHKERRQ(ierr);
518   /* Distribute remote ranks to leaves */
519   ierr = PetscObjectSetName((PetscObject) leafSection, "Leaf Section");CHKERRQ(ierr);
520   ierr = DMPlexDistributeFieldIS(dm, sfPoint, rootSection, *rootrank, leafSection, leafrank);CHKERRQ(ierr);
521   PetscFunctionReturn(0);
522 }
523 
524 #undef __FUNCT__
525 #define __FUNCT__ "DMPlexCreateOverlap"
526 /*@C
527   DMPlexCreateOverlap - Compute owner information for shared points. This basically gets two-sided for an SF.
528 
529   Collective on DM
530 
531   Input Parameters:
532 + dm          - The DM
533 . levels      - Number of overlap levels
534 . rootSection - The number of leaves for a given root point
535 . rootrank    - The rank of each edge into the root point
536 . leafSection - The number of processes sharing a given leaf point
537 - leafrank    - The rank of each process sharing a leaf point
538 
539   Output Parameters:
540 + ovLabel     - DMLabel containing remote overlap contributions as point/rank pairings
541 
542   Level: developer
543 
544 .seealso: DMPlexDistributeOwnership(), DMPlexDistribute()
545 @*/
546 PetscErrorCode DMPlexCreateOverlap(DM dm, PetscInt levels, PetscSection rootSection, IS rootrank, PetscSection leafSection, IS leafrank, DMLabel *ovLabel)
547 {
548   MPI_Comm           comm;
549   DMLabel            ovAdjByRank; /* A DMLabel containing all points adjacent to shared points, separated by rank (value in label) */
550   PetscSF            sfPoint, sfProc;
551   const PetscSFNode *remote;
552   const PetscInt    *local;
553   const PetscInt    *nrank, *rrank;
554   PetscInt          *adj = NULL;
555   PetscInt           pStart, pEnd, p, sStart, sEnd, nleaves, l;
556   PetscMPIInt        rank, numProcs;
557   PetscBool          useCone, useClosure, flg;
558   PetscErrorCode     ierr;
559 
560   PetscFunctionBegin;
561   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
562   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
563   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
564   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
565   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
566   ierr = PetscSectionGetChart(leafSection, &sStart, &sEnd);CHKERRQ(ierr);
567   ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr);
568   ierr = DMLabelCreate("Overlap adjacency", &ovAdjByRank);CHKERRQ(ierr);
569   /* Handle leaves: shared with the root point */
570   for (l = 0; l < nleaves; ++l) {
571     PetscInt adjSize = PETSC_DETERMINE, a;
572 
573     ierr = DMPlexGetAdjacency(dm, local[l], &adjSize, &adj);CHKERRQ(ierr);
574     for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remote[l].rank);CHKERRQ(ierr);}
575   }
576   ierr = ISGetIndices(rootrank, &rrank);CHKERRQ(ierr);
577   ierr = ISGetIndices(leafrank, &nrank);CHKERRQ(ierr);
578   /* Handle roots */
579   for (p = pStart; p < pEnd; ++p) {
580     PetscInt adjSize = PETSC_DETERMINE, neighbors = 0, noff, n, a;
581 
582     if ((p >= sStart) && (p < sEnd)) {
583       /* Some leaves share a root with other leaves on different processes */
584       ierr = PetscSectionGetDof(leafSection, p, &neighbors);CHKERRQ(ierr);
585       if (neighbors) {
586         ierr = PetscSectionGetOffset(leafSection, p, &noff);CHKERRQ(ierr);
587         ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
588         for (n = 0; n < neighbors; ++n) {
589           const PetscInt remoteRank = nrank[noff+n];
590 
591           if (remoteRank == rank) continue;
592           for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);}
593         }
594       }
595     }
596     /* Roots are shared with leaves */
597     ierr = PetscSectionGetDof(rootSection, p, &neighbors);CHKERRQ(ierr);
598     if (!neighbors) continue;
599     ierr = PetscSectionGetOffset(rootSection, p, &noff);CHKERRQ(ierr);
600     ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
601     for (n = 0; n < neighbors; ++n) {
602       const PetscInt remoteRank = rrank[noff+n];
603 
604       if (remoteRank == rank) continue;
605       for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);}
606     }
607   }
608   ierr = PetscFree(adj);CHKERRQ(ierr);
609   ierr = ISRestoreIndices(rootrank, &rrank);CHKERRQ(ierr);
610   ierr = ISRestoreIndices(leafrank, &nrank);CHKERRQ(ierr);
611   /* Add additional overlap levels */
612   for (l = 1; l < levels; l++) {ierr = DMPlexPartitionLabelAdjacency(dm, ovAdjByRank);CHKERRQ(ierr);}
613   /* We require the closure in the overlap */
614   ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
615   ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr);
616   if (useCone || !useClosure) {
617     ierr = DMPlexPartitionLabelClosure(dm, ovAdjByRank);CHKERRQ(ierr);
618   }
619   ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-overlap_view", &flg);CHKERRQ(ierr);
620   if (flg) {
621     ierr = DMLabelView(ovAdjByRank, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
622   }
623   /* Make process SF and invert sender to receiver label */
624   ierr = DMPlexCreateTwoSidedProcessSF(dm, sfPoint, rootSection, rootrank, leafSection, leafrank, NULL, &sfProc);CHKERRQ(ierr);
625   ierr = DMLabelCreate("Overlap label", ovLabel);CHKERRQ(ierr);
626   ierr = DMPlexPartitionLabelInvert(dm, ovAdjByRank, sfProc, *ovLabel);CHKERRQ(ierr);
627   /* Add owned points, except for shared local points */
628   for (p = pStart; p < pEnd; ++p) {ierr = DMLabelSetValue(*ovLabel, p, rank);CHKERRQ(ierr);}
629   for (l = 0; l < nleaves; ++l) {
630     ierr = DMLabelClearValue(*ovLabel, local[l], rank);CHKERRQ(ierr);
631     ierr = DMLabelSetValue(*ovLabel, remote[l].index, remote[l].rank);CHKERRQ(ierr);
632   }
633   /* Clean up */
634   ierr = DMLabelDestroy(&ovAdjByRank);CHKERRQ(ierr);
635   ierr = PetscSFDestroy(&sfProc);CHKERRQ(ierr);
636   PetscFunctionReturn(0);
637 }
638 
639 #undef __FUNCT__
640 #define __FUNCT__ "DMPlexCreateOverlapMigrationSF"
641 PetscErrorCode DMPlexCreateOverlapMigrationSF(DM dm, PetscSF overlapSF, PetscSF *migrationSF)
642 {
643   MPI_Comm           comm;
644   PetscMPIInt        rank, numProcs;
645   PetscInt           d, dim, p, pStart, pEnd, nroots, nleaves, newLeaves, point, numSharedPoints;
646   PetscInt          *pointDepths, *remoteDepths, *ilocal;
647   PetscInt          *depthRecv, *depthShift, *depthIdx;
648   PetscSFNode       *iremote;
649   PetscSF            pointSF;
650   const PetscInt    *sharedLocal;
651   const PetscSFNode *overlapRemote, *sharedRemote;
652   PetscErrorCode     ierr;
653 
654   PetscFunctionBegin;
655   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
656   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
657   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
658   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
659   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
660 
661   /* Before building the migration SF we need to know the new stratum offsets */
662   ierr = PetscSFGetGraph(overlapSF, &nroots, &nleaves, NULL, &overlapRemote);CHKERRQ(ierr);
663   ierr = PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);CHKERRQ(ierr);
664   for (d=0; d<dim+1; d++) {
665     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
666     for (p=pStart; p<pEnd; p++) pointDepths[p] = d;
667   }
668   for (p=0; p<nleaves; p++) remoteDepths[p] = -1;
669   ierr = PetscSFBcastBegin(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
670   ierr = PetscSFBcastEnd(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
671 
672   /* Count recevied points in each stratum and compute the internal strata shift */
673   ierr = PetscMalloc3(dim+1, &depthRecv, dim+1, &depthShift, dim+1, &depthIdx);CHKERRQ(ierr);
674   for (d=0; d<dim+1; d++) depthRecv[d]=0;
675   for (p=0; p<nleaves; p++) depthRecv[remoteDepths[p]]++;
676   depthShift[dim] = 0;
677   for (d=0; d<dim; d++) depthShift[d] = depthRecv[dim];
678   for (d=1; d<dim; d++) depthShift[d] += depthRecv[0];
679   for (d=dim-2; d>0; d--) depthShift[d] += depthRecv[d+1];
680   for (d=0; d<dim+1; d++) {
681     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
682     depthIdx[d] = pStart + depthShift[d];
683   }
684 
685   /* Form the overlap SF build an SF that describes the full overlap migration SF */
686   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
687   newLeaves = pEnd - pStart + nleaves;
688   ierr = PetscMalloc1(newLeaves, &ilocal);CHKERRQ(ierr);
689   ierr = PetscMalloc1(newLeaves, &iremote);CHKERRQ(ierr);
690   /* First map local points to themselves */
691   for (d=0; d<dim+1; d++) {
692     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
693     for (p=pStart; p<pEnd; p++) {
694       point = p + depthShift[d];
695       ilocal[point] = point;
696       iremote[point].index = p;
697       iremote[point].rank = rank;
698       depthIdx[d]++;
699     }
700   }
701 
702   /* Add in the remote roots for currently shared points */
703   ierr = DMGetPointSF(dm, &pointSF);CHKERRQ(ierr);
704   ierr = PetscSFGetGraph(pointSF, NULL, &numSharedPoints, &sharedLocal, &sharedRemote);CHKERRQ(ierr);
705   for (d=0; d<dim+1; d++) {
706     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
707     for (p=0; p<numSharedPoints; p++) {
708       if (pStart <= sharedLocal[p] && sharedLocal[p] < pEnd) {
709         point = sharedLocal[p] + depthShift[d];
710         iremote[point].index = sharedRemote[p].index;
711         iremote[point].rank = sharedRemote[p].rank;
712       }
713     }
714   }
715 
716   /* Now add the incoming overlap points */
717   for (p=0; p<nleaves; p++) {
718     point = depthIdx[remoteDepths[p]];
719     ilocal[point] = point;
720     iremote[point].index = overlapRemote[p].index;
721     iremote[point].rank = overlapRemote[p].rank;
722     depthIdx[remoteDepths[p]]++;
723   }
724   ierr = PetscFree2(pointDepths,remoteDepths);CHKERRQ(ierr);
725 
726   ierr = PetscSFCreate(comm, migrationSF);CHKERRQ(ierr);
727   ierr = PetscObjectSetName((PetscObject) *migrationSF, "Overlap Migration SF");CHKERRQ(ierr);
728   ierr = PetscSFSetFromOptions(*migrationSF);CHKERRQ(ierr);
729   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
730   ierr = PetscSFSetGraph(*migrationSF, pEnd-pStart, newLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
731 
732   ierr = PetscFree3(depthRecv, depthShift, depthIdx);CHKERRQ(ierr);
733   PetscFunctionReturn(0);
734 }
735 
736 #undef __FUNCT__
737 #define __FUNCT__ "DMPlexStratifyMigrationSF"
738 /*@
739   DMPlexStratifyMigrationSF - Add partition overlap to a distributed non-overlapping DM.
740 
741   Input Parameter:
742 + dm          - The DM
743 - sf          - A star forest with non-ordered leaves, usually defining a DM point migration
744 
745   Output Parameter:
746 . migrationSF - A star forest with added leaf indirection that ensures the resulting DM is stratified
747 
748   Level: developer
749 
750 .seealso: DMPlexPartitionLabelCreateSF(), DMPlexDistribute(), DMPlexDistributeOverlap()
751 @*/
752 PetscErrorCode DMPlexStratifyMigrationSF(DM dm, PetscSF sf, PetscSF *migrationSF)
753 {
754   MPI_Comm           comm;
755   PetscMPIInt        rank, numProcs;
756   PetscInt           d, ldepth, depth, p, pStart, pEnd, nroots, nleaves;
757   PetscInt          *pointDepths, *remoteDepths, *ilocal;
758   PetscInt          *depthRecv, *depthShift, *depthIdx;
759   const PetscSFNode *iremote;
760   PetscErrorCode     ierr;
761 
762   PetscFunctionBegin;
763   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
764   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
765   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
766   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
767   ierr = DMPlexGetDepth(dm, &ldepth);CHKERRQ(ierr);
768   ierr = MPIU_Allreduce(&ldepth, &depth, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
769   if ((ldepth >= 0) && (depth != ldepth)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent Plex depth %d != %d", ldepth, depth);
770 
771   /* Before building the migration SF we need to know the new stratum offsets */
772   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, &iremote);CHKERRQ(ierr);
773   ierr = PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);CHKERRQ(ierr);
774   for (d = 0; d < depth+1; ++d) {
775     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
776     for (p = pStart; p < pEnd; ++p) pointDepths[p] = d;
777   }
778   for (p = 0; p < nleaves; ++p) remoteDepths[p] = -1;
779   ierr = PetscSFBcastBegin(sf, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
780   ierr = PetscSFBcastEnd(sf, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
781   /* Count recevied points in each stratum and compute the internal strata shift */
782   ierr = PetscMalloc3(depth+1, &depthRecv, depth+1, &depthShift, depth+1, &depthIdx);CHKERRQ(ierr);
783   for (d = 0; d < depth+1; ++d) depthRecv[d] = 0;
784   for (p = 0; p < nleaves; ++p) depthRecv[remoteDepths[p]]++;
785   depthShift[depth] = 0;
786   for (d = 0; d < depth; ++d) depthShift[d] = depthRecv[depth];
787   for (d = 1; d < depth; ++d) depthShift[d] += depthRecv[0];
788   for (d = depth-2; d > 0; --d) depthShift[d] += depthRecv[d+1];
789   for (d = 0; d < depth+1; ++d) {depthIdx[d] = 0;}
790   /* Derive a new local permutation based on stratified indices */
791   ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr);
792   for (p = 0; p < nleaves; ++p) {
793     const PetscInt dep = remoteDepths[p];
794 
795     ilocal[p] = depthShift[dep] + depthIdx[dep];
796     depthIdx[dep]++;
797   }
798   ierr = PetscSFCreate(comm, migrationSF);CHKERRQ(ierr);
799   ierr = PetscObjectSetName((PetscObject) *migrationSF, "Migration SF");CHKERRQ(ierr);
800   ierr = PetscSFSetGraph(*migrationSF, nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_COPY_VALUES);CHKERRQ(ierr);
801   ierr = PetscFree2(pointDepths,remoteDepths);CHKERRQ(ierr);
802   ierr = PetscFree3(depthRecv, depthShift, depthIdx);CHKERRQ(ierr);
803   PetscFunctionReturn(0);
804 }
805 
806 #undef __FUNCT__
807 #define __FUNCT__ "DMPlexDistributeField"
808 /*@
809   DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
810 
811   Collective on DM
812 
813   Input Parameters:
814 + dm - The DMPlex object
815 . pointSF - The PetscSF describing the communication pattern
816 . originalSection - The PetscSection for existing data layout
817 - originalVec - The existing data
818 
819   Output Parameters:
820 + newSection - The PetscSF describing the new data layout
821 - newVec - The new data
822 
823   Level: developer
824 
825 .seealso: DMPlexDistribute(), DMPlexDistributeFieldIS(), DMPlexDistributeData()
826 @*/
827 PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec)
828 {
829   PetscSF        fieldSF;
830   PetscInt      *remoteOffsets, fieldSize;
831   PetscScalar   *originalValues, *newValues;
832   PetscErrorCode ierr;
833 
834   PetscFunctionBegin;
835   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
836   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
837 
838   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
839   ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr);
840   ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr);
841 
842   ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr);
843   ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr);
844   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
845   ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
846   ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
847   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
848   ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr);
849   ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr);
850   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
851   PetscFunctionReturn(0);
852 }
853 
854 #undef __FUNCT__
855 #define __FUNCT__ "DMPlexDistributeFieldIS"
856 /*@
857   DMPlexDistributeFieldIS - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
858 
859   Collective on DM
860 
861   Input Parameters:
862 + dm - The DMPlex object
863 . pointSF - The PetscSF describing the communication pattern
864 . originalSection - The PetscSection for existing data layout
865 - originalIS - The existing data
866 
867   Output Parameters:
868 + newSection - The PetscSF describing the new data layout
869 - newIS - The new data
870 
871   Level: developer
872 
873 .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexDistributeData()
874 @*/
875 PetscErrorCode DMPlexDistributeFieldIS(DM dm, PetscSF pointSF, PetscSection originalSection, IS originalIS, PetscSection newSection, IS *newIS)
876 {
877   PetscSF         fieldSF;
878   PetscInt       *newValues, *remoteOffsets, fieldSize;
879   const PetscInt *originalValues;
880   PetscErrorCode  ierr;
881 
882   PetscFunctionBegin;
883   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
884   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
885 
886   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
887   ierr = PetscMalloc1(fieldSize, &newValues);CHKERRQ(ierr);
888 
889   ierr = ISGetIndices(originalIS, &originalValues);CHKERRQ(ierr);
890   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
891   ierr = PetscSFBcastBegin(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr);
892   ierr = PetscSFBcastEnd(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr);
893   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
894   ierr = ISRestoreIndices(originalIS, &originalValues);CHKERRQ(ierr);
895   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) pointSF), fieldSize, newValues, PETSC_OWN_POINTER, newIS);CHKERRQ(ierr);
896   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
897   PetscFunctionReturn(0);
898 }
899 
900 #undef __FUNCT__
901 #define __FUNCT__ "DMPlexDistributeData"
902 /*@
903   DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
904 
905   Collective on DM
906 
907   Input Parameters:
908 + dm - The DMPlex object
909 . pointSF - The PetscSF describing the communication pattern
910 . originalSection - The PetscSection for existing data layout
911 . datatype - The type of data
912 - originalData - The existing data
913 
914   Output Parameters:
915 + newSection - The PetscSection describing the new data layout
916 - newData - The new data
917 
918   Level: developer
919 
920 .seealso: DMPlexDistribute(), DMPlexDistributeField()
921 @*/
922 PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData)
923 {
924   PetscSF        fieldSF;
925   PetscInt      *remoteOffsets, fieldSize;
926   PetscMPIInt    dataSize;
927   PetscErrorCode ierr;
928 
929   PetscFunctionBegin;
930   ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
931   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
932 
933   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
934   ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr);
935   ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr);
936 
937   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
938   ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
939   ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
940   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
941   ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
942   PetscFunctionReturn(0);
943 }
944 
945 #undef __FUNCT__
946 #define __FUNCT__ "DMPlexDistributeCones"
947 PetscErrorCode DMPlexDistributeCones(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping original, ISLocalToGlobalMapping renumbering, DM dmParallel)
948 {
949   DM_Plex               *mesh  = (DM_Plex*) dm->data;
950   DM_Plex               *pmesh = (DM_Plex*) (dmParallel)->data;
951   MPI_Comm               comm;
952   PetscSF                coneSF;
953   PetscSection           originalConeSection, newConeSection;
954   PetscInt              *remoteOffsets, *cones, *globCones, *newCones, newConesSize;
955   PetscBool              flg;
956   PetscErrorCode         ierr;
957 
958   PetscFunctionBegin;
959   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
960   PetscValidPointer(dmParallel,4);
961   ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
962 
963   /* Distribute cone section */
964   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
965   ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr);
966   ierr = DMPlexGetConeSection(dmParallel, &newConeSection);CHKERRQ(ierr);
967   ierr = PetscSFDistributeSection(migrationSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr);
968   ierr = DMSetUp(dmParallel);CHKERRQ(ierr);
969   {
970     PetscInt pStart, pEnd, p;
971 
972     ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr);
973     for (p = pStart; p < pEnd; ++p) {
974       PetscInt coneSize;
975       ierr               = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr);
976       pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize);
977     }
978   }
979   /* Communicate and renumber cones */
980   ierr = PetscSFCreateSectionSF(migrationSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr);
981   ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
982   if (original) {
983     PetscInt numCones;
984 
985     ierr = PetscSectionGetStorageSize(originalConeSection,&numCones);CHKERRQ(ierr); ierr = PetscMalloc1(numCones,&globCones);CHKERRQ(ierr);
986     ierr = ISLocalToGlobalMappingApplyBlock(original, numCones, cones, globCones);CHKERRQ(ierr);
987   }
988   else {
989     globCones = cones;
990   }
991   ierr = DMPlexGetCones(dmParallel, &newCones);CHKERRQ(ierr);
992   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, globCones, newCones);CHKERRQ(ierr);
993   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, globCones, newCones);CHKERRQ(ierr);
994   if (original) {
995     ierr = PetscFree(globCones);CHKERRQ(ierr);
996   }
997   ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr);
998   ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr);
999 #if PETSC_USE_DEBUG
1000   {
1001     PetscInt  p;
1002     PetscBool valid = PETSC_TRUE;
1003     for (p = 0; p < newConesSize; ++p) {
1004       if (newCones[p] < 0) {valid = PETSC_FALSE; ierr = PetscPrintf(PETSC_COMM_SELF, "Point %d not in overlap SF\n", p);CHKERRQ(ierr);}
1005     }
1006     if (!valid) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map");
1007   }
1008 #endif
1009   ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr);
1010   if (flg) {
1011     ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr);
1012     ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1013     ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr);
1014     ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1015     ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr);
1016   }
1017   ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr);
1018   ierr = DMPlexGetConeOrientations(dmParallel, &newCones);CHKERRQ(ierr);
1019   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
1020   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
1021   ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr);
1022   ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
1023   /* Create supports and stratify sieve */
1024   {
1025     PetscInt pStart, pEnd;
1026 
1027     ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr);
1028     ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr);
1029   }
1030   ierr = DMPlexSymmetrize(dmParallel);CHKERRQ(ierr);
1031   ierr = DMPlexStratify(dmParallel);CHKERRQ(ierr);
1032   pmesh->useCone    = mesh->useCone;
1033   pmesh->useClosure = mesh->useClosure;
1034   PetscFunctionReturn(0);
1035 }
1036 
1037 #undef __FUNCT__
1038 #define __FUNCT__ "DMPlexDistributeCoordinates"
1039 PetscErrorCode DMPlexDistributeCoordinates(DM dm, PetscSF migrationSF, DM dmParallel)
1040 {
1041   MPI_Comm         comm;
1042   PetscSection     originalCoordSection, newCoordSection;
1043   Vec              originalCoordinates, newCoordinates;
1044   PetscInt         bs;
1045   const char      *name;
1046   const PetscReal *maxCell, *L;
1047   const DMBoundaryType *bd;
1048   PetscErrorCode   ierr;
1049 
1050   PetscFunctionBegin;
1051   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1052   PetscValidPointer(dmParallel, 3);
1053 
1054   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
1055   ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr);
1056   ierr = DMGetCoordinateSection(dmParallel, &newCoordSection);CHKERRQ(ierr);
1057   ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr);
1058   if (originalCoordinates) {
1059     ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr);
1060     ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr);
1061     ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr);
1062 
1063     ierr = DMPlexDistributeField(dm, migrationSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr);
1064     ierr = DMSetCoordinatesLocal(dmParallel, newCoordinates);CHKERRQ(ierr);
1065     ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr);
1066     ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr);
1067     ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr);
1068   }
1069   ierr = DMGetPeriodicity(dm, &maxCell, &L, &bd);CHKERRQ(ierr);
1070   if (L) {ierr = DMSetPeriodicity(dmParallel, maxCell, L, bd);CHKERRQ(ierr);}
1071   PetscFunctionReturn(0);
1072 }
1073 
1074 #undef __FUNCT__
1075 #define __FUNCT__ "DMPlexDistributeLabels"
1076 /* Here we are assuming that process 0 always has everything */
1077 PetscErrorCode DMPlexDistributeLabels(DM dm, PetscSF migrationSF, DM dmParallel)
1078 {
1079   DM_Plex         *mesh = (DM_Plex*) dm->data;
1080   MPI_Comm         comm;
1081   DMLabel          depth;
1082   PetscMPIInt      rank;
1083   PetscInt         numLabels, numLocalLabels, l;
1084   PetscBool        hasLabels = PETSC_FALSE, lsendDepth, sendDepth;
1085   PetscObjectState depthState = -1;
1086   PetscErrorCode   ierr;
1087 
1088   PetscFunctionBegin;
1089   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1090   PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
1091   ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
1092   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
1093   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1094 
1095   /* If the user has changed the depth label, communicate it instead */
1096   ierr = DMPlexGetDepthLabel(dm, &depth);CHKERRQ(ierr);
1097   if (depth) {ierr = DMLabelGetState(depth, &depthState);CHKERRQ(ierr);}
1098   lsendDepth = mesh->depthState != depthState ? PETSC_TRUE : PETSC_FALSE;
1099   ierr = MPIU_Allreduce(&lsendDepth, &sendDepth, 1, MPIU_BOOL, MPI_LOR, comm);CHKERRQ(ierr);
1100   if (sendDepth) {
1101     ierr = DMPlexRemoveLabel(dmParallel, "depth", &depth);CHKERRQ(ierr);
1102     ierr = DMLabelDestroy(&depth);CHKERRQ(ierr);
1103   }
1104   /* Everyone must have either the same number of labels, or none */
1105   ierr = DMPlexGetNumLabels(dm, &numLocalLabels);CHKERRQ(ierr);
1106   numLabels = numLocalLabels;
1107   ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1108   if (numLabels == numLocalLabels) hasLabels = PETSC_TRUE;
1109   for (l = numLabels-1; l >= 0; --l) {
1110     DMLabel     label = NULL, labelNew = NULL;
1111     PetscBool   isdepth;
1112 
1113     if (hasLabels) {
1114       ierr = DMPlexGetLabelByNum(dm, l, &label);CHKERRQ(ierr);
1115       /* Skip "depth" because it is recreated */
1116       ierr = PetscStrcmp(label->name, "depth", &isdepth);CHKERRQ(ierr);
1117     }
1118     ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
1119     if (isdepth && !sendDepth) continue;
1120     ierr = DMLabelDistribute(label, migrationSF, &labelNew);CHKERRQ(ierr);
1121     ierr = DMPlexAddLabel(dmParallel, labelNew);CHKERRQ(ierr);
1122   }
1123   ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
1124   PetscFunctionReturn(0);
1125 }
1126 
1127 #undef __FUNCT__
1128 #define __FUNCT__ "DMPlexDistributeSetupHybrid"
1129 PetscErrorCode DMPlexDistributeSetupHybrid(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
1130 {
1131   DM_Plex        *mesh  = (DM_Plex*) dm->data;
1132   DM_Plex        *pmesh = (DM_Plex*) (dmParallel)->data;
1133   MPI_Comm        comm;
1134   const PetscInt *gpoints;
1135   PetscInt        dim, depth, n, d;
1136   PetscErrorCode  ierr;
1137 
1138   PetscFunctionBegin;
1139   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1140   PetscValidPointer(dmParallel, 4);
1141 
1142   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
1143   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1144 
1145   /* Setup hybrid structure */
1146   for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];}
1147   ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1148   ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr);
1149   ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr);
1150   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1151   for (d = 0; d <= dim; ++d) {
1152     PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p;
1153 
1154     if (pmax < 0) continue;
1155     ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr);
1156     ierr = DMPlexGetDepthStratum(dmParallel, d, NULL, &pEnd);CHKERRQ(ierr);
1157     ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr);
1158     for (p = 0; p < n; ++p) {
1159       const PetscInt point = gpoints[p];
1160 
1161       if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax;
1162     }
1163     if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax;
1164     else            pmesh->hybridPointMax[d] = -1;
1165   }
1166   ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr);
1167   PetscFunctionReturn(0);
1168 }
1169 
1170 #undef __FUNCT__
1171 #define __FUNCT__ "DMPlexDistributeSetupTree"
1172 PetscErrorCode DMPlexDistributeSetupTree(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping original, ISLocalToGlobalMapping renumbering, DM dmParallel)
1173 {
1174   DM_Plex        *mesh  = (DM_Plex*) dm->data;
1175   DM_Plex        *pmesh = (DM_Plex*) (dmParallel)->data;
1176   MPI_Comm        comm;
1177   DM              refTree;
1178   PetscSection    origParentSection, newParentSection;
1179   PetscInt        *origParents, *origChildIDs;
1180   PetscBool       flg;
1181   PetscErrorCode  ierr;
1182 
1183   PetscFunctionBegin;
1184   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1185   PetscValidHeaderSpecific(dm, DM_CLASSID, 4);
1186   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
1187 
1188   /* Set up tree */
1189   ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr);
1190   ierr = DMPlexSetReferenceTree(dmParallel,refTree);CHKERRQ(ierr);
1191   ierr = DMPlexGetTree(dm,&origParentSection,&origParents,&origChildIDs,NULL,NULL);CHKERRQ(ierr);
1192   if (origParentSection) {
1193     PetscInt        pStart, pEnd;
1194     PetscInt        *newParents, *newChildIDs, *globParents;
1195     PetscInt        *remoteOffsetsParents, newParentSize;
1196     PetscSF         parentSF;
1197 
1198     ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
1199     ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dmParallel),&newParentSection);CHKERRQ(ierr);
1200     ierr = PetscSectionSetChart(newParentSection,pStart,pEnd);CHKERRQ(ierr);
1201     ierr = PetscSFDistributeSection(migrationSF, origParentSection, &remoteOffsetsParents, newParentSection);CHKERRQ(ierr);
1202     ierr = PetscSFCreateSectionSF(migrationSF, origParentSection, remoteOffsetsParents, newParentSection, &parentSF);CHKERRQ(ierr);
1203     ierr = PetscSectionGetStorageSize(newParentSection,&newParentSize);CHKERRQ(ierr);
1204     ierr = PetscMalloc2(newParentSize,&newParents,newParentSize,&newChildIDs);CHKERRQ(ierr);
1205     if (original) {
1206       PetscInt numParents;
1207 
1208       ierr = PetscSectionGetStorageSize(origParentSection,&numParents);CHKERRQ(ierr);
1209       ierr = PetscMalloc1(numParents,&globParents);CHKERRQ(ierr);
1210       ierr = ISLocalToGlobalMappingApplyBlock(original, numParents, origParents, globParents);CHKERRQ(ierr);
1211     }
1212     else {
1213       globParents = origParents;
1214     }
1215     ierr = PetscSFBcastBegin(parentSF, MPIU_INT, globParents, newParents);CHKERRQ(ierr);
1216     ierr = PetscSFBcastEnd(parentSF, MPIU_INT, globParents, newParents);CHKERRQ(ierr);
1217     if (original) {
1218       ierr = PetscFree(globParents);CHKERRQ(ierr);
1219     }
1220     ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1221     ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1222     ierr = ISGlobalToLocalMappingApplyBlock(renumbering,IS_GTOLM_MASK, newParentSize, newParents, NULL, newParents);CHKERRQ(ierr);
1223 #if PETSC_USE_DEBUG
1224     {
1225       PetscInt  p;
1226       PetscBool valid = PETSC_TRUE;
1227       for (p = 0; p < newParentSize; ++p) {
1228         if (newParents[p] < 0) {valid = PETSC_FALSE; ierr = PetscPrintf(PETSC_COMM_SELF, "Point %d not in overlap SF\n", p);CHKERRQ(ierr);}
1229       }
1230       if (!valid) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map");
1231     }
1232 #endif
1233     ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-parents_view", &flg);CHKERRQ(ierr);
1234     if (flg) {
1235       ierr = PetscPrintf(comm, "Serial Parent Section: \n");CHKERRQ(ierr);
1236       ierr = PetscSectionView(origParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1237       ierr = PetscPrintf(comm, "Parallel Parent Section: \n");CHKERRQ(ierr);
1238       ierr = PetscSectionView(newParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1239       ierr = PetscSFView(parentSF, NULL);CHKERRQ(ierr);
1240     }
1241     ierr = DMPlexSetTree(dmParallel,newParentSection,newParents,newChildIDs);CHKERRQ(ierr);
1242     ierr = PetscSectionDestroy(&newParentSection);CHKERRQ(ierr);
1243     ierr = PetscFree2(newParents,newChildIDs);CHKERRQ(ierr);
1244     ierr = PetscSFDestroy(&parentSF);CHKERRQ(ierr);
1245   }
1246   pmesh->useAnchors = mesh->useAnchors;
1247   PetscFunctionReturn(0);
1248 }
1249 
1250 #undef __FUNCT__
1251 #define __FUNCT__ "DMPlexDistributeSF"
1252 PetscErrorCode DMPlexDistributeSF(DM dm, PetscSF migrationSF, DM dmParallel)
1253 {
1254   DM_Plex               *mesh  = (DM_Plex*) dm->data;
1255   DM_Plex               *pmesh = (DM_Plex*) (dmParallel)->data;
1256   PetscMPIInt            rank, numProcs;
1257   MPI_Comm               comm;
1258   PetscErrorCode         ierr;
1259 
1260   PetscFunctionBegin;
1261   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1262   PetscValidPointer(dmParallel,7);
1263 
1264   /* Create point SF for parallel mesh */
1265   ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
1266   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
1267   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1268   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
1269   {
1270     const PetscInt *leaves;
1271     PetscSFNode    *remotePoints, *rowners, *lowners;
1272     PetscInt        numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints;
1273     PetscInt        pStart, pEnd;
1274 
1275     ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
1276     ierr = PetscSFGetGraph(migrationSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr);
1277     ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr);
1278     for (p=0; p<numRoots; p++) {
1279       rowners[p].rank  = -1;
1280       rowners[p].index = -1;
1281     }
1282     ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
1283     ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
1284     for (p = 0; p < numLeaves; ++p) {
1285       if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */
1286         lowners[p].rank  = rank;
1287         lowners[p].index = leaves ? leaves[p] : p;
1288       } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */
1289         lowners[p].rank  = -2;
1290         lowners[p].index = -2;
1291       }
1292     }
1293     for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */
1294       rowners[p].rank  = -3;
1295       rowners[p].index = -3;
1296     }
1297     ierr = PetscSFReduceBegin(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
1298     ierr = PetscSFReduceEnd(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
1299     ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
1300     ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
1301     for (p = 0; p < numLeaves; ++p) {
1302       if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed");
1303       if (lowners[p].rank != rank) ++numGhostPoints;
1304     }
1305     ierr = PetscMalloc1(numGhostPoints, &ghostPoints);CHKERRQ(ierr);
1306     ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr);
1307     for (p = 0, gp = 0; p < numLeaves; ++p) {
1308       if (lowners[p].rank != rank) {
1309         ghostPoints[gp]        = leaves ? leaves[p] : p;
1310         remotePoints[gp].rank  = lowners[p].rank;
1311         remotePoints[gp].index = lowners[p].index;
1312         ++gp;
1313       }
1314     }
1315     ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr);
1316     ierr = PetscSFSetGraph((dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
1317     ierr = PetscSFSetFromOptions((dmParallel)->sf);CHKERRQ(ierr);
1318   }
1319   pmesh->useCone    = mesh->useCone;
1320   pmesh->useClosure = mesh->useClosure;
1321   ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
1322   PetscFunctionReturn(0);
1323 }
1324 
1325 #undef __FUNCT__
1326 #define __FUNCT__ "DMPlexCreatePointSF"
1327 /*@C
1328   DMPlexDerivePointSF - Build a point SF from an SF describing a point migration
1329 
1330   Input Parameter:
1331 + dm          - The source DMPlex object
1332 . migrationSF - The star forest that describes the parallel point remapping
1333 . ownership   - Flag causing a vote to determine point ownership
1334 
1335   Output Parameter:
1336 - pointSF     - The star forest describing the point overlap in the remapped DM
1337 
1338   Level: developer
1339 
1340 .seealso: DMPlexDistribute(), DMPlexDistributeOverlap()
1341 @*/
1342 PetscErrorCode DMPlexCreatePointSF(DM dm, PetscSF migrationSF, PetscBool ownership, PetscSF *pointSF)
1343 {
1344   PetscMPIInt        rank;
1345   PetscInt           p, nroots, nleaves, idx, npointLeaves;
1346   PetscInt          *pointLocal;
1347   const PetscInt    *leaves;
1348   const PetscSFNode *roots;
1349   PetscSFNode       *rootNodes, *leafNodes, *pointRemote;
1350   PetscErrorCode     ierr;
1351 
1352   PetscFunctionBegin;
1353   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1354   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
1355 
1356   ierr = PetscSFGetGraph(migrationSF, &nroots, &nleaves, &leaves, &roots);CHKERRQ(ierr);
1357   ierr = PetscMalloc2(nroots, &rootNodes, nleaves, &leafNodes);CHKERRQ(ierr);
1358   if (ownership) {
1359     /* Point ownership vote: Process with highest rank ownes shared points */
1360     for (p = 0; p < nleaves; ++p) {
1361       /* Either put in a bid or we know we own it */
1362       leafNodes[p].rank  = rank;
1363       leafNodes[p].index = p;
1364     }
1365     for (p = 0; p < nroots; p++) {
1366       /* Root must not participate in the reduction, flag so that MAXLOC does not use */
1367       rootNodes[p].rank  = -3;
1368       rootNodes[p].index = -3;
1369     }
1370     ierr = PetscSFReduceBegin(migrationSF, MPIU_2INT, leafNodes, rootNodes, MPI_MAXLOC);CHKERRQ(ierr);
1371     ierr = PetscSFReduceEnd(migrationSF, MPIU_2INT, leafNodes, rootNodes, MPI_MAXLOC);CHKERRQ(ierr);
1372   } else {
1373     for (p = 0; p < nroots; p++) {
1374       rootNodes[p].index = -1;
1375       rootNodes[p].rank = rank;
1376     };
1377     for (p = 0; p < nleaves; p++) {
1378       /* Write new local id into old location */
1379       if (roots[p].rank == rank) {
1380         rootNodes[roots[p].index].index = leaves[p];
1381       }
1382     }
1383   }
1384   ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rootNodes, leafNodes);CHKERRQ(ierr);
1385   ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rootNodes, leafNodes);CHKERRQ(ierr);
1386 
1387   for (npointLeaves = 0, p = 0; p < nleaves; p++) {if (leafNodes[p].rank != rank) npointLeaves++;}
1388   ierr = PetscMalloc1(npointLeaves, &pointLocal);CHKERRQ(ierr);
1389   ierr = PetscMalloc1(npointLeaves, &pointRemote);CHKERRQ(ierr);
1390   for (idx = 0, p = 0; p < nleaves; p++) {
1391     if (leafNodes[p].rank != rank) {
1392       pointLocal[idx] = p;
1393       pointRemote[idx] = leafNodes[p];
1394       idx++;
1395     }
1396   }
1397   ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), pointSF);CHKERRQ(ierr);
1398   ierr = PetscSFSetFromOptions(*pointSF);CHKERRQ(ierr);
1399   ierr = PetscSFSetGraph(*pointSF, nleaves, npointLeaves, pointLocal, PETSC_OWN_POINTER, pointRemote, PETSC_OWN_POINTER);CHKERRQ(ierr);
1400   ierr = PetscFree2(rootNodes, leafNodes);CHKERRQ(ierr);
1401   PetscFunctionReturn(0);
1402 }
1403 
1404 #undef __FUNCT__
1405 #define __FUNCT__ "DMPlexMigrate"
1406 /*@C
1407   DMPlexMigrate  - Migrates internal DM data over the supplied star forest
1408 
1409   Input Parameter:
1410 + dm       - The source DMPlex object
1411 . sf       - The star forest communication context describing the migration pattern
1412 
1413   Output Parameter:
1414 - targetDM - The target DMPlex object
1415 
1416   Level: intermediate
1417 
1418 .seealso: DMPlexDistribute(), DMPlexDistributeOverlap()
1419 @*/
1420 PetscErrorCode DMPlexMigrate(DM dm, PetscSF sf, DM targetDM)
1421 {
1422   MPI_Comm               comm;
1423   PetscInt               dim, nroots;
1424   PetscSF                sfPoint;
1425   ISLocalToGlobalMapping ltogMigration;
1426   ISLocalToGlobalMapping ltogOriginal = NULL;
1427   PetscBool              flg;
1428   PetscErrorCode         ierr;
1429 
1430   PetscFunctionBegin;
1431   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1432   ierr = PetscLogEventBegin(DMPLEX_Migrate, dm, 0, 0, 0);CHKERRQ(ierr);
1433   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
1434   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1435   ierr = DMSetDimension(targetDM, dim);CHKERRQ(ierr);
1436 
1437   /* Check for a one-to-all distribution pattern */
1438   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
1439   ierr = PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1440   if (nroots >= 0) {
1441     IS                     isOriginal;
1442     PetscInt               n, size, nleaves;
1443     PetscInt              *numbering_orig, *numbering_new;
1444     /* Get the original point numbering */
1445     ierr = DMPlexCreatePointNumbering(dm, &isOriginal);CHKERRQ(ierr);
1446     ierr = ISLocalToGlobalMappingCreateIS(isOriginal, &ltogOriginal);CHKERRQ(ierr);
1447     ierr = ISLocalToGlobalMappingGetSize(ltogOriginal, &size);CHKERRQ(ierr);
1448     ierr = ISLocalToGlobalMappingGetBlockIndices(ltogOriginal, (const PetscInt**)&numbering_orig);CHKERRQ(ierr);
1449     /* Convert to positive global numbers */
1450     for (n=0; n<size; n++) {if (numbering_orig[n] < 0) numbering_orig[n] = -(numbering_orig[n]+1);}
1451     /* Derive the new local-to-global mapping from the old one */
1452     ierr = PetscSFGetGraph(sf, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr);
1453     ierr = PetscMalloc1(nleaves, &numbering_new);CHKERRQ(ierr);
1454     ierr = PetscSFBcastBegin(sf, MPIU_INT, (PetscInt *) numbering_orig, numbering_new);CHKERRQ(ierr);
1455     ierr = PetscSFBcastEnd(sf, MPIU_INT, (PetscInt *) numbering_orig, numbering_new);CHKERRQ(ierr);
1456     ierr = ISLocalToGlobalMappingCreate(comm, 1, nleaves, (const PetscInt*) numbering_new, PETSC_OWN_POINTER, &ltogMigration);CHKERRQ(ierr);
1457     ierr = ISLocalToGlobalMappingRestoreIndices(ltogOriginal, (const PetscInt**)&numbering_orig);CHKERRQ(ierr);
1458     ierr = ISDestroy(&isOriginal);CHKERRQ(ierr);
1459   } else {
1460     /* One-to-all distribution pattern: We can derive LToG from SF */
1461     ierr = ISLocalToGlobalMappingCreateSF(sf, 0, &ltogMigration);CHKERRQ(ierr);
1462   }
1463   ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
1464   if (flg) {
1465     ierr = PetscPrintf(comm, "Point renumbering for DM migration:\n");CHKERRQ(ierr);
1466     ierr = ISLocalToGlobalMappingView(ltogMigration, NULL);CHKERRQ(ierr);
1467   }
1468   /* Migrate DM data to target DM */
1469   ierr = DMPlexDistributeCones(dm, sf, ltogOriginal, ltogMigration, targetDM);CHKERRQ(ierr);
1470   ierr = DMPlexDistributeCoordinates(dm, sf, targetDM);CHKERRQ(ierr);
1471   ierr = DMPlexDistributeLabels(dm, sf, targetDM);CHKERRQ(ierr);
1472   ierr = DMPlexDistributeSetupHybrid(dm, sf, ltogMigration, targetDM);CHKERRQ(ierr);
1473   ierr = DMPlexDistributeSetupTree(dm, sf, ltogOriginal, ltogMigration, targetDM);CHKERRQ(ierr);
1474   ierr = ISLocalToGlobalMappingDestroy(&ltogOriginal);CHKERRQ(ierr);
1475   ierr = ISLocalToGlobalMappingDestroy(&ltogMigration);CHKERRQ(ierr);
1476   ierr = PetscLogEventEnd(DMPLEX_Migrate, dm, 0, 0, 0);CHKERRQ(ierr);
1477   PetscFunctionReturn(0);
1478 }
1479 
1480 #undef __FUNCT__
1481 #define __FUNCT__ "DMPlexDistribute"
1482 /*@C
1483   DMPlexDistribute - Distributes the mesh and any associated sections.
1484 
1485   Not Collective
1486 
1487   Input Parameter:
1488 + dm  - The original DMPlex object
1489 - overlap - The overlap of partitions, 0 is the default
1490 
1491   Output Parameter:
1492 + sf - The PetscSF used for point distribution
1493 - parallelMesh - The distributed DMPlex object, or NULL
1494 
1495   Note: If the mesh was not distributed, the return value is NULL.
1496 
1497   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
1498   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
1499   representation on the mesh.
1500 
1501   Level: intermediate
1502 
1503 .keywords: mesh, elements
1504 .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
1505 @*/
1506 PetscErrorCode DMPlexDistribute(DM dm, PetscInt overlap, PetscSF *sf, DM *dmParallel)
1507 {
1508   MPI_Comm               comm;
1509   PetscPartitioner       partitioner;
1510   IS                     cellPart;
1511   PetscSection           cellPartSection;
1512   DM                     dmCoord;
1513   DMLabel                lblPartition, lblMigration;
1514   PetscSF                sfProcess, sfMigration, sfStratified, sfPoint;
1515   PetscBool              flg;
1516   PetscMPIInt            rank, numProcs, p;
1517   PetscErrorCode         ierr;
1518 
1519   PetscFunctionBegin;
1520   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1521   if (sf) PetscValidPointer(sf,4);
1522   PetscValidPointer(dmParallel,5);
1523 
1524   ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
1525   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1526   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1527   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
1528 
1529   *dmParallel = NULL;
1530   if (numProcs == 1) PetscFunctionReturn(0);
1531 
1532   /* Create cell partition */
1533   ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1534   ierr = PetscSectionCreate(comm, &cellPartSection);CHKERRQ(ierr);
1535   ierr = DMPlexGetPartitioner(dm, &partitioner);CHKERRQ(ierr);
1536   ierr = PetscPartitionerPartition(partitioner, dm, cellPartSection, &cellPart);CHKERRQ(ierr);
1537   {
1538     /* Convert partition to DMLabel */
1539     PetscInt proc, pStart, pEnd, npoints, poffset;
1540     const PetscInt *points;
1541     ierr = DMLabelCreate("Point Partition", &lblPartition);CHKERRQ(ierr);
1542     ierr = ISGetIndices(cellPart, &points);CHKERRQ(ierr);
1543     ierr = PetscSectionGetChart(cellPartSection, &pStart, &pEnd);CHKERRQ(ierr);
1544     for (proc = pStart; proc < pEnd; proc++) {
1545       ierr = PetscSectionGetDof(cellPartSection, proc, &npoints);CHKERRQ(ierr);
1546       ierr = PetscSectionGetOffset(cellPartSection, proc, &poffset);CHKERRQ(ierr);
1547       for (p = poffset; p < poffset+npoints; p++) {
1548         ierr = DMLabelSetValue(lblPartition, points[p], proc);CHKERRQ(ierr);
1549       }
1550     }
1551     ierr = ISRestoreIndices(cellPart, &points);CHKERRQ(ierr);
1552   }
1553   ierr = DMPlexPartitionLabelClosure(dm, lblPartition);CHKERRQ(ierr);
1554   {
1555     /* Build a global process SF */
1556     PetscSFNode *remoteProc;
1557     ierr = PetscMalloc1(numProcs, &remoteProc);CHKERRQ(ierr);
1558     for (p = 0; p < numProcs; ++p) {
1559       remoteProc[p].rank  = p;
1560       remoteProc[p].index = rank;
1561     }
1562     ierr = PetscSFCreate(comm, &sfProcess);CHKERRQ(ierr);
1563     ierr = PetscObjectSetName((PetscObject) sfProcess, "Process SF");CHKERRQ(ierr);
1564     ierr = PetscSFSetGraph(sfProcess, numProcs, numProcs, NULL, PETSC_OWN_POINTER, remoteProc, PETSC_OWN_POINTER);CHKERRQ(ierr);
1565   }
1566   ierr = DMLabelCreate("Point migration", &lblMigration);CHKERRQ(ierr);
1567   ierr = DMPlexPartitionLabelInvert(dm, lblPartition, sfProcess, lblMigration);CHKERRQ(ierr);
1568   ierr = DMPlexPartitionLabelCreateSF(dm, lblMigration, &sfMigration);CHKERRQ(ierr);
1569   /* Stratify the SF in case we are migrating an already parallel plex */
1570   ierr = DMPlexStratifyMigrationSF(dm, sfMigration, &sfStratified);CHKERRQ(ierr);
1571   ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr);
1572   sfMigration = sfStratified;
1573   ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1574   ierr = PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
1575   if (flg) {
1576     ierr = DMLabelView(lblPartition, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1577     ierr = PetscSFView(sfMigration, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1578   }
1579 
1580   /* Create non-overlapping parallel DM and migrate internal data */
1581   ierr = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr);
1582   ierr = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr);
1583   ierr = DMPlexMigrate(dm, sfMigration, *dmParallel);CHKERRQ(ierr);
1584 
1585   /* Build the point SF without overlap */
1586   ierr = DMPlexCreatePointSF(*dmParallel, sfMigration, PETSC_TRUE, &sfPoint);CHKERRQ(ierr);
1587   ierr = DMSetPointSF(*dmParallel, sfPoint);CHKERRQ(ierr);
1588   ierr = DMGetCoordinateDM(*dmParallel, &dmCoord);CHKERRQ(ierr);
1589   if (dmCoord) {ierr = DMSetPointSF(dmCoord, sfPoint);CHKERRQ(ierr);}
1590   if (flg) {ierr = PetscSFView(sfPoint, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
1591 
1592   if (overlap > 0) {
1593     DM                 dmOverlap;
1594     PetscInt           nroots, nleaves;
1595     PetscSFNode       *newRemote;
1596     const PetscSFNode *oldRemote;
1597     PetscSF            sfOverlap, sfOverlapPoint;
1598     /* Add the partition overlap to the distributed DM */
1599     ierr = DMPlexDistributeOverlap(*dmParallel, overlap, &sfOverlap, &dmOverlap);CHKERRQ(ierr);
1600     ierr = DMDestroy(dmParallel);CHKERRQ(ierr);
1601     *dmParallel = dmOverlap;
1602     if (flg) {
1603       ierr = PetscPrintf(comm, "Overlap Migration SF:\n");CHKERRQ(ierr);
1604       ierr = PetscSFView(sfOverlap, NULL);CHKERRQ(ierr);
1605     }
1606 
1607     /* Re-map the migration SF to establish the full migration pattern */
1608     ierr = PetscSFGetGraph(sfMigration, &nroots, NULL, NULL, &oldRemote);CHKERRQ(ierr);
1609     ierr = PetscSFGetGraph(sfOverlap, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr);
1610     ierr = PetscMalloc1(nleaves, &newRemote);CHKERRQ(ierr);
1611     ierr = PetscSFBcastBegin(sfOverlap, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr);
1612     ierr = PetscSFBcastEnd(sfOverlap, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr);
1613     ierr = PetscSFCreate(comm, &sfOverlapPoint);CHKERRQ(ierr);
1614     ierr = PetscSFSetGraph(sfOverlapPoint, nroots, nleaves, NULL, PETSC_OWN_POINTER, newRemote, PETSC_OWN_POINTER);CHKERRQ(ierr);
1615     ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr);
1616     ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr);
1617     sfMigration = sfOverlapPoint;
1618   }
1619   /* Cleanup Partition */
1620   ierr = PetscSFDestroy(&sfProcess);CHKERRQ(ierr);
1621   ierr = DMLabelDestroy(&lblPartition);CHKERRQ(ierr);
1622   ierr = DMLabelDestroy(&lblMigration);CHKERRQ(ierr);
1623   ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr);
1624   ierr = ISDestroy(&cellPart);CHKERRQ(ierr);
1625   /* Copy BC */
1626   ierr = DMPlexCopyBoundary(dm, *dmParallel);CHKERRQ(ierr);
1627   /* Create sfNatural */
1628   if (dm->useNatural) {
1629     PetscSection section;
1630 
1631     ierr = DMGetDefaultSection(dm, &section);CHKERRQ(ierr);
1632     ierr = DMPlexCreateGlobalToNaturalSF(*dmParallel, section, sfMigration, &(*dmParallel)->sfNatural);CHKERRQ(ierr);
1633   }
1634   /* Cleanup */
1635   if (sf) {*sf = sfMigration;}
1636   else    {ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr);}
1637   ierr = PetscSFDestroy(&sfPoint);CHKERRQ(ierr);
1638   ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
1639   PetscFunctionReturn(0);
1640 }
1641 
1642 #undef __FUNCT__
1643 #define __FUNCT__ "DMPlexDistributeOverlap"
1644 /*@C
1645   DMPlexDistribute - Add partition overlap to a distributed non-overlapping DM.
1646 
1647   Not Collective
1648 
1649   Input Parameter:
1650 + dm  - The non-overlapping distrbuted DMPlex object
1651 - overlap - The overlap of partitions, 0 is the default
1652 
1653   Output Parameter:
1654 + sf - The PetscSF used for point distribution
1655 - dmOverlap - The overlapping distributed DMPlex object, or NULL
1656 
1657   Note: If the mesh was not distributed, the return value is NULL.
1658 
1659   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
1660   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
1661   representation on the mesh.
1662 
1663   Level: intermediate
1664 
1665 .keywords: mesh, elements
1666 .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
1667 @*/
1668 PetscErrorCode DMPlexDistributeOverlap(DM dm, PetscInt overlap, PetscSF *sf, DM *dmOverlap)
1669 {
1670   MPI_Comm               comm;
1671   PetscMPIInt            rank;
1672   PetscSection           rootSection, leafSection;
1673   IS                     rootrank, leafrank;
1674   DM                     dmCoord;
1675   DMLabel                lblOverlap;
1676   PetscSF                sfOverlap, sfStratified, sfPoint;
1677   PetscErrorCode         ierr;
1678 
1679   PetscFunctionBegin;
1680   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1681   if (sf) PetscValidPointer(sf, 3);
1682   PetscValidPointer(dmOverlap, 4);
1683 
1684   ierr = PetscLogEventBegin(DMPLEX_DistributeOverlap, dm, 0, 0, 0);CHKERRQ(ierr);
1685   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1686   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1687 
1688   /* Compute point overlap with neighbouring processes on the distributed DM */
1689   ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1690   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
1691   ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr);
1692   ierr = DMPlexDistributeOwnership(dm, rootSection, &rootrank, leafSection, &leafrank);CHKERRQ(ierr);
1693   ierr = DMPlexCreateOverlap(dm, overlap, rootSection, rootrank, leafSection, leafrank, &lblOverlap);CHKERRQ(ierr);
1694   /* Convert overlap label to stratified migration SF */
1695   ierr = DMPlexPartitionLabelCreateSF(dm, lblOverlap, &sfOverlap);CHKERRQ(ierr);
1696   ierr = DMPlexStratifyMigrationSF(dm, sfOverlap, &sfStratified);CHKERRQ(ierr);
1697   ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr);
1698   sfOverlap = sfStratified;
1699   ierr = PetscObjectSetName((PetscObject) sfOverlap, "Overlap SF");CHKERRQ(ierr);
1700   ierr = PetscSFSetFromOptions(sfOverlap);CHKERRQ(ierr);
1701 
1702   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
1703   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
1704   ierr = ISDestroy(&rootrank);CHKERRQ(ierr);
1705   ierr = ISDestroy(&leafrank);CHKERRQ(ierr);
1706   ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1707 
1708   /* Build the overlapping DM */
1709   ierr = DMPlexCreate(comm, dmOverlap);CHKERRQ(ierr);
1710   ierr = PetscObjectSetName((PetscObject) *dmOverlap, "Parallel Mesh");CHKERRQ(ierr);
1711   ierr = DMPlexMigrate(dm, sfOverlap, *dmOverlap);CHKERRQ(ierr);
1712   /* Build the new point SF */
1713   ierr = DMPlexCreatePointSF(*dmOverlap, sfOverlap, PETSC_FALSE, &sfPoint);CHKERRQ(ierr);
1714   ierr = DMSetPointSF(*dmOverlap, sfPoint);CHKERRQ(ierr);
1715   ierr = DMGetCoordinateDM(*dmOverlap, &dmCoord);CHKERRQ(ierr);
1716   if (dmCoord) {ierr = DMSetPointSF(dmCoord, sfPoint);CHKERRQ(ierr);}
1717   ierr = PetscSFDestroy(&sfPoint);CHKERRQ(ierr);
1718   /* Cleanup overlap partition */
1719   ierr = DMLabelDestroy(&lblOverlap);CHKERRQ(ierr);
1720   if (sf) *sf = sfOverlap;
1721   else    {ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr);}
1722   ierr = PetscLogEventEnd(DMPLEX_DistributeOverlap, dm, 0, 0, 0);CHKERRQ(ierr);
1723   PetscFunctionReturn(0);
1724 }
1725