xref: /petsc/src/dm/impls/plex/plexdistribute.c (revision f5476a2b443e0eedb22d31cb1b2d85b9b918f73f)
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 
374 #undef __FUNCT__
375 #define __FUNCT__ "DMPlexDistributeField"
376 /*@
377   DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
378 
379   Collective on DM
380 
381   Input Parameters:
382 + dm - The DMPlex object
383 . pointSF - The PetscSF describing the communication pattern
384 . originalSection - The PetscSection for existing data layout
385 - originalVec - The existing data
386 
387   Output Parameters:
388 + newSection - The PetscSF describing the new data layout
389 - newVec - The new data
390 
391   Level: developer
392 
393 .seealso: DMPlexDistribute(), DMPlexDistributeData()
394 @*/
395 PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec)
396 {
397   PetscSF        fieldSF;
398   PetscInt      *remoteOffsets, fieldSize;
399   PetscScalar   *originalValues, *newValues;
400   PetscErrorCode ierr;
401 
402   PetscFunctionBegin;
403   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
404   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
405 
406   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
407   ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr);
408   ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr);
409 
410   ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr);
411   ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr);
412   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
413   ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
414   ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
415   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
416   ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr);
417   ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr);
418   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
419   PetscFunctionReturn(0);
420 }
421 
422 #undef __FUNCT__
423 #define __FUNCT__ "DMPlexDistributeData"
424 /*@
425   DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
426 
427   Collective on DM
428 
429   Input Parameters:
430 + dm - The DMPlex object
431 . pointSF - The PetscSF describing the communication pattern
432 . originalSection - The PetscSection for existing data layout
433 . datatype - The type of data
434 - originalData - The existing data
435 
436   Output Parameters:
437 + newSection - The PetscSF describing the new data layout
438 - newData - The new data
439 
440   Level: developer
441 
442 .seealso: DMPlexDistribute(), DMPlexDistributeField()
443 @*/
444 PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData)
445 {
446   PetscSF        fieldSF;
447   PetscInt      *remoteOffsets, fieldSize;
448   PetscMPIInt    dataSize;
449   PetscErrorCode ierr;
450 
451   PetscFunctionBegin;
452   ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
453   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
454 
455   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
456   ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr);
457   ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr);
458 
459   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
460   ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
461   ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
462   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
463   ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
464   PetscFunctionReturn(0);
465 }
466 
467 #undef __FUNCT__
468 #define __FUNCT__ "DMPlexDistribute"
469 /*@C
470   DMPlexDistribute - Distributes the mesh and any associated sections.
471 
472   Not Collective
473 
474   Input Parameter:
475 + dm  - The original DMPlex object
476 . partitioner - The partitioning package, or NULL for the default
477 - overlap - The overlap of partitions, 0 is the default
478 
479   Output Parameter:
480 + sf - The PetscSF used for point distribution
481 - parallelMesh - The distributed DMPlex object, or NULL
482 
483   Note: If the mesh was not distributed, the return value is NULL.
484 
485   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
486   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
487   representation on the mesh.
488 
489   Level: intermediate
490 
491 .keywords: mesh, elements
492 .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
493 @*/
494 PetscErrorCode DMPlexDistribute(DM dm, const char partitioner[], PetscInt overlap, PetscSF *sf, DM *dmParallel)
495 {
496   DM_Plex               *mesh   = (DM_Plex*) dm->data, *pmesh;
497   MPI_Comm               comm;
498   const PetscInt         height = 0;
499   PetscInt               dim, numRemoteRanks;
500   IS                     origCellPart,        origPart,        cellPart,        part;
501   PetscSection           origCellPartSection, origPartSection, cellPartSection, partSection;
502   PetscSFNode           *remoteRanks;
503   PetscSF                partSF, pointSF, coneSF;
504   ISLocalToGlobalMapping renumbering;
505   PetscSection           originalConeSection, newConeSection;
506   PetscInt              *remoteOffsets;
507   PetscInt              *cones, *newCones, newConesSize;
508   PetscBool              flg;
509   PetscMPIInt            rank, numProcs, p;
510   PetscErrorCode         ierr;
511 
512   PetscFunctionBegin;
513   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
514   if (sf) PetscValidPointer(sf,4);
515   PetscValidPointer(dmParallel,5);
516 
517   ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
518   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
519   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
520   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
521 
522   *dmParallel = NULL;
523   if (numProcs == 1) PetscFunctionReturn(0);
524 
525   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
526   /* Create cell partition - We need to rewrite to use IS, use the MatPartition stuff */
527   ierr = PetscLogEventBegin(DMPLEX_Partition,dm,0,0,0);CHKERRQ(ierr);
528   if (overlap > 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Overlap > 1 not yet implemented");
529   ierr = DMPlexCreatePartition(dm, partitioner, height, overlap > 0 ? PETSC_TRUE : PETSC_FALSE, &cellPartSection, &cellPart, &origCellPartSection, &origCellPart);CHKERRQ(ierr);
530   /* Create SF assuming a serial partition for all processes: Could check for IS length here */
531   if (!rank) numRemoteRanks = numProcs;
532   else       numRemoteRanks = 0;
533   ierr = PetscMalloc1(numRemoteRanks, &remoteRanks);CHKERRQ(ierr);
534   for (p = 0; p < numRemoteRanks; ++p) {
535     remoteRanks[p].rank  = p;
536     remoteRanks[p].index = 0;
537   }
538   ierr = PetscSFCreate(comm, &partSF);CHKERRQ(ierr);
539   ierr = PetscSFSetGraph(partSF, 1, numRemoteRanks, NULL, PETSC_OWN_POINTER, remoteRanks, PETSC_OWN_POINTER);CHKERRQ(ierr);
540   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
541   if (flg) {
542     ierr = PetscPrintf(comm, "Cell Partition:\n");CHKERRQ(ierr);
543     ierr = PetscSectionView(cellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
544     ierr = ISView(cellPart, NULL);CHKERRQ(ierr);
545     if (origCellPart) {
546       ierr = PetscPrintf(comm, "Original Cell Partition:\n");CHKERRQ(ierr);
547       ierr = PetscSectionView(origCellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
548       ierr = ISView(origCellPart, NULL);CHKERRQ(ierr);
549     }
550     ierr = PetscSFView(partSF, NULL);CHKERRQ(ierr);
551   }
552   /* Close the partition over the mesh */
553   ierr = DMPlexCreatePartitionClosure(dm, cellPartSection, cellPart, &partSection, &part);CHKERRQ(ierr);
554   ierr = ISDestroy(&cellPart);CHKERRQ(ierr);
555   ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr);
556   /* Create new mesh */
557   ierr  = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr);
558   ierr  = DMSetDimension(*dmParallel, dim);CHKERRQ(ierr);
559   ierr  = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr);
560   pmesh = (DM_Plex*) (*dmParallel)->data;
561   /* Distribute sieve points and the global point numbering (replaces creating remote bases) */
562   ierr = PetscSFConvertPartition(partSF, partSection, part, &renumbering, &pointSF);CHKERRQ(ierr);
563   if (flg) {
564     ierr = PetscPrintf(comm, "Point Partition:\n");CHKERRQ(ierr);
565     ierr = PetscSectionView(partSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
566     ierr = ISView(part, NULL);CHKERRQ(ierr);
567     ierr = PetscSFView(pointSF, NULL);CHKERRQ(ierr);
568     ierr = PetscPrintf(comm, "Point Renumbering after partition:\n");CHKERRQ(ierr);
569     ierr = ISLocalToGlobalMappingView(renumbering, NULL);CHKERRQ(ierr);
570   }
571   ierr = PetscLogEventEnd(DMPLEX_Partition,dm,0,0,0);CHKERRQ(ierr);
572   ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
573   /* Distribute cone section */
574   ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr);
575   ierr = DMPlexGetConeSection(*dmParallel, &newConeSection);CHKERRQ(ierr);
576   ierr = PetscSFDistributeSection(pointSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr);
577   ierr = DMSetUp(*dmParallel);CHKERRQ(ierr);
578   {
579     PetscInt pStart, pEnd, p;
580 
581     ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr);
582     for (p = pStart; p < pEnd; ++p) {
583       PetscInt coneSize;
584       ierr               = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr);
585       pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize);
586     }
587   }
588   /* Communicate and renumber cones */
589   ierr = PetscSFCreateSectionSF(pointSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr);
590   ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
591   ierr = DMPlexGetCones(*dmParallel, &newCones);CHKERRQ(ierr);
592   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
593   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
594   ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr);
595   ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr);
596   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr);
597   if (flg) {
598     ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr);
599     ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
600     ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr);
601     ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
602     ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr);
603   }
604   ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr);
605   ierr = DMPlexGetConeOrientations(*dmParallel, &newCones);CHKERRQ(ierr);
606   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
607   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
608   ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr);
609   ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
610   /* Create supports and stratify sieve */
611   {
612     PetscInt pStart, pEnd;
613 
614     ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr);
615     ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr);
616   }
617   ierr = DMPlexSymmetrize(*dmParallel);CHKERRQ(ierr);
618   ierr = DMPlexStratify(*dmParallel);CHKERRQ(ierr);
619   /* Create point SF for parallel mesh */
620   ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
621   {
622     const PetscInt *leaves;
623     PetscSFNode    *remotePoints, *rowners, *lowners;
624     PetscInt        numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints;
625     PetscInt        pStart, pEnd;
626 
627     ierr = DMPlexGetChart(*dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
628     ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr);
629     ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr);
630     for (p=0; p<numRoots; p++) {
631       rowners[p].rank  = -1;
632       rowners[p].index = -1;
633     }
634     if (origCellPart) {
635       /* Make sure points in the original partition are not assigned to other procs */
636       const PetscInt *origPoints;
637 
638       ierr = DMPlexCreatePartitionClosure(dm, origCellPartSection, origCellPart, &origPartSection, &origPart);CHKERRQ(ierr);
639       ierr = ISGetIndices(origPart, &origPoints);CHKERRQ(ierr);
640       for (p = 0; p < numProcs; ++p) {
641         PetscInt dof, off, d;
642 
643         ierr = PetscSectionGetDof(origPartSection, p, &dof);CHKERRQ(ierr);
644         ierr = PetscSectionGetOffset(origPartSection, p, &off);CHKERRQ(ierr);
645         for (d = off; d < off+dof; ++d) {
646           rowners[origPoints[d]].rank = p;
647         }
648       }
649       ierr = ISRestoreIndices(origPart, &origPoints);CHKERRQ(ierr);
650       ierr = ISDestroy(&origPart);CHKERRQ(ierr);
651       ierr = PetscSectionDestroy(&origPartSection);CHKERRQ(ierr);
652     }
653     ierr = ISDestroy(&origCellPart);CHKERRQ(ierr);
654     ierr = PetscSectionDestroy(&origCellPartSection);CHKERRQ(ierr);
655 
656     ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
657     ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
658     for (p = 0; p < numLeaves; ++p) {
659       if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */
660         lowners[p].rank  = rank;
661         lowners[p].index = leaves ? leaves[p] : p;
662       } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */
663         lowners[p].rank  = -2;
664         lowners[p].index = -2;
665       }
666     }
667     for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */
668       rowners[p].rank  = -3;
669       rowners[p].index = -3;
670     }
671     ierr = PetscSFReduceBegin(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
672     ierr = PetscSFReduceEnd(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
673     ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
674     ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
675     for (p = 0; p < numLeaves; ++p) {
676       if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed");
677       if (lowners[p].rank != rank) ++numGhostPoints;
678     }
679     ierr = PetscMalloc1(numGhostPoints,    &ghostPoints);CHKERRQ(ierr);
680     ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr);
681     for (p = 0, gp = 0; p < numLeaves; ++p) {
682       if (lowners[p].rank != rank) {
683         ghostPoints[gp]        = leaves ? leaves[p] : p;
684         remotePoints[gp].rank  = lowners[p].rank;
685         remotePoints[gp].index = lowners[p].index;
686         ++gp;
687       }
688     }
689     ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr);
690     ierr = PetscSFSetGraph((*dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
691     ierr = PetscSFSetFromOptions((*dmParallel)->sf);CHKERRQ(ierr);
692   }
693   pmesh->useCone    = mesh->useCone;
694   pmesh->useClosure = mesh->useClosure;
695   pmesh->useAnchors = mesh->useAnchors;
696   ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
697   /* Distribute Coordinates */
698   {
699     PetscSection     originalCoordSection, newCoordSection;
700     Vec              originalCoordinates, newCoordinates;
701     PetscInt         bs;
702     const char      *name;
703     const PetscReal *maxCell, *L;
704 
705     ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr);
706     ierr = DMGetCoordinateSection(*dmParallel, &newCoordSection);CHKERRQ(ierr);
707     ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr);
708     ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr);
709     ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr);
710     ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr);
711 
712     ierr = DMPlexDistributeField(dm, pointSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr);
713     ierr = DMSetCoordinatesLocal(*dmParallel, newCoordinates);CHKERRQ(ierr);
714     ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr);
715     ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr);
716     ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr);
717     ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr);
718     if (L) {ierr = DMSetPeriodicity(*dmParallel, maxCell, L);CHKERRQ(ierr);}
719   }
720   /* Distribute labels */
721   ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
722   {
723     DMLabel  next      = mesh->labels, newNext = pmesh->labels;
724     PetscInt numLabels = 0, l;
725 
726     /* Bcast number of labels */
727     while (next) {++numLabels; next = next->next;}
728     ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
729     next = mesh->labels;
730     for (l = 0; l < numLabels; ++l) {
731       DMLabel   labelNew;
732       PetscBool isdepth;
733 
734       /* Skip "depth" because it is recreated */
735       if (!rank) {ierr = PetscStrcmp(next->name, "depth", &isdepth);CHKERRQ(ierr);}
736       ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
737       if (isdepth) {if (!rank) next = next->next; continue;}
738       ierr = DMLabelDistribute(next, partSection, part, renumbering, &labelNew);CHKERRQ(ierr);
739       /* Insert into list */
740       if (newNext) newNext->next = labelNew;
741       else         pmesh->labels = labelNew;
742       newNext = labelNew;
743       if (!rank) next = next->next;
744     }
745   }
746   ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
747   /* Setup hybrid structure */
748   {
749     const PetscInt *gpoints;
750     PetscInt        depth, n, d;
751 
752     for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];}
753     ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr);
754     ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr);
755     ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr);
756     ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
757     for (d = 0; d <= dim; ++d) {
758       PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p;
759 
760       if (pmax < 0) continue;
761       ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr);
762       ierr = DMPlexGetDepthStratum(*dmParallel, d, NULL, &pEnd);CHKERRQ(ierr);
763       ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr);
764       for (p = 0; p < n; ++p) {
765         const PetscInt point = gpoints[p];
766 
767         if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax;
768       }
769       if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax;
770       else            pmesh->hybridPointMax[d] = -1;
771     }
772     ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr);
773   }
774   /* Set up tree */
775   {
776     DM              refTree;
777     PetscSection    origParentSection, newParentSection;
778     PetscInt        *origParents, *origChildIDs;
779 
780     ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr);
781     ierr = DMPlexSetReferenceTree(*dmParallel,refTree);CHKERRQ(ierr);
782     ierr = DMPlexGetTree(dm,&origParentSection,&origParents,&origChildIDs,NULL,NULL);CHKERRQ(ierr);
783     if (origParentSection) {
784       PetscInt        pStart, pEnd;
785       PetscInt        *newParents, *newChildIDs;
786       PetscInt        *remoteOffsetsParents, newParentSize;
787       PetscSF         parentSF;
788 
789       ierr = DMPlexGetChart(*dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
790       ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dmParallel),&newParentSection);CHKERRQ(ierr);
791       ierr = PetscSectionSetChart(newParentSection,pStart,pEnd);CHKERRQ(ierr);
792       ierr = PetscSFDistributeSection(pointSF, origParentSection, &remoteOffsetsParents, newParentSection);CHKERRQ(ierr);
793       ierr = PetscSFCreateSectionSF(pointSF, origParentSection, remoteOffsetsParents, newParentSection, &parentSF);CHKERRQ(ierr);
794       ierr = PetscSectionGetStorageSize(newParentSection,&newParentSize);CHKERRQ(ierr);
795       ierr = PetscMalloc2(newParentSize,&newParents,newParentSize,&newChildIDs);CHKERRQ(ierr);
796       ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
797       ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
798       ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
799       ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
800       ierr = ISGlobalToLocalMappingApplyBlock(renumbering,IS_GTOLM_MASK, newParentSize, newParents, NULL, newParents);CHKERRQ(ierr);
801       ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-parents_view", &flg);CHKERRQ(ierr);
802       if (flg) {
803         ierr = PetscPrintf(comm, "Serial Parent Section: \n");CHKERRQ(ierr);
804         ierr = PetscSectionView(origParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
805         ierr = PetscPrintf(comm, "Parallel Parent Section: \n");CHKERRQ(ierr);
806         ierr = PetscSectionView(newParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
807         ierr = PetscSFView(parentSF, NULL);CHKERRQ(ierr);
808       }
809       ierr = DMPlexSetTree(*dmParallel,newParentSection,newParents,newChildIDs);CHKERRQ(ierr);
810       ierr = PetscSectionDestroy(&newParentSection);CHKERRQ(ierr);
811       ierr = PetscFree2(newParents,newChildIDs);CHKERRQ(ierr);
812       ierr = PetscSFDestroy(&parentSF);CHKERRQ(ierr);
813     }
814   }
815   /* Cleanup Partition */
816   ierr = ISLocalToGlobalMappingDestroy(&renumbering);CHKERRQ(ierr);
817   ierr = PetscSFDestroy(&partSF);CHKERRQ(ierr);
818   ierr = PetscSectionDestroy(&partSection);CHKERRQ(ierr);
819   ierr = ISDestroy(&part);CHKERRQ(ierr);
820   /* Copy BC */
821   ierr = DMPlexCopyBoundary(dm, *dmParallel);CHKERRQ(ierr);
822   /* Cleanup */
823   if (sf) {*sf = pointSF;}
824   else    {ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr);}
825   ierr = DMSetFromOptions(*dmParallel);CHKERRQ(ierr);
826   ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
827   PetscFunctionReturn(0);
828 }
829