xref: /petsc/src/dm/impls/plex/plexdistribute.c (revision 8b0b4c70077becf5cb5aa61295a890b450cc9c47)
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__ "DMPlexSetAdjacencyUseConstraints"
124 /*@
125   DMPlexSetAdjacencyUseConstraints - Define adjacency in the mesh using the point-to-point constraints.
126 
127   Input Parameters:
128 + dm      - The DM object
129 - useConstraints - 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(), DMPlexSetConstraints()
134 @*/
135 PetscErrorCode DMPlexSetAdjacencyUseConstraints(DM dm, PetscBool useConstraints)
136 {
137   DM_Plex *mesh = (DM_Plex *) dm->data;
138 
139   PetscFunctionBegin;
140   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
141   mesh->useConstraints = useConstraints;
142   PetscFunctionReturn(0);
143 }
144 
145 #undef __FUNCT__
146 #define __FUNCT__ "DMPlexGetAdjacencyUseConstraints"
147 /*@
148   DMPlexGetAdjacencyUseConstraints - 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 . useConstraints - 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: DMPlexSetAdjacencyUseConstraints(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetConstraints()
159 @*/
160 PetscErrorCode DMPlexGetAdjacencyUseConstraints(DM dm, PetscBool *useConstraints)
161 {
162   DM_Plex *mesh = (DM_Plex *) dm->data;
163 
164   PetscFunctionBegin;
165   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
166   PetscValidIntPointer(useConstraints, 2);
167   *useConstraints = mesh->useConstraints;
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 useConstraints, 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 (useConstraints) {
270     ierr = DMPlexGetConstraints(dm,&aSec,&aIS);CHKERRQ(ierr);
271     if (aSec) {
272       ierr = PetscSectionGetMaxDof(aSec,&maxAnchors);CHKERRQ(ierr);
273       ierr = PetscSectionGetChart(aSec,&aStart,&aEnd);CHKERRQ(ierr);
274       ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr);
275     }
276   }
277   if (!*adj) {
278     PetscInt depth, maxConeSize, maxSupportSize;
279 
280     ierr  = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
281     ierr  = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
282     asiz  = PetscPowInt(maxConeSize, depth+1) * PetscPowInt(maxSupportSize, depth+1) + 1;
283     asiz *= maxAnchors;
284     ierr  = PetscMalloc1(asiz,adj);CHKERRQ(ierr);
285   }
286   if (*adjSize < 0) *adjSize = asiz;
287   maxAdjSize = *adjSize;
288   if (useTransitiveClosure) {
289     ierr = DMPlexGetAdjacency_Transitive_Internal(dm, p, useCone, adjSize, *adj);CHKERRQ(ierr);
290   } else if (useCone) {
291     ierr = DMPlexGetAdjacency_Cone_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
292   } else {
293     ierr = DMPlexGetAdjacency_Support_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
294   }
295   if (useConstraints && aSec) {
296     PetscInt origSize = *adjSize;
297     PetscInt numAdj = origSize;
298     PetscInt i = 0, j;
299     PetscInt *orig = *adj;
300 
301     while (i < origSize) {
302       PetscInt p = orig[i];
303       PetscInt aDof = 0;
304 
305       if (p >= aStart && p < aEnd) {
306         ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr);
307       }
308       if (aDof) {
309         PetscInt aOff;
310         PetscInt s, q;
311 
312         for (j = i + 1; j < numAdj; j++) {
313           orig[j - 1] = orig[j];
314         }
315         origSize--;
316         numAdj--;
317         ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr);
318         for (s = 0; s < aDof; ++s) {
319           for (q = 0; q < numAdj || (orig[numAdj++] = anchors[aOff+s],0); ++q) {
320             if (anchors[aOff+s] == orig[q]) break;
321           }
322           if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
323         }
324       }
325       else {
326         i++;
327       }
328     }
329     *adjSize = numAdj;
330     ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr);
331   }
332   PetscFunctionReturn(0);
333 }
334 
335 #undef __FUNCT__
336 #define __FUNCT__ "DMPlexGetAdjacency"
337 /*@
338   DMPlexGetAdjacency - Return all points adjacent to the given point
339 
340   Input Parameters:
341 + dm - The DM object
342 . p  - The point
343 . adjSize - The maximum size of adj if it is non-NULL, or PETSC_DETERMINE
344 - adj - Either NULL so that the array is allocated, or an existing array with size adjSize
345 
346   Output Parameters:
347 + adjSize - The number of adjacent points
348 - adj - The adjacent points
349 
350   Level: advanced
351 
352   Notes: The user must PetscFree the adj array if it was not passed in.
353 
354 .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexDistribute(), DMCreateMatrix(), DMPlexPreallocateOperator()
355 @*/
356 PetscErrorCode DMPlexGetAdjacency(DM dm, PetscInt p, PetscInt *adjSize, PetscInt *adj[])
357 {
358   DM_Plex       *mesh = (DM_Plex *) dm->data;
359   PetscErrorCode ierr;
360 
361   PetscFunctionBeginHot;
362   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
363   PetscValidPointer(adjSize,3);
364   PetscValidPointer(adj,4);
365   ierr = DMPlexGetAdjacency_Internal(dm, p, mesh->useCone, mesh->useClosure, mesh->useConstraints, adjSize, adj);CHKERRQ(ierr);
366   PetscFunctionReturn(0);
367 }
368 
369 #undef __FUNCT__
370 #define __FUNCT__ "DMPlexDistributeField"
371 /*@
372   DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
373 
374   Collective on DM
375 
376   Input Parameters:
377 + dm - The DMPlex object
378 . pointSF - The PetscSF describing the communication pattern
379 . originalSection - The PetscSection for existing data layout
380 - originalVec - The existing data
381 
382   Output Parameters:
383 + newSection - The PetscSF describing the new data layout
384 - newVec - The new data
385 
386   Level: developer
387 
388 .seealso: DMPlexDistribute(), DMPlexDistributeData()
389 @*/
390 PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec)
391 {
392   PetscSF        fieldSF;
393   PetscInt      *remoteOffsets, fieldSize;
394   PetscScalar   *originalValues, *newValues;
395   PetscErrorCode ierr;
396 
397   PetscFunctionBegin;
398   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
399   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
400 
401   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
402   ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr);
403   ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr);
404 
405   ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr);
406   ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr);
407   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
408   ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
409   ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
410   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
411   ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr);
412   ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr);
413   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
414   PetscFunctionReturn(0);
415 }
416 
417 #undef __FUNCT__
418 #define __FUNCT__ "DMPlexDistributeData"
419 /*@
420   DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
421 
422   Collective on DM
423 
424   Input Parameters:
425 + dm - The DMPlex object
426 . pointSF - The PetscSF describing the communication pattern
427 . originalSection - The PetscSection for existing data layout
428 . datatype - The type of data
429 - originalData - The existing data
430 
431   Output Parameters:
432 + newSection - The PetscSF describing the new data layout
433 - newData - The new data
434 
435   Level: developer
436 
437 .seealso: DMPlexDistribute(), DMPlexDistributeField()
438 @*/
439 PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData)
440 {
441   PetscSF        fieldSF;
442   PetscInt      *remoteOffsets, fieldSize;
443   PetscMPIInt    dataSize;
444   PetscErrorCode ierr;
445 
446   PetscFunctionBegin;
447   ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
448   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
449 
450   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
451   ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr);
452   ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr);
453 
454   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
455   ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
456   ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
457   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
458   ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
459   PetscFunctionReturn(0);
460 }
461 
462 #undef __FUNCT__
463 #define __FUNCT__ "DMPlexDistribute"
464 /*@C
465   DMPlexDistribute - Distributes the mesh and any associated sections.
466 
467   Not Collective
468 
469   Input Parameter:
470 + dm  - The original DMPlex object
471 . partitioner - The partitioning package, or NULL for the default
472 - overlap - The overlap of partitions, 0 is the default
473 
474   Output Parameter:
475 + sf - The PetscSF used for point distribution
476 - parallelMesh - The distributed DMPlex object, or NULL
477 
478   Note: If the mesh was not distributed, the return value is NULL.
479 
480   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
481   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
482   representation on the mesh.
483 
484   Level: intermediate
485 
486 .keywords: mesh, elements
487 .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
488 @*/
489 PetscErrorCode DMPlexDistribute(DM dm, const char partitioner[], PetscInt overlap, PetscSF *sf, DM *dmParallel)
490 {
491   DM_Plex               *mesh   = (DM_Plex*) dm->data, *pmesh;
492   MPI_Comm               comm;
493   const PetscInt         height = 0;
494   PetscInt               dim, numRemoteRanks;
495   IS                     origCellPart,        origPart,        cellPart,        part;
496   PetscSection           origCellPartSection, origPartSection, cellPartSection, partSection;
497   PetscSFNode           *remoteRanks;
498   PetscSF                partSF, pointSF, coneSF;
499   ISLocalToGlobalMapping renumbering;
500   PetscSection           originalConeSection, newConeSection;
501   PetscInt              *remoteOffsets;
502   PetscInt              *cones, *newCones, newConesSize;
503   PetscBool              flg;
504   PetscMPIInt            rank, numProcs, p;
505   PetscErrorCode         ierr;
506 
507   PetscFunctionBegin;
508   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
509   if (sf) PetscValidPointer(sf,4);
510   PetscValidPointer(dmParallel,5);
511 
512   ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
513   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
514   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
515   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
516 
517   *dmParallel = NULL;
518   if (numProcs == 1) PetscFunctionReturn(0);
519 
520   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
521   /* Create cell partition - We need to rewrite to use IS, use the MatPartition stuff */
522   ierr = PetscLogEventBegin(DMPLEX_Partition,dm,0,0,0);CHKERRQ(ierr);
523   if (overlap > 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Overlap > 1 not yet implemented");
524   ierr = DMPlexCreatePartition(dm, partitioner, height, overlap > 0 ? PETSC_TRUE : PETSC_FALSE, &cellPartSection, &cellPart, &origCellPartSection, &origCellPart);CHKERRQ(ierr);
525   /* Create SF assuming a serial partition for all processes: Could check for IS length here */
526   if (!rank) numRemoteRanks = numProcs;
527   else       numRemoteRanks = 0;
528   ierr = PetscMalloc1(numRemoteRanks, &remoteRanks);CHKERRQ(ierr);
529   for (p = 0; p < numRemoteRanks; ++p) {
530     remoteRanks[p].rank  = p;
531     remoteRanks[p].index = 0;
532   }
533   ierr = PetscSFCreate(comm, &partSF);CHKERRQ(ierr);
534   ierr = PetscSFSetGraph(partSF, 1, numRemoteRanks, NULL, PETSC_OWN_POINTER, remoteRanks, PETSC_OWN_POINTER);CHKERRQ(ierr);
535   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
536   if (flg) {
537     ierr = PetscPrintf(comm, "Cell Partition:\n");CHKERRQ(ierr);
538     ierr = PetscSectionView(cellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
539     ierr = ISView(cellPart, NULL);CHKERRQ(ierr);
540     if (origCellPart) {
541       ierr = PetscPrintf(comm, "Original Cell Partition:\n");CHKERRQ(ierr);
542       ierr = PetscSectionView(origCellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
543       ierr = ISView(origCellPart, NULL);CHKERRQ(ierr);
544     }
545     ierr = PetscSFView(partSF, NULL);CHKERRQ(ierr);
546   }
547   /* Close the partition over the mesh */
548   ierr = DMPlexCreatePartitionClosure(dm, cellPartSection, cellPart, &partSection, &part);CHKERRQ(ierr);
549   ierr = ISDestroy(&cellPart);CHKERRQ(ierr);
550   ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr);
551   /* Create new mesh */
552   ierr  = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr);
553   ierr  = DMSetDimension(*dmParallel, dim);CHKERRQ(ierr);
554   ierr  = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr);
555   pmesh = (DM_Plex*) (*dmParallel)->data;
556   /* Distribute sieve points and the global point numbering (replaces creating remote bases) */
557   ierr = PetscSFConvertPartition(partSF, partSection, part, &renumbering, &pointSF);CHKERRQ(ierr);
558   if (flg) {
559     ierr = PetscPrintf(comm, "Point Partition:\n");CHKERRQ(ierr);
560     ierr = PetscSectionView(partSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
561     ierr = ISView(part, NULL);CHKERRQ(ierr);
562     ierr = PetscSFView(pointSF, NULL);CHKERRQ(ierr);
563     ierr = PetscPrintf(comm, "Point Renumbering after partition:\n");CHKERRQ(ierr);
564     ierr = ISLocalToGlobalMappingView(renumbering, NULL);CHKERRQ(ierr);
565   }
566   ierr = PetscLogEventEnd(DMPLEX_Partition,dm,0,0,0);CHKERRQ(ierr);
567   ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
568   /* Distribute cone section */
569   ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr);
570   ierr = DMPlexGetConeSection(*dmParallel, &newConeSection);CHKERRQ(ierr);
571   ierr = PetscSFDistributeSection(pointSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr);
572   ierr = DMSetUp(*dmParallel);CHKERRQ(ierr);
573   {
574     PetscInt pStart, pEnd, p;
575 
576     ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr);
577     for (p = pStart; p < pEnd; ++p) {
578       PetscInt coneSize;
579       ierr               = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr);
580       pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize);
581     }
582   }
583   /* Communicate and renumber cones */
584   ierr = PetscSFCreateSectionSF(pointSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr);
585   ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
586   ierr = DMPlexGetCones(*dmParallel, &newCones);CHKERRQ(ierr);
587   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
588   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
589   ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr);
590   ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr);
591   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr);
592   if (flg) {
593     ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr);
594     ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
595     ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr);
596     ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
597     ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr);
598   }
599   ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr);
600   ierr = DMPlexGetConeOrientations(*dmParallel, &newCones);CHKERRQ(ierr);
601   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
602   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
603   ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr);
604   ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
605   /* Create supports and stratify sieve */
606   {
607     PetscInt pStart, pEnd;
608 
609     ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr);
610     ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr);
611   }
612   ierr = DMPlexSymmetrize(*dmParallel);CHKERRQ(ierr);
613   ierr = DMPlexStratify(*dmParallel);CHKERRQ(ierr);
614   /* Create point SF for parallel mesh */
615   ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
616   {
617     const PetscInt *leaves;
618     PetscSFNode    *remotePoints, *rowners, *lowners;
619     PetscInt        numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints;
620     PetscInt        pStart, pEnd;
621 
622     ierr = DMPlexGetChart(*dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
623     ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr);
624     ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr);
625     for (p=0; p<numRoots; p++) {
626       rowners[p].rank  = -1;
627       rowners[p].index = -1;
628     }
629     if (origCellPart) {
630       /* Make sure points in the original partition are not assigned to other procs */
631       const PetscInt *origPoints;
632 
633       ierr = DMPlexCreatePartitionClosure(dm, origCellPartSection, origCellPart, &origPartSection, &origPart);CHKERRQ(ierr);
634       ierr = ISGetIndices(origPart, &origPoints);CHKERRQ(ierr);
635       for (p = 0; p < numProcs; ++p) {
636         PetscInt dof, off, d;
637 
638         ierr = PetscSectionGetDof(origPartSection, p, &dof);CHKERRQ(ierr);
639         ierr = PetscSectionGetOffset(origPartSection, p, &off);CHKERRQ(ierr);
640         for (d = off; d < off+dof; ++d) {
641           rowners[origPoints[d]].rank = p;
642         }
643       }
644       ierr = ISRestoreIndices(origPart, &origPoints);CHKERRQ(ierr);
645       ierr = ISDestroy(&origPart);CHKERRQ(ierr);
646       ierr = PetscSectionDestroy(&origPartSection);CHKERRQ(ierr);
647     }
648     ierr = ISDestroy(&origCellPart);CHKERRQ(ierr);
649     ierr = PetscSectionDestroy(&origCellPartSection);CHKERRQ(ierr);
650 
651     ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
652     ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
653     for (p = 0; p < numLeaves; ++p) {
654       if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */
655         lowners[p].rank  = rank;
656         lowners[p].index = leaves ? leaves[p] : p;
657       } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */
658         lowners[p].rank  = -2;
659         lowners[p].index = -2;
660       }
661     }
662     for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */
663       rowners[p].rank  = -3;
664       rowners[p].index = -3;
665     }
666     ierr = PetscSFReduceBegin(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
667     ierr = PetscSFReduceEnd(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
668     ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
669     ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
670     for (p = 0; p < numLeaves; ++p) {
671       if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed");
672       if (lowners[p].rank != rank) ++numGhostPoints;
673     }
674     ierr = PetscMalloc1(numGhostPoints,    &ghostPoints);CHKERRQ(ierr);
675     ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr);
676     for (p = 0, gp = 0; p < numLeaves; ++p) {
677       if (lowners[p].rank != rank) {
678         ghostPoints[gp]        = leaves ? leaves[p] : p;
679         remotePoints[gp].rank  = lowners[p].rank;
680         remotePoints[gp].index = lowners[p].index;
681         ++gp;
682       }
683     }
684     ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr);
685     ierr = PetscSFSetGraph((*dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
686     ierr = PetscSFSetFromOptions((*dmParallel)->sf);CHKERRQ(ierr);
687   }
688   pmesh->useCone    = mesh->useCone;
689   pmesh->useClosure = mesh->useClosure;
690   ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
691   /* Distribute Coordinates */
692   {
693     PetscSection     originalCoordSection, newCoordSection;
694     Vec              originalCoordinates, newCoordinates;
695     PetscInt         bs;
696     const char      *name;
697     const PetscReal *maxCell, *L;
698 
699     ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr);
700     ierr = DMGetCoordinateSection(*dmParallel, &newCoordSection);CHKERRQ(ierr);
701     ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr);
702     ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr);
703     ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr);
704     ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr);
705 
706     ierr = DMPlexDistributeField(dm, pointSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr);
707     ierr = DMSetCoordinatesLocal(*dmParallel, newCoordinates);CHKERRQ(ierr);
708     ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr);
709     ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr);
710     ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr);
711     ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr);
712     if (L) {ierr = DMSetPeriodicity(*dmParallel, maxCell, L);CHKERRQ(ierr);}
713   }
714   /* Distribute labels */
715   ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
716   {
717     DMLabel  next      = mesh->labels, newNext = pmesh->labels;
718     PetscInt numLabels = 0, l;
719 
720     /* Bcast number of labels */
721     while (next) {++numLabels; next = next->next;}
722     ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
723     next = mesh->labels;
724     for (l = 0; l < numLabels; ++l) {
725       DMLabel   labelNew;
726       PetscBool isdepth;
727 
728       /* Skip "depth" because it is recreated */
729       if (!rank) {ierr = PetscStrcmp(next->name, "depth", &isdepth);CHKERRQ(ierr);}
730       ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
731       if (isdepth) {if (!rank) next = next->next; continue;}
732       ierr = DMLabelDistribute(next, partSection, part, renumbering, &labelNew);CHKERRQ(ierr);
733       /* Insert into list */
734       if (newNext) newNext->next = labelNew;
735       else         pmesh->labels = labelNew;
736       newNext = labelNew;
737       if (!rank) next = next->next;
738     }
739   }
740   ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
741   /* Setup hybrid structure */
742   {
743     const PetscInt *gpoints;
744     PetscInt        depth, n, d;
745 
746     for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];}
747     ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr);
748     ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr);
749     ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr);
750     ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
751     for (d = 0; d <= dim; ++d) {
752       PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p;
753 
754       if (pmax < 0) continue;
755       ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr);
756       ierr = DMPlexGetDepthStratum(*dmParallel, d, NULL, &pEnd);CHKERRQ(ierr);
757       ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr);
758       for (p = 0; p < n; ++p) {
759         const PetscInt point = gpoints[p];
760 
761         if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax;
762       }
763       if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax;
764       else            pmesh->hybridPointMax[d] = -1;
765     }
766     ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr);
767   }
768   /* Cleanup Partition */
769   ierr = ISLocalToGlobalMappingDestroy(&renumbering);CHKERRQ(ierr);
770   ierr = PetscSFDestroy(&partSF);CHKERRQ(ierr);
771   ierr = PetscSectionDestroy(&partSection);CHKERRQ(ierr);
772   ierr = ISDestroy(&part);CHKERRQ(ierr);
773   /* Copy BC */
774   ierr = DMPlexCopyBoundary(dm, *dmParallel);CHKERRQ(ierr);
775   /* Cleanup */
776   if (sf) {*sf = pointSF;}
777   else    {ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr);}
778   ierr = DMSetFromOptions(*dmParallel);CHKERRQ(ierr);
779   ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
780   PetscFunctionReturn(0);
781 }
782