xref: /petsc/src/dm/impls/plex/plexdistribute.c (revision cc71bff1308b5c843157f5fa098d26a9289e3ffb)
1 #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 #include <petscsf.h>
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "DMPlexSetAdjacencyUseCone"
6 /*@
7   DMPlexSetAdjacencyUseCone - Define adjacency in the mesh using either the cone or the support first
8 
9   Input Parameters:
10 + dm      - The DM object
11 - useCone - Flag to use the cone first
12 
13   Level: intermediate
14 
15   Notes:
16 $     FEM:   Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE
17 $     FVM:   Two points p and q are adjacent if q \in star(cone(p)),    useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
18 $     FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
19 
20 .seealso: DMPlexGetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
21 @*/
22 PetscErrorCode DMPlexSetAdjacencyUseCone(DM dm, PetscBool useCone)
23 {
24   DM_Plex *mesh = (DM_Plex *) dm->data;
25 
26   PetscFunctionBegin;
27   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
28   mesh->useCone = useCone;
29   PetscFunctionReturn(0);
30 }
31 
32 #undef __FUNCT__
33 #define __FUNCT__ "DMPlexGetAdjacencyUseCone"
34 /*@
35   DMPlexGetAdjacencyUseCone - Query whether adjacency in the mesh uses the cone or the support first
36 
37   Input Parameter:
38 . dm      - The DM object
39 
40   Output Parameter:
41 . useCone - Flag to use the cone first
42 
43   Level: intermediate
44 
45   Notes:
46 $     FEM:   Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE
47 $     FVM:   Two points p and q are adjacent if q \in star(cone(p)),    useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
48 $     FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
49 
50 .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
51 @*/
52 PetscErrorCode DMPlexGetAdjacencyUseCone(DM dm, PetscBool *useCone)
53 {
54   DM_Plex *mesh = (DM_Plex *) dm->data;
55 
56   PetscFunctionBegin;
57   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
58   PetscValidIntPointer(useCone, 2);
59   *useCone = mesh->useCone;
60   PetscFunctionReturn(0);
61 }
62 
63 #undef __FUNCT__
64 #define __FUNCT__ "DMPlexSetAdjacencyUseClosure"
65 /*@
66   DMPlexSetAdjacencyUseClosure - Define adjacency in the mesh using the transitive closure
67 
68   Input Parameters:
69 + dm      - The DM object
70 - useClosure - Flag to use the closure
71 
72   Level: intermediate
73 
74   Notes:
75 $     FEM:   Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE
76 $     FVM:   Two points p and q are adjacent if q \in star(cone(p)),    useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
77 $     FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
78 
79 .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
80 @*/
81 PetscErrorCode DMPlexSetAdjacencyUseClosure(DM dm, PetscBool useClosure)
82 {
83   DM_Plex *mesh = (DM_Plex *) dm->data;
84 
85   PetscFunctionBegin;
86   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
87   mesh->useClosure = useClosure;
88   PetscFunctionReturn(0);
89 }
90 
91 #undef __FUNCT__
92 #define __FUNCT__ "DMPlexGetAdjacencyUseClosure"
93 /*@
94   DMPlexGetAdjacencyUseClosure - Query whether adjacency in the mesh uses the transitive closure
95 
96   Input Parameter:
97 . dm      - The DM object
98 
99   Output Parameter:
100 . useClosure - Flag to use the closure
101 
102   Level: intermediate
103 
104   Notes:
105 $     FEM:   Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE
106 $     FVM:   Two points p and q are adjacent if q \in star(cone(p)),    useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
107 $     FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
108 
109 .seealso: DMPlexSetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
110 @*/
111 PetscErrorCode DMPlexGetAdjacencyUseClosure(DM dm, PetscBool *useClosure)
112 {
113   DM_Plex *mesh = (DM_Plex *) dm->data;
114 
115   PetscFunctionBegin;
116   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
117   PetscValidIntPointer(useClosure, 2);
118   *useClosure = mesh->useClosure;
119   PetscFunctionReturn(0);
120 }
121 
122 #undef __FUNCT__
123 #define __FUNCT__ "DMPlexSetAdjacencyUseAnchors"
124 /*@
125   DMPlexSetAdjacencyUseAnchors - Define adjacency in the mesh using the point-to-point constraints.
126 
127   Input Parameters:
128 + dm      - The DM object
129 - useAnchors - Flag to use the constraints.  If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place.
130 
131   Level: intermediate
132 
133 .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
134 @*/
135 PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM dm, PetscBool useAnchors)
136 {
137   DM_Plex *mesh = (DM_Plex *) dm->data;
138 
139   PetscFunctionBegin;
140   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
141   mesh->useAnchors = useAnchors;
142   PetscFunctionReturn(0);
143 }
144 
145 #undef __FUNCT__
146 #define __FUNCT__ "DMPlexGetAdjacencyUseAnchors"
147 /*@
148   DMPlexGetAdjacencyUseAnchors - Query whether adjacency in the mesh uses the point-to-point constraints.
149 
150   Input Parameter:
151 . dm      - The DM object
152 
153   Output Parameter:
154 . useAnchors - Flag to use the closure.  If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place.
155 
156   Level: intermediate
157 
158 .seealso: DMPlexSetAdjacencyUseAnchors(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
159 @*/
160 PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM dm, PetscBool *useAnchors)
161 {
162   DM_Plex *mesh = (DM_Plex *) dm->data;
163 
164   PetscFunctionBegin;
165   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
166   PetscValidIntPointer(useAnchors, 2);
167   *useAnchors = mesh->useAnchors;
168   PetscFunctionReturn(0);
169 }
170 
171 #undef __FUNCT__
172 #define __FUNCT__ "DMPlexGetAdjacency_Cone_Internal"
173 static PetscErrorCode DMPlexGetAdjacency_Cone_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
174 {
175   const PetscInt *cone = NULL;
176   PetscInt        numAdj = 0, maxAdjSize = *adjSize, coneSize, c;
177   PetscErrorCode  ierr;
178 
179   PetscFunctionBeginHot;
180   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
181   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
182   for (c = 0; c < coneSize; ++c) {
183     const PetscInt *support = NULL;
184     PetscInt        supportSize, s, q;
185 
186     ierr = DMPlexGetSupportSize(dm, cone[c], &supportSize);CHKERRQ(ierr);
187     ierr = DMPlexGetSupport(dm, cone[c], &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 *cone = NULL;
212     PetscInt        coneSize, c, q;
213 
214     ierr = DMPlexGetConeSize(dm, support[s], &coneSize);CHKERRQ(ierr);
215     ierr = DMPlexGetCone(dm, support[s], &cone);CHKERRQ(ierr);
216     for (c = 0; c < coneSize; ++c) {
217       for (q = 0; q < numAdj || (adj[numAdj++] = cone[c],0); ++q) {
218         if (cone[c] == adj[q]) break;
219       }
220       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
221     }
222   }
223   *adjSize = numAdj;
224   PetscFunctionReturn(0);
225 }
226 
227 #undef __FUNCT__
228 #define __FUNCT__ "DMPlexGetAdjacency_Transitive_Internal"
229 static PetscErrorCode DMPlexGetAdjacency_Transitive_Internal(DM dm, PetscInt p, PetscBool useClosure, PetscInt *adjSize, PetscInt adj[])
230 {
231   PetscInt      *star = NULL;
232   PetscInt       numAdj = 0, maxAdjSize = *adjSize, starSize, s;
233   PetscErrorCode ierr;
234 
235   PetscFunctionBeginHot;
236   ierr = DMPlexGetTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr);
237   for (s = 0; s < starSize*2; s += 2) {
238     const PetscInt *closure = NULL;
239     PetscInt        closureSize, c, q;
240 
241     ierr = DMPlexGetTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
242     for (c = 0; c < closureSize*2; c += 2) {
243       for (q = 0; q < numAdj || (adj[numAdj++] = closure[c],0); ++q) {
244         if (closure[c] == adj[q]) break;
245       }
246       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
247     }
248     ierr = DMPlexRestoreTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
249   }
250   ierr = DMPlexRestoreTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr);
251   *adjSize = numAdj;
252   PetscFunctionReturn(0);
253 }
254 
255 #undef __FUNCT__
256 #define __FUNCT__ "DMPlexGetAdjacency_Internal"
257 PetscErrorCode DMPlexGetAdjacency_Internal(DM dm, PetscInt p, PetscBool useCone, PetscBool useTransitiveClosure, PetscBool useAnchors, PetscInt *adjSize, PetscInt *adj[])
258 {
259   static PetscInt asiz = 0;
260   PetscInt maxAnchors = 1;
261   PetscInt aStart = -1, aEnd = -1;
262   PetscInt maxAdjSize;
263   PetscSection aSec = NULL;
264   IS aIS = NULL;
265   const PetscInt *anchors;
266   PetscErrorCode  ierr;
267 
268   PetscFunctionBeginHot;
269   if (useAnchors) {
270     ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr);
271     if (aSec) {
272       ierr = PetscSectionGetMaxDof(aSec,&maxAnchors);CHKERRQ(ierr);
273       maxAnchors = PetscMax(1,maxAnchors);
274       ierr = PetscSectionGetChart(aSec,&aStart,&aEnd);CHKERRQ(ierr);
275       ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr);
276     }
277   }
278   if (!*adj) {
279     PetscInt depth, coneSeries, supportSeries, maxC, maxS, pStart, pEnd;
280 
281     ierr  = DMPlexGetChart(dm, &pStart,&pEnd);CHKERRQ(ierr);
282     ierr  = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
283     ierr  = DMPlexGetMaxSizes(dm, &maxC, &maxS);CHKERRQ(ierr);
284     coneSeries    = (maxC > 1) ? ((PetscPowInt(maxC,depth+1)-1)/(maxC-1)) : depth+1;
285     supportSeries = (maxS > 1) ? ((PetscPowInt(maxS,depth+1)-1)/(maxS-1)) : depth+1;
286     asiz  = PetscMax(PetscPowInt(maxS,depth)*coneSeries,PetscPowInt(maxC,depth)*supportSeries);
287     asiz *= maxAnchors;
288     asiz  = PetscMin(asiz,pEnd-pStart);
289     ierr  = PetscMalloc1(asiz,adj);CHKERRQ(ierr);
290   }
291   if (*adjSize < 0) *adjSize = asiz;
292   maxAdjSize = *adjSize;
293   if (useTransitiveClosure) {
294     ierr = DMPlexGetAdjacency_Transitive_Internal(dm, p, useCone, adjSize, *adj);CHKERRQ(ierr);
295   } else if (useCone) {
296     ierr = DMPlexGetAdjacency_Cone_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
297   } else {
298     ierr = DMPlexGetAdjacency_Support_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
299   }
300   if (useAnchors && aSec) {
301     PetscInt origSize = *adjSize;
302     PetscInt numAdj = origSize;
303     PetscInt i = 0, j;
304     PetscInt *orig = *adj;
305 
306     while (i < origSize) {
307       PetscInt p = orig[i];
308       PetscInt aDof = 0;
309 
310       if (p >= aStart && p < aEnd) {
311         ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr);
312       }
313       if (aDof) {
314         PetscInt aOff;
315         PetscInt s, q;
316 
317         for (j = i + 1; j < numAdj; j++) {
318           orig[j - 1] = orig[j];
319         }
320         origSize--;
321         numAdj--;
322         ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr);
323         for (s = 0; s < aDof; ++s) {
324           for (q = 0; q < numAdj || (orig[numAdj++] = anchors[aOff+s],0); ++q) {
325             if (anchors[aOff+s] == orig[q]) break;
326           }
327           if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
328         }
329       }
330       else {
331         i++;
332       }
333     }
334     *adjSize = numAdj;
335     ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr);
336   }
337   PetscFunctionReturn(0);
338 }
339 
340 #undef __FUNCT__
341 #define __FUNCT__ "DMPlexGetAdjacency"
342 /*@
343   DMPlexGetAdjacency - Return all points adjacent to the given point
344 
345   Input Parameters:
346 + dm - The DM object
347 . p  - The point
348 . adjSize - The maximum size of adj if it is non-NULL, or PETSC_DETERMINE
349 - adj - Either NULL so that the array is allocated, or an existing array with size adjSize
350 
351   Output Parameters:
352 + adjSize - The number of adjacent points
353 - adj - The adjacent points
354 
355   Level: advanced
356 
357   Notes: The user must PetscFree the adj array if it was not passed in.
358 
359 .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexDistribute(), DMCreateMatrix(), DMPlexPreallocateOperator()
360 @*/
361 PetscErrorCode DMPlexGetAdjacency(DM dm, PetscInt p, PetscInt *adjSize, PetscInt *adj[])
362 {
363   DM_Plex       *mesh = (DM_Plex *) dm->data;
364   PetscErrorCode ierr;
365 
366   PetscFunctionBeginHot;
367   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
368   PetscValidPointer(adjSize,3);
369   PetscValidPointer(adj,4);
370   ierr = DMPlexGetAdjacency_Internal(dm, p, mesh->useCone, mesh->useClosure, mesh->useAnchors, adjSize, adj);CHKERRQ(ierr);
371   PetscFunctionReturn(0);
372 }
373 #undef __FUNCT__
374 #define __FUNCT__ "DMPlexCreateTwoSidedProcessSF"
375 /*@
376   DMPlexCreateTwoSidedProcessSF - Create an SF which just has process connectivity
377 
378   Collective on DM
379 
380   Input Parameters:
381 + dm      - The DM
382 - sfPoint - The PetscSF which encodes point connectivity
383 
384   Output Parameters:
385 + processRanks - A list of process neighbors, or NULL
386 - sfProcess    - An SF encoding the two-sided process connectivity, or NULL
387 
388   Level: developer
389 
390 .seealso: PetscSFCreate(), DMPlexCreateProcessSF()
391 @*/
392 PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM dm, PetscSF sfPoint, PetscSection rootRankSection, IS rootRanks, PetscSection leafRankSection, IS leafRanks, IS *processRanks, PetscSF *sfProcess)
393 {
394   const PetscSFNode *remotePoints;
395   PetscInt          *localPointsNew;
396   PetscSFNode       *remotePointsNew;
397   const PetscInt    *nranks;
398   PetscInt          *ranksNew;
399   PetscBT            neighbors;
400   PetscInt           pStart, pEnd, p, numLeaves, l, numNeighbors, n;
401   PetscMPIInt        numProcs, proc, rank;
402   PetscErrorCode     ierr;
403 
404   PetscFunctionBegin;
405   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
406   PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2);
407   if (processRanks) {PetscValidPointer(processRanks, 3);}
408   if (sfProcess)    {PetscValidPointer(sfProcess, 4);}
409   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr);
410   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
411   ierr = PetscSFGetGraph(sfPoint, NULL, &numLeaves, NULL, &remotePoints);CHKERRQ(ierr);
412   ierr = PetscBTCreate(numProcs, &neighbors);CHKERRQ(ierr);
413   ierr = PetscBTMemzero(numProcs, neighbors);CHKERRQ(ierr);
414   /* Compute root-to-leaf process connectivity */
415   ierr = PetscSectionGetChart(rootRankSection, &pStart, &pEnd);CHKERRQ(ierr);
416   ierr = ISGetIndices(rootRanks, &nranks);CHKERRQ(ierr);
417   for (p = pStart; p < pEnd; ++p) {
418     PetscInt ndof, noff, n;
419 
420     ierr = PetscSectionGetDof(rootRankSection, p, &ndof);CHKERRQ(ierr);
421     ierr = PetscSectionGetOffset(rootRankSection, p, &noff);CHKERRQ(ierr);
422     for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);}
423   }
424   ierr = ISRestoreIndices(rootRanks, &nranks);CHKERRQ(ierr);
425   /* Compute leaf-to-neighbor process connectivity */
426   ierr = PetscSectionGetChart(leafRankSection, &pStart, &pEnd);CHKERRQ(ierr);
427   ierr = ISGetIndices(leafRanks, &nranks);CHKERRQ(ierr);
428   for (p = pStart; p < pEnd; ++p) {
429     PetscInt ndof, noff, n;
430 
431     ierr = PetscSectionGetDof(leafRankSection, p, &ndof);CHKERRQ(ierr);
432     ierr = PetscSectionGetOffset(leafRankSection, p, &noff);CHKERRQ(ierr);
433     for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);}
434   }
435   ierr = ISRestoreIndices(leafRanks, &nranks);CHKERRQ(ierr);
436   /* Compute leaf-to-root process connectivity */
437   for (l = 0; l < numLeaves; ++l) {PetscBTSet(neighbors, remotePoints[l].rank);}
438   /* Calculate edges */
439   PetscBTClear(neighbors, rank);
440   for(proc = 0, numNeighbors = 0; proc < numProcs; ++proc) {if (PetscBTLookup(neighbors, proc)) ++numNeighbors;}
441   ierr = PetscMalloc1(numNeighbors, &ranksNew);CHKERRQ(ierr);
442   ierr = PetscMalloc1(numNeighbors, &localPointsNew);CHKERRQ(ierr);
443   ierr = PetscMalloc1(numNeighbors, &remotePointsNew);CHKERRQ(ierr);
444   for(proc = 0, n = 0; proc < numProcs; ++proc) {
445     if (PetscBTLookup(neighbors, proc)) {
446       ranksNew[n]              = proc;
447       localPointsNew[n]        = proc;
448       remotePointsNew[n].index = rank;
449       remotePointsNew[n].rank  = proc;
450       ++n;
451     }
452   }
453   ierr = PetscBTDestroy(&neighbors);CHKERRQ(ierr);
454   if (processRanks) {ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), numNeighbors, ranksNew, PETSC_OWN_POINTER, processRanks);CHKERRQ(ierr);}
455   else              {ierr = PetscFree(ranksNew);CHKERRQ(ierr);}
456   if (sfProcess) {
457     ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess);CHKERRQ(ierr);
458     ierr = PetscObjectSetName((PetscObject) *sfProcess, "Two-Sided Process SF");CHKERRQ(ierr);
459     ierr = PetscSFSetFromOptions(*sfProcess);CHKERRQ(ierr);
460     ierr = PetscSFSetGraph(*sfProcess, numProcs, numNeighbors, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr);
461   }
462   PetscFunctionReturn(0);
463 }
464 
465 #undef __FUNCT__
466 #define __FUNCT__ "DMPlexDistributeOwnership"
467 /*@
468   DMPlexDistributeOwnership - Compute owner information for shared points. This basically gets two-sided for an SF.
469 
470   Collective on DM
471 
472   Input Parameter:
473 . dm - The DM
474 
475   Output Parameters:
476 + rootSection - The number of leaves for a given root point
477 . rootrank    - The rank of each edge into the root point
478 . leafSection - The number of processes sharing a given leaf point
479 - leafrank    - The rank of each process sharing a leaf point
480 
481   Level: developer
482 
483 .seealso: DMPlexCreateOverlap()
484 @*/
485 PetscErrorCode DMPlexDistributeOwnership(DM dm, PetscSection rootSection, IS *rootrank, PetscSection leafSection, IS *leafrank)
486 {
487   MPI_Comm        comm;
488   PetscSF         sfPoint;
489   const PetscInt *rootdegree;
490   PetscInt       *myrank, *remoterank;
491   PetscInt        pStart, pEnd, p, nedges;
492   PetscMPIInt     rank;
493   PetscErrorCode  ierr;
494 
495   PetscFunctionBegin;
496   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
497   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
498   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
499   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
500   /* Compute number of leaves for each root */
501   ierr = PetscObjectSetName((PetscObject) rootSection, "Root Section");CHKERRQ(ierr);
502   ierr = PetscSectionSetChart(rootSection, pStart, pEnd);CHKERRQ(ierr);
503   ierr = PetscSFComputeDegreeBegin(sfPoint, &rootdegree);CHKERRQ(ierr);
504   ierr = PetscSFComputeDegreeEnd(sfPoint, &rootdegree);CHKERRQ(ierr);
505   for (p = pStart; p < pEnd; ++p) {ierr = PetscSectionSetDof(rootSection, p, rootdegree[p-pStart]);CHKERRQ(ierr);}
506   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
507   /* Gather rank of each leaf to root */
508   ierr = PetscSectionGetStorageSize(rootSection, &nedges);CHKERRQ(ierr);
509   ierr = PetscMalloc1(pEnd-pStart, &myrank);CHKERRQ(ierr);
510   ierr = PetscMalloc1(nedges,  &remoterank);CHKERRQ(ierr);
511   for (p = 0; p < pEnd-pStart; ++p) myrank[p] = rank;
512   ierr = PetscSFGatherBegin(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr);
513   ierr = PetscSFGatherEnd(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr);
514   ierr = PetscFree(myrank);CHKERRQ(ierr);
515   ierr = ISCreateGeneral(comm, nedges, remoterank, PETSC_OWN_POINTER, rootrank);CHKERRQ(ierr);
516   /* Distribute remote ranks to leaves */
517   ierr = PetscObjectSetName((PetscObject) leafSection, "Leaf Section");CHKERRQ(ierr);
518   ierr = DMPlexDistributeFieldIS(dm, sfPoint, rootSection, *rootrank, leafSection, leafrank);CHKERRQ(ierr);
519   PetscFunctionReturn(0);
520 }
521 
522 #undef __FUNCT__
523 #define __FUNCT__ "DMPlexCreateOverlap"
524 /*@C
525   DMPlexCreateOverlap - Compute owner information for shared points. This basically gets two-sided for an SF.
526 
527   Collective on DM
528 
529   Input Parameters:
530 + dm          - The DM
531 . rootSection - The number of leaves for a given root point
532 . rootrank    - The rank of each edge into the root point
533 . leafSection - The number of processes sharing a given leaf point
534 - leafrank    - The rank of each process sharing a leaf point
535 
536   Output Parameters:
537 + ovRootSection - The number of new overlap points for each neighboring process
538 . ovRootPoints  - The new overlap points for each neighboring process
539 . ovLeafSection - The number of new overlap points from each neighboring process
540 - ovLeafPoints  - The new overlap points from each neighboring process
541 
542   Level: developer
543 
544 .seealso: DMPlexDistributeOwnership()
545 @*/
546 PetscErrorCode DMPlexCreateOverlap(DM dm, PetscSection rootSection, IS rootrank, PetscSection leafSection, IS leafrank, PetscSection ovRootSection, PetscSFNode **ovRootPoints, PetscSection ovLeafSection, PetscSFNode **ovLeafPoints)
547 {
548   DMLabel            ovAdjByRank; /* A DMLabel containing all points adjacent to shared points, separated by rank (value in label) */
549   PetscSF            sfPoint, sfProc;
550   IS                 valueIS;
551   const PetscSFNode *remote;
552   const PetscInt    *local;
553   const PetscInt    *nrank, *rrank, *neighbors;
554   PetscInt          *adj = NULL;
555   PetscInt           pStart, pEnd, p, sStart, sEnd, nleaves, l, numNeighbors, n, ovSize;
556   PetscMPIInt        rank, numProcs;
557   PetscErrorCode     ierr;
558 
559   PetscFunctionBegin;
560   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr);
561   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
562   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
563   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
564   ierr = PetscSectionGetChart(leafSection, &sStart, &sEnd);CHKERRQ(ierr);
565   ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr);
566   ierr = DMLabelCreate("Overlap adjacency", &ovAdjByRank);CHKERRQ(ierr);
567   /* Handle leaves: shared with the root point */
568   for (l = 0; l < nleaves; ++l) {
569     PetscInt adjSize = PETSC_DETERMINE, a;
570 
571     ierr = DMPlexGetAdjacency(dm, local[l], &adjSize, &adj);CHKERRQ(ierr);
572     for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remote[l].rank);CHKERRQ(ierr);}
573   }
574   ierr = ISGetIndices(rootrank, &rrank);CHKERRQ(ierr);
575   ierr = ISGetIndices(leafrank, &nrank);CHKERRQ(ierr);
576   /* Handle roots */
577   for (p = pStart; p < pEnd; ++p) {
578     PetscInt adjSize = PETSC_DETERMINE, neighbors = 0, noff, n, a;
579 
580     if ((p >= sStart) && (p < sEnd)) {
581       /* Some leaves share a root with other leaves on different processes */
582       ierr = PetscSectionGetDof(leafSection, p, &neighbors);CHKERRQ(ierr);
583       if (neighbors) {
584         ierr = PetscSectionGetOffset(leafSection, p, &noff);CHKERRQ(ierr);
585         ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
586         for (n = 0; n < neighbors; ++n) {
587           const PetscInt remoteRank = nrank[noff+n];
588 
589           if (remoteRank == rank) continue;
590           for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);}
591         }
592       }
593     }
594     /* Roots are shared with leaves */
595     ierr = PetscSectionGetDof(rootSection, p, &neighbors);CHKERRQ(ierr);
596     if (!neighbors) continue;
597     ierr = PetscSectionGetOffset(rootSection, p, &noff);CHKERRQ(ierr);
598     ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
599     for (n = 0; n < neighbors; ++n) {
600       const PetscInt remoteRank = rrank[noff+n];
601 
602       if (remoteRank == rank) continue;
603       for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);}
604     }
605   }
606   ierr = PetscFree(adj);CHKERRQ(ierr);
607   ierr = ISRestoreIndices(rootrank, &rrank);CHKERRQ(ierr);
608   ierr = ISRestoreIndices(leafrank, &nrank);CHKERRQ(ierr);
609   {
610     ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr);
611     ierr = DMLabelView(ovAdjByRank, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
612   }
613   /* Convert to (point, rank) and use actual owners */
614   ierr = PetscSectionSetChart(ovRootSection, 0, numProcs);CHKERRQ(ierr);
615   ierr = DMLabelGetValueIS(ovAdjByRank, &valueIS);CHKERRQ(ierr);
616   ierr = ISGetLocalSize(valueIS, &numNeighbors);CHKERRQ(ierr);
617   ierr = ISGetIndices(valueIS, &neighbors);CHKERRQ(ierr);
618   for (n = 0; n < numNeighbors; ++n) {
619     PetscInt numPoints;
620 
621     ierr = DMLabelGetStratumSize(ovAdjByRank, neighbors[n], &numPoints);CHKERRQ(ierr);
622     ierr = PetscSectionAddDof(ovRootSection, neighbors[n], numPoints);CHKERRQ(ierr);
623   }
624   ierr = PetscSectionSetUp(ovRootSection);CHKERRQ(ierr);
625   ierr = PetscSectionGetStorageSize(ovRootSection, &ovSize);CHKERRQ(ierr);
626   ierr = PetscMalloc1(ovSize, ovRootPoints);CHKERRQ(ierr);
627   ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr);
628   for (n = 0; n < numNeighbors; ++n) {
629     IS              pointIS;
630     const PetscInt *points;
631     PetscInt        off, numPoints, p;
632 
633     ierr = PetscSectionGetOffset(ovRootSection, neighbors[n], &off);CHKERRQ(ierr);
634     ierr = DMLabelGetStratumIS(ovAdjByRank, neighbors[n], &pointIS);CHKERRQ(ierr);
635     ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
636     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
637     for (p = 0; p < numPoints; ++p) {
638       ierr = PetscFindInt(points[p], nleaves, local, &l);CHKERRQ(ierr);
639       if (l >= 0) {(*ovRootPoints)[off+p] = remote[l];}
640       else        {(*ovRootPoints)[off+p].index = points[p]; (*ovRootPoints)[off+p].rank = rank;}
641     }
642     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
643     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
644   }
645   ierr = ISRestoreIndices(valueIS, &neighbors);CHKERRQ(ierr);
646   ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
647   ierr = DMLabelDestroy(&ovAdjByRank);CHKERRQ(ierr);
648   /* Make process SF */
649   ierr = DMPlexCreateTwoSidedProcessSF(dm, sfPoint, rootSection, rootrank, leafSection, leafrank, NULL, &sfProc);CHKERRQ(ierr);
650   {
651     ierr = PetscSFView(sfProc, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
652   }
653   /* Communicate overlap */
654   ierr = DMPlexDistributeData(dm, sfProc, ovRootSection, MPIU_2INT, (void *) *ovRootPoints, ovLeafSection, (void **) ovLeafPoints);CHKERRQ(ierr);
655   ierr = PetscSFDestroy(&sfProc);CHKERRQ(ierr);
656   PetscFunctionReturn(0);
657 }
658 
659 #undef __FUNCT__
660 #define __FUNCT__ "DMPlexDistributeField"
661 /*@
662   DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
663 
664   Collective on DM
665 
666   Input Parameters:
667 + dm - The DMPlex object
668 . pointSF - The PetscSF describing the communication pattern
669 . originalSection - The PetscSection for existing data layout
670 - originalVec - The existing data
671 
672   Output Parameters:
673 + newSection - The PetscSF describing the new data layout
674 - newVec - The new data
675 
676   Level: developer
677 
678 .seealso: DMPlexDistribute(), DMPlexDistributeFieldIS(), DMPlexDistributeData()
679 @*/
680 PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec)
681 {
682   PetscSF        fieldSF;
683   PetscInt      *remoteOffsets, fieldSize;
684   PetscScalar   *originalValues, *newValues;
685   PetscErrorCode ierr;
686 
687   PetscFunctionBegin;
688   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
689   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
690 
691   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
692   ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr);
693   ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr);
694 
695   ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr);
696   ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr);
697   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
698   ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
699   ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
700   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
701   ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr);
702   ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr);
703   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
704   PetscFunctionReturn(0);
705 }
706 
707 #undef __FUNCT__
708 #define __FUNCT__ "DMPlexDistributeFieldIS"
709 /*@
710   DMPlexDistributeFieldIS - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
711 
712   Collective on DM
713 
714   Input Parameters:
715 + dm - The DMPlex object
716 . pointSF - The PetscSF describing the communication pattern
717 . originalSection - The PetscSection for existing data layout
718 - originalIS - The existing data
719 
720   Output Parameters:
721 + newSection - The PetscSF describing the new data layout
722 - newIS - The new data
723 
724   Level: developer
725 
726 .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexDistributeData()
727 @*/
728 PetscErrorCode DMPlexDistributeFieldIS(DM dm, PetscSF pointSF, PetscSection originalSection, IS originalIS, PetscSection newSection, IS *newIS)
729 {
730   PetscSF         fieldSF;
731   PetscInt       *newValues, *remoteOffsets, fieldSize;
732   const PetscInt *originalValues;
733   PetscErrorCode  ierr;
734 
735   PetscFunctionBegin;
736   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
737   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
738 
739   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
740   ierr = PetscMalloc(fieldSize * sizeof(PetscInt), &newValues);CHKERRQ(ierr);
741 
742   ierr = ISGetIndices(originalIS, &originalValues);CHKERRQ(ierr);
743   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
744   ierr = PetscSFBcastBegin(fieldSF, MPIU_INT, originalValues, newValues);CHKERRQ(ierr);
745   ierr = PetscSFBcastEnd(fieldSF, MPIU_INT, originalValues, newValues);CHKERRQ(ierr);
746   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
747   ierr = ISRestoreIndices(originalIS, &originalValues);CHKERRQ(ierr);
748   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) pointSF), fieldSize, newValues, PETSC_OWN_POINTER, newIS);CHKERRQ(ierr);
749   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
750   PetscFunctionReturn(0);
751 }
752 
753 #undef __FUNCT__
754 #define __FUNCT__ "DMPlexDistributeData"
755 /*@
756   DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
757 
758   Collective on DM
759 
760   Input Parameters:
761 + dm - The DMPlex object
762 . pointSF - The PetscSF describing the communication pattern
763 . originalSection - The PetscSection for existing data layout
764 . datatype - The type of data
765 - originalData - The existing data
766 
767   Output Parameters:
768 + newSection - The PetscSection describing the new data layout
769 - newData - The new data
770 
771   Level: developer
772 
773 .seealso: DMPlexDistribute(), DMPlexDistributeField()
774 @*/
775 PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData)
776 {
777   PetscSF        fieldSF;
778   PetscInt      *remoteOffsets, fieldSize;
779   PetscMPIInt    dataSize;
780   PetscErrorCode ierr;
781 
782   PetscFunctionBegin;
783   ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
784   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
785 
786   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
787   ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr);
788   ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr);
789 
790   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
791   ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
792   ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
793   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
794   ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
795   PetscFunctionReturn(0);
796 }
797 
798 #undef __FUNCT__
799 #define __FUNCT__ "DMPlexDistributeCones"
800 PetscErrorCode DMPlexDistributeCones(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
801 {
802   DM_Plex               *pmesh   = (DM_Plex*) (dmParallel)->data;
803   MPI_Comm               comm;
804   PetscSF                coneSF;
805   PetscSection           originalConeSection, newConeSection;
806   PetscInt              *remoteOffsets, *cones, *newCones, newConesSize;
807   PetscBool              flg;
808   PetscErrorCode         ierr;
809 
810   PetscFunctionBegin;
811   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
812   PetscValidPointer(dmParallel,4);
813   ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
814 
815   /* Distribute cone section */
816   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
817   ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr);
818   ierr = DMPlexGetConeSection(dmParallel, &newConeSection);CHKERRQ(ierr);
819   ierr = PetscSFDistributeSection(migrationSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr);
820   ierr = DMSetUp(dmParallel);CHKERRQ(ierr);
821   {
822     PetscInt pStart, pEnd, p;
823 
824     ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr);
825     for (p = pStart; p < pEnd; ++p) {
826       PetscInt coneSize;
827       ierr               = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr);
828       pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize);
829     }
830   }
831   /* Communicate and renumber cones */
832   ierr = PetscSFCreateSectionSF(migrationSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr);
833   ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
834   ierr = DMPlexGetCones(dmParallel, &newCones);CHKERRQ(ierr);
835   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
836   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
837   ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr);
838   ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr);
839   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr);
840   if (flg) {
841     ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr);
842     ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
843     ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr);
844     ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
845     ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr);
846   }
847   ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr);
848   ierr = DMPlexGetConeOrientations(dmParallel, &newCones);CHKERRQ(ierr);
849   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
850   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
851   ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr);
852   ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
853   /* Create supports and stratify sieve */
854   {
855     PetscInt pStart, pEnd;
856 
857     ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr);
858     ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr);
859   }
860   ierr = DMPlexSymmetrize(dmParallel);CHKERRQ(ierr);
861   ierr = DMPlexStratify(dmParallel);CHKERRQ(ierr);
862   PetscFunctionReturn(0);
863 }
864 
865 #undef __FUNCT__
866 #define __FUNCT__ "DMPlexDistribute"
867 /*@C
868   DMPlexDistribute - Distributes the mesh and any associated sections.
869 
870   Not Collective
871 
872   Input Parameter:
873 + dm  - The original DMPlex object
874 . partitioner - The partitioning package, or NULL for the default
875 - overlap - The overlap of partitions, 0 is the default
876 
877   Output Parameter:
878 + sf - The PetscSF used for point distribution
879 - parallelMesh - The distributed DMPlex object, or NULL
880 
881   Note: If the mesh was not distributed, the return value is NULL.
882 
883   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
884   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
885   representation on the mesh.
886 
887   Level: intermediate
888 
889 .keywords: mesh, elements
890 .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
891 @*/
892 PetscErrorCode DMPlexDistribute(DM dm, const char partitioner[], PetscInt overlap, PetscSF *sf, DM *dmParallel)
893 {
894   DM_Plex               *mesh   = (DM_Plex*) dm->data, *pmesh;
895   MPI_Comm               comm;
896   PetscInt               dim, numRemoteRanks;
897   IS                     origCellPart,        origPart,        cellPart,        part;
898   PetscSection           origCellPartSection, origPartSection, cellPartSection, partSection;
899   PetscSFNode           *remoteRanks;
900   PetscSF                partSF, pointSF;
901   ISLocalToGlobalMapping renumbering;
902   PetscBool              flg;
903   PetscMPIInt            rank, numProcs, p;
904   PetscErrorCode         ierr;
905 
906   PetscFunctionBegin;
907   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
908   if (sf) PetscValidPointer(sf,4);
909   PetscValidPointer(dmParallel,5);
910 
911   ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
912   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
913   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
914   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
915 
916   *dmParallel = NULL;
917   if (numProcs == 1) PetscFunctionReturn(0);
918 
919   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
920   /* Create cell partition - We need to rewrite to use IS, use the MatPartition stuff */
921   ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
922   if (overlap > 1) SETERRQ(comm, PETSC_ERR_SUP, "Overlap > 1 not yet implemented");
923   ierr = PetscSectionCreate(comm, &cellPartSection);CHKERRQ(ierr);
924   ierr = PetscSectionCreate(comm, &origCellPartSection);CHKERRQ(ierr);
925   ierr = PetscPartitionerPartition(mesh->partitioner, dm, overlap > 0 ? PETSC_TRUE : PETSC_FALSE, cellPartSection, &cellPart, origCellPartSection, &origCellPart);CHKERRQ(ierr);
926   /* Create SF assuming a serial partition for all processes: Could check for IS length here */
927   if (!rank) numRemoteRanks = numProcs;
928   else       numRemoteRanks = 0;
929   ierr = PetscMalloc1(numRemoteRanks, &remoteRanks);CHKERRQ(ierr);
930   for (p = 0; p < numRemoteRanks; ++p) {
931     remoteRanks[p].rank  = p;
932     remoteRanks[p].index = 0;
933   }
934   ierr = PetscSFCreate(comm, &partSF);CHKERRQ(ierr);
935   ierr = PetscSFSetGraph(partSF, 1, numRemoteRanks, NULL, PETSC_OWN_POINTER, remoteRanks, PETSC_OWN_POINTER);CHKERRQ(ierr);
936   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
937   if (flg) {
938     ierr = PetscPrintf(comm, "Cell Partition:\n");CHKERRQ(ierr);
939     ierr = PetscSectionView(cellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
940     ierr = ISView(cellPart, NULL);CHKERRQ(ierr);
941     if (origCellPart) {
942       ierr = PetscPrintf(comm, "Original Cell Partition:\n");CHKERRQ(ierr);
943       ierr = PetscSectionView(origCellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
944       ierr = ISView(origCellPart, NULL);CHKERRQ(ierr);
945     }
946     ierr = PetscSFView(partSF, NULL);CHKERRQ(ierr);
947   }
948   /* Close the partition over the mesh */
949   ierr = DMPlexCreatePartitionClosure(dm, cellPartSection, cellPart, &partSection, &part);CHKERRQ(ierr);
950   ierr = ISDestroy(&cellPart);CHKERRQ(ierr);
951   ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr);
952   /* Create new mesh */
953   ierr  = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr);
954   ierr  = DMSetDimension(*dmParallel, dim);CHKERRQ(ierr);
955   ierr  = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr);
956   pmesh = (DM_Plex*) (*dmParallel)->data;
957   /* Distribute sieve points and the global point numbering (replaces creating remote bases) */
958   ierr = PetscSFConvertPartition(partSF, partSection, part, &renumbering, &pointSF);CHKERRQ(ierr);
959   if (flg) {
960     ierr = PetscPrintf(comm, "Point Partition:\n");CHKERRQ(ierr);
961     ierr = PetscSectionView(partSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
962     ierr = ISView(part, NULL);CHKERRQ(ierr);
963     ierr = PetscSFView(pointSF, NULL);CHKERRQ(ierr);
964     ierr = PetscPrintf(comm, "Point Renumbering after partition:\n");CHKERRQ(ierr);
965     ierr = ISLocalToGlobalMappingView(renumbering, NULL);CHKERRQ(ierr);
966   }
967   ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
968 
969   ierr = DMPlexDistributeCones(dm, pointSF, renumbering, *dmParallel);CHKERRQ(ierr);
970 
971   /* Create point SF for parallel mesh */
972   ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
973   {
974     const PetscInt *leaves;
975     PetscSFNode    *remotePoints, *rowners, *lowners;
976     PetscInt        numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints;
977     PetscInt        pStart, pEnd;
978 
979     ierr = DMPlexGetChart(*dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
980     ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr);
981     ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr);
982     for (p=0; p<numRoots; p++) {
983       rowners[p].rank  = -1;
984       rowners[p].index = -1;
985     }
986     if (origCellPart) {
987       /* Make sure points in the original partition are not assigned to other procs */
988       const PetscInt *origPoints;
989 
990       ierr = DMPlexCreatePartitionClosure(dm, origCellPartSection, origCellPart, &origPartSection, &origPart);CHKERRQ(ierr);
991       ierr = ISGetIndices(origPart, &origPoints);CHKERRQ(ierr);
992       for (p = 0; p < numProcs; ++p) {
993         PetscInt dof, off, d;
994 
995         ierr = PetscSectionGetDof(origPartSection, p, &dof);CHKERRQ(ierr);
996         ierr = PetscSectionGetOffset(origPartSection, p, &off);CHKERRQ(ierr);
997         for (d = off; d < off+dof; ++d) {
998           rowners[origPoints[d]].rank = p;
999         }
1000       }
1001       ierr = ISRestoreIndices(origPart, &origPoints);CHKERRQ(ierr);
1002       ierr = ISDestroy(&origPart);CHKERRQ(ierr);
1003       ierr = PetscSectionDestroy(&origPartSection);CHKERRQ(ierr);
1004     }
1005     ierr = ISDestroy(&origCellPart);CHKERRQ(ierr);
1006     ierr = PetscSectionDestroy(&origCellPartSection);CHKERRQ(ierr);
1007 
1008     ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
1009     ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
1010     for (p = 0; p < numLeaves; ++p) {
1011       if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */
1012         lowners[p].rank  = rank;
1013         lowners[p].index = leaves ? leaves[p] : p;
1014       } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */
1015         lowners[p].rank  = -2;
1016         lowners[p].index = -2;
1017       }
1018     }
1019     for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */
1020       rowners[p].rank  = -3;
1021       rowners[p].index = -3;
1022     }
1023     ierr = PetscSFReduceBegin(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
1024     ierr = PetscSFReduceEnd(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
1025     ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
1026     ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
1027     for (p = 0; p < numLeaves; ++p) {
1028       if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed");
1029       if (lowners[p].rank != rank) ++numGhostPoints;
1030     }
1031     ierr = PetscMalloc1(numGhostPoints,    &ghostPoints);CHKERRQ(ierr);
1032     ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr);
1033     for (p = 0, gp = 0; p < numLeaves; ++p) {
1034       if (lowners[p].rank != rank) {
1035         ghostPoints[gp]        = leaves ? leaves[p] : p;
1036         remotePoints[gp].rank  = lowners[p].rank;
1037         remotePoints[gp].index = lowners[p].index;
1038         ++gp;
1039       }
1040     }
1041     ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr);
1042     ierr = PetscSFSetGraph((*dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
1043     ierr = PetscSFSetFromOptions((*dmParallel)->sf);CHKERRQ(ierr);
1044   }
1045   pmesh->useCone    = mesh->useCone;
1046   pmesh->useClosure = mesh->useClosure;
1047   pmesh->useAnchors = mesh->useAnchors;
1048   ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
1049   /* Distribute Coordinates */
1050   {
1051     PetscSection     originalCoordSection, newCoordSection;
1052     Vec              originalCoordinates, newCoordinates;
1053     PetscInt         bs;
1054     const char      *name;
1055     const PetscReal *maxCell, *L;
1056 
1057     ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr);
1058     ierr = DMGetCoordinateSection(*dmParallel, &newCoordSection);CHKERRQ(ierr);
1059     ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr);
1060     if (originalCoordinates) {
1061       ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr);
1062       ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr);
1063       ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr);
1064 
1065       ierr = DMPlexDistributeField(dm, pointSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr);
1066       ierr = DMSetCoordinatesLocal(*dmParallel, newCoordinates);CHKERRQ(ierr);
1067       ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr);
1068       ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr);
1069       ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr);
1070     }
1071     ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr);
1072     if (L) {ierr = DMSetPeriodicity(*dmParallel, maxCell, L);CHKERRQ(ierr);}
1073   }
1074   /* Distribute labels */
1075   ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
1076   {
1077     DMLabel  next      = mesh->labels, newNext = pmesh->labels;
1078     PetscInt numLabels = 0, l;
1079 
1080     /* Bcast number of labels */
1081     while (next) {++numLabels; next = next->next;}
1082     ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1083     next = mesh->labels;
1084     for (l = 0; l < numLabels; ++l) {
1085       DMLabel   labelNew;
1086       PetscBool isdepth;
1087 
1088       /* Skip "depth" because it is recreated */
1089       if (!rank) {ierr = PetscStrcmp(next->name, "depth", &isdepth);CHKERRQ(ierr);}
1090       ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
1091       if (isdepth) {if (!rank) next = next->next; continue;}
1092       ierr = DMLabelDistribute(next, partSection, part, renumbering, &labelNew);CHKERRQ(ierr);
1093       /* Insert into list */
1094       if (newNext) newNext->next = labelNew;
1095       else         pmesh->labels = labelNew;
1096       newNext = labelNew;
1097       if (!rank) next = next->next;
1098     }
1099   }
1100   ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
1101   /* Setup hybrid structure */
1102   {
1103     const PetscInt *gpoints;
1104     PetscInt        depth, n, d;
1105 
1106     for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];}
1107     ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1108     ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr);
1109     ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr);
1110     ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
1111     for (d = 0; d <= dim; ++d) {
1112       PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p;
1113 
1114       if (pmax < 0) continue;
1115       ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr);
1116       ierr = DMPlexGetDepthStratum(*dmParallel, d, NULL, &pEnd);CHKERRQ(ierr);
1117       ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr);
1118       for (p = 0; p < n; ++p) {
1119         const PetscInt point = gpoints[p];
1120 
1121         if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax;
1122       }
1123       if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax;
1124       else            pmesh->hybridPointMax[d] = -1;
1125     }
1126     ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr);
1127   }
1128   /* Set up tree */
1129   {
1130     DM              refTree;
1131     PetscSection    origParentSection, newParentSection;
1132     PetscInt        *origParents, *origChildIDs;
1133 
1134     ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr);
1135     ierr = DMPlexSetReferenceTree(*dmParallel,refTree);CHKERRQ(ierr);
1136     ierr = DMPlexGetTree(dm,&origParentSection,&origParents,&origChildIDs,NULL,NULL);CHKERRQ(ierr);
1137     if (origParentSection) {
1138       PetscInt        pStart, pEnd;
1139       PetscInt        *newParents, *newChildIDs;
1140       PetscInt        *remoteOffsetsParents, newParentSize;
1141       PetscSF         parentSF;
1142 
1143       ierr = DMPlexGetChart(*dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
1144       ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dmParallel),&newParentSection);CHKERRQ(ierr);
1145       ierr = PetscSectionSetChart(newParentSection,pStart,pEnd);CHKERRQ(ierr);
1146       ierr = PetscSFDistributeSection(pointSF, origParentSection, &remoteOffsetsParents, newParentSection);CHKERRQ(ierr);
1147       ierr = PetscSFCreateSectionSF(pointSF, origParentSection, remoteOffsetsParents, newParentSection, &parentSF);CHKERRQ(ierr);
1148       ierr = PetscSectionGetStorageSize(newParentSection,&newParentSize);CHKERRQ(ierr);
1149       ierr = PetscMalloc2(newParentSize,&newParents,newParentSize,&newChildIDs);CHKERRQ(ierr);
1150       ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
1151       ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
1152       ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1153       ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1154       ierr = ISGlobalToLocalMappingApplyBlock(renumbering,IS_GTOLM_MASK, newParentSize, newParents, NULL, newParents);CHKERRQ(ierr);
1155       ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-parents_view", &flg);CHKERRQ(ierr);
1156       if (flg) {
1157         ierr = PetscPrintf(comm, "Serial Parent Section: \n");CHKERRQ(ierr);
1158         ierr = PetscSectionView(origParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1159         ierr = PetscPrintf(comm, "Parallel Parent Section: \n");CHKERRQ(ierr);
1160         ierr = PetscSectionView(newParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1161         ierr = PetscSFView(parentSF, NULL);CHKERRQ(ierr);
1162       }
1163       ierr = DMPlexSetTree(*dmParallel,newParentSection,newParents,newChildIDs);CHKERRQ(ierr);
1164       ierr = PetscSectionDestroy(&newParentSection);CHKERRQ(ierr);
1165       ierr = PetscFree2(newParents,newChildIDs);CHKERRQ(ierr);
1166       ierr = PetscSFDestroy(&parentSF);CHKERRQ(ierr);
1167     }
1168   }
1169   /* Cleanup Partition */
1170   ierr = ISLocalToGlobalMappingDestroy(&renumbering);CHKERRQ(ierr);
1171   ierr = PetscSFDestroy(&partSF);CHKERRQ(ierr);
1172   ierr = PetscSectionDestroy(&partSection);CHKERRQ(ierr);
1173   ierr = ISDestroy(&part);CHKERRQ(ierr);
1174   /* Copy BC */
1175   ierr = DMPlexCopyBoundary(dm, *dmParallel);CHKERRQ(ierr);
1176   /* Cleanup */
1177   if (sf) {*sf = pointSF;}
1178   else    {ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr);}
1179   ierr = DMSetFromOptions(*dmParallel);CHKERRQ(ierr);
1180   ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
1181   PetscFunctionReturn(0);
1182 }
1183