Actual source code: plexdistribute.c

petsc-master 2017-01-18
Report Typos and Errors
  1:  #include <petsc/private/dmpleximpl.h>
  2:  #include <petsc/private/dmlabelimpl.h>

  4: /*@
  5:   DMPlexSetAdjacencyUseCone - Define adjacency in the mesh using either the cone or the support first

  7:   Input Parameters:
  8: + dm      - The DM object
  9: - useCone - Flag to use the cone first

 11:   Level: intermediate

 13:   Notes:
 14: $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
 15: $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
 16: $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE

 18: .seealso: DMPlexGetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
 19: @*/
 20: PetscErrorCode DMPlexSetAdjacencyUseCone(DM dm, PetscBool useCone)
 21: {
 22:   DM_Plex *mesh = (DM_Plex *) dm->data;

 26:   mesh->useCone = useCone;
 27:   return(0);
 28: }

 30: /*@
 31:   DMPlexGetAdjacencyUseCone - Query whether adjacency in the mesh uses the cone or the support first

 33:   Input Parameter:
 34: . dm      - The DM object

 36:   Output Parameter:
 37: . useCone - Flag to use the cone first

 39:   Level: intermediate

 41:   Notes:
 42: $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
 43: $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
 44: $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE

 46: .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
 47: @*/
 48: PetscErrorCode DMPlexGetAdjacencyUseCone(DM dm, PetscBool *useCone)
 49: {
 50:   DM_Plex *mesh = (DM_Plex *) dm->data;

 55:   *useCone = mesh->useCone;
 56:   return(0);
 57: }

 59: /*@
 60:   DMPlexSetAdjacencyUseClosure - Define adjacency in the mesh using the transitive closure

 62:   Input Parameters:
 63: + dm      - The DM object
 64: - useClosure - Flag to use the closure

 66:   Level: intermediate

 68:   Notes:
 69: $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
 70: $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
 71: $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE

 73: .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
 74: @*/
 75: PetscErrorCode DMPlexSetAdjacencyUseClosure(DM dm, PetscBool useClosure)
 76: {
 77:   DM_Plex *mesh = (DM_Plex *) dm->data;

 81:   mesh->useClosure = useClosure;
 82:   return(0);
 83: }

 85: /*@
 86:   DMPlexGetAdjacencyUseClosure - Query whether adjacency in the mesh uses the transitive closure

 88:   Input Parameter:
 89: . dm      - The DM object

 91:   Output Parameter:
 92: . useClosure - Flag to use the closure

 94:   Level: intermediate

 96:   Notes:
 97: $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
 98: $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
 99: $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE

101: .seealso: DMPlexSetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
102: @*/
103: PetscErrorCode DMPlexGetAdjacencyUseClosure(DM dm, PetscBool *useClosure)
104: {
105:   DM_Plex *mesh = (DM_Plex *) dm->data;

110:   *useClosure = mesh->useClosure;
111:   return(0);
112: }

114: /*@
115:   DMPlexSetAdjacencyUseAnchors - Define adjacency in the mesh using the point-to-point constraints.

117:   Input Parameters:
118: + dm      - The DM object
119: - useAnchors - Flag to use the constraints.  If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place.

121:   Level: intermediate

123: .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
124: @*/
125: PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM dm, PetscBool useAnchors)
126: {
127:   DM_Plex *mesh = (DM_Plex *) dm->data;

131:   mesh->useAnchors = useAnchors;
132:   return(0);
133: }

135: /*@
136:   DMPlexGetAdjacencyUseAnchors - Query whether adjacency in the mesh uses the point-to-point constraints.

138:   Input Parameter:
139: . dm      - The DM object

141:   Output Parameter:
142: . useAnchors - Flag to use the closure.  If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place.

144:   Level: intermediate

146: .seealso: DMPlexSetAdjacencyUseAnchors(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
147: @*/
148: PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM dm, PetscBool *useAnchors)
149: {
150:   DM_Plex *mesh = (DM_Plex *) dm->data;

155:   *useAnchors = mesh->useAnchors;
156:   return(0);
157: }

159: static PetscErrorCode DMPlexGetAdjacency_Cone_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
160: {
161:   const PetscInt *cone = NULL;
162:   PetscInt        numAdj = 0, maxAdjSize = *adjSize, coneSize, c;
163:   PetscErrorCode  ierr;

166:   DMPlexGetConeSize(dm, p, &coneSize);
167:   DMPlexGetCone(dm, p, &cone);
168:   for (c = 0; c <= coneSize; ++c) {
169:     const PetscInt  point   = !c ? p : cone[c-1];
170:     const PetscInt *support = NULL;
171:     PetscInt        supportSize, s, q;

173:     DMPlexGetSupportSize(dm, point, &supportSize);
174:     DMPlexGetSupport(dm, point, &support);
175:     for (s = 0; s < supportSize; ++s) {
176:       for (q = 0; q < numAdj || ((void)(adj[numAdj++] = support[s]),0); ++q) {
177:         if (support[s] == adj[q]) break;
178:       }
179:       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
180:     }
181:   }
182:   *adjSize = numAdj;
183:   return(0);
184: }

186: static PetscErrorCode DMPlexGetAdjacency_Support_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
187: {
188:   const PetscInt *support = NULL;
189:   PetscInt        numAdj   = 0, maxAdjSize = *adjSize, supportSize, s;
190:   PetscErrorCode  ierr;

193:   DMPlexGetSupportSize(dm, p, &supportSize);
194:   DMPlexGetSupport(dm, p, &support);
195:   for (s = 0; s <= supportSize; ++s) {
196:     const PetscInt  point = !s ? p : support[s-1];
197:     const PetscInt *cone  = NULL;
198:     PetscInt        coneSize, c, q;

200:     DMPlexGetConeSize(dm, point, &coneSize);
201:     DMPlexGetCone(dm, point, &cone);
202:     for (c = 0; c < coneSize; ++c) {
203:       for (q = 0; q < numAdj || ((void)(adj[numAdj++] = cone[c]),0); ++q) {
204:         if (cone[c] == adj[q]) break;
205:       }
206:       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
207:     }
208:   }
209:   *adjSize = numAdj;
210:   return(0);
211: }

213: static PetscErrorCode DMPlexGetAdjacency_Transitive_Internal(DM dm, PetscInt p, PetscBool useClosure, PetscInt *adjSize, PetscInt adj[])
214: {
215:   PetscInt      *star = NULL;
216:   PetscInt       numAdj = 0, maxAdjSize = *adjSize, starSize, s;

220:   DMPlexGetTransitiveClosure(dm, p, useClosure, &starSize, &star);
221:   for (s = 0; s < starSize*2; s += 2) {
222:     const PetscInt *closure = NULL;
223:     PetscInt        closureSize, c, q;

225:     DMPlexGetTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);
226:     for (c = 0; c < closureSize*2; c += 2) {
227:       for (q = 0; q < numAdj || ((void)(adj[numAdj++] = closure[c]),0); ++q) {
228:         if (closure[c] == adj[q]) break;
229:       }
230:       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
231:     }
232:     DMPlexRestoreTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);
233:   }
234:   DMPlexRestoreTransitiveClosure(dm, p, useClosure, &starSize, &star);
235:   *adjSize = numAdj;
236:   return(0);
237: }

239: PetscErrorCode DMPlexGetAdjacency_Internal(DM dm, PetscInt p, PetscBool useCone, PetscBool useTransitiveClosure, PetscBool useAnchors, PetscInt *adjSize, PetscInt *adj[])
240: {
241:   static PetscInt asiz = 0;
242:   PetscInt maxAnchors = 1;
243:   PetscInt aStart = -1, aEnd = -1;
244:   PetscInt maxAdjSize;
245:   PetscSection aSec = NULL;
246:   IS aIS = NULL;
247:   const PetscInt *anchors;
248:   PetscErrorCode  ierr;

251:   if (useAnchors) {
252:     DMPlexGetAnchors(dm,&aSec,&aIS);
253:     if (aSec) {
254:       PetscSectionGetMaxDof(aSec,&maxAnchors);
255:       maxAnchors = PetscMax(1,maxAnchors);
256:       PetscSectionGetChart(aSec,&aStart,&aEnd);
257:       ISGetIndices(aIS,&anchors);
258:     }
259:   }
260:   if (!*adj) {
261:     PetscInt depth, coneSeries, supportSeries, maxC, maxS, pStart, pEnd;

263:     DMPlexGetChart(dm, &pStart,&pEnd);
264:     DMPlexGetDepth(dm, &depth);
265:     DMPlexGetMaxSizes(dm, &maxC, &maxS);
266:     coneSeries    = (maxC > 1) ? ((PetscPowInt(maxC,depth+1)-1)/(maxC-1)) : depth+1;
267:     supportSeries = (maxS > 1) ? ((PetscPowInt(maxS,depth+1)-1)/(maxS-1)) : depth+1;
268:     asiz  = PetscMax(PetscPowInt(maxS,depth)*coneSeries,PetscPowInt(maxC,depth)*supportSeries);
269:     asiz *= maxAnchors;
270:     asiz  = PetscMin(asiz,pEnd-pStart);
271:     PetscMalloc1(asiz,adj);
272:   }
273:   if (*adjSize < 0) *adjSize = asiz;
274:   maxAdjSize = *adjSize;
275:   if (useTransitiveClosure) {
276:     DMPlexGetAdjacency_Transitive_Internal(dm, p, useCone, adjSize, *adj);
277:   } else if (useCone) {
278:     DMPlexGetAdjacency_Cone_Internal(dm, p, adjSize, *adj);
279:   } else {
280:     DMPlexGetAdjacency_Support_Internal(dm, p, adjSize, *adj);
281:   }
282:   if (useAnchors && aSec) {
283:     PetscInt origSize = *adjSize;
284:     PetscInt numAdj = origSize;
285:     PetscInt i = 0, j;
286:     PetscInt *orig = *adj;

288:     while (i < origSize) {
289:       PetscInt p = orig[i];
290:       PetscInt aDof = 0;

292:       if (p >= aStart && p < aEnd) {
293:         PetscSectionGetDof(aSec,p,&aDof);
294:       }
295:       if (aDof) {
296:         PetscInt aOff;
297:         PetscInt s, q;

299:         for (j = i + 1; j < numAdj; j++) {
300:           orig[j - 1] = orig[j];
301:         }
302:         origSize--;
303:         numAdj--;
304:         PetscSectionGetOffset(aSec,p,&aOff);
305:         for (s = 0; s < aDof; ++s) {
306:           for (q = 0; q < numAdj || ((void)(orig[numAdj++] = anchors[aOff+s]),0); ++q) {
307:             if (anchors[aOff+s] == orig[q]) break;
308:           }
309:           if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
310:         }
311:       }
312:       else {
313:         i++;
314:       }
315:     }
316:     *adjSize = numAdj;
317:     ISRestoreIndices(aIS,&anchors);
318:   }
319:   return(0);
320: }

322: /*@
323:   DMPlexGetAdjacency - Return all points adjacent to the given point

325:   Input Parameters:
326: + dm - The DM object
327: . p  - The point
328: . adjSize - The maximum size of adj if it is non-NULL, or PETSC_DETERMINE
329: - adj - Either NULL so that the array is allocated, or an existing array with size adjSize

331:   Output Parameters:
332: + adjSize - The number of adjacent points
333: - adj - The adjacent points

335:   Level: advanced

337:   Notes: The user must PetscFree the adj array if it was not passed in.

339: .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexDistribute(), DMCreateMatrix(), DMPlexPreallocateOperator()
340: @*/
341: PetscErrorCode DMPlexGetAdjacency(DM dm, PetscInt p, PetscInt *adjSize, PetscInt *adj[])
342: {
343:   DM_Plex       *mesh = (DM_Plex *) dm->data;

350:   DMPlexGetAdjacency_Internal(dm, p, mesh->useCone, mesh->useClosure, mesh->useAnchors, adjSize, adj);
351:   return(0);
352: }

354: /*@
355:   DMPlexCreateTwoSidedProcessSF - Create an SF which just has process connectivity

357:   Collective on DM

359:   Input Parameters:
360: + dm      - The DM
361: - sfPoint - The PetscSF which encodes point connectivity

363:   Output Parameters:
364: + processRanks - A list of process neighbors, or NULL
365: - sfProcess    - An SF encoding the two-sided process connectivity, or NULL

367:   Level: developer

369: .seealso: PetscSFCreate(), DMPlexCreateProcessSF()
370: @*/
371: PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM dm, PetscSF sfPoint, PetscSection rootRankSection, IS rootRanks, PetscSection leafRankSection, IS leafRanks, IS *processRanks, PetscSF *sfProcess)
372: {
373:   const PetscSFNode *remotePoints;
374:   PetscInt          *localPointsNew;
375:   PetscSFNode       *remotePointsNew;
376:   const PetscInt    *nranks;
377:   PetscInt          *ranksNew;
378:   PetscBT            neighbors;
379:   PetscInt           pStart, pEnd, p, numLeaves, l, numNeighbors, n;
380:   PetscMPIInt        numProcs, proc, rank;
381:   PetscErrorCode     ierr;

388:   MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);
389:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);
390:   PetscSFGetGraph(sfPoint, NULL, &numLeaves, NULL, &remotePoints);
391:   PetscBTCreate(numProcs, &neighbors);
392:   PetscBTMemzero(numProcs, neighbors);
393:   /* Compute root-to-leaf process connectivity */
394:   PetscSectionGetChart(rootRankSection, &pStart, &pEnd);
395:   ISGetIndices(rootRanks, &nranks);
396:   for (p = pStart; p < pEnd; ++p) {
397:     PetscInt ndof, noff, n;

399:     PetscSectionGetDof(rootRankSection, p, &ndof);
400:     PetscSectionGetOffset(rootRankSection, p, &noff);
401:     for (n = 0; n < ndof; ++n) {PetscBTSet(neighbors, nranks[noff+n]);}
402:   }
403:   ISRestoreIndices(rootRanks, &nranks);
404:   /* Compute leaf-to-neighbor process connectivity */
405:   PetscSectionGetChart(leafRankSection, &pStart, &pEnd);
406:   ISGetIndices(leafRanks, &nranks);
407:   for (p = pStart; p < pEnd; ++p) {
408:     PetscInt ndof, noff, n;

410:     PetscSectionGetDof(leafRankSection, p, &ndof);
411:     PetscSectionGetOffset(leafRankSection, p, &noff);
412:     for (n = 0; n < ndof; ++n) {PetscBTSet(neighbors, nranks[noff+n]);}
413:   }
414:   ISRestoreIndices(leafRanks, &nranks);
415:   /* Compute leaf-to-root process connectivity */
416:   for (l = 0; l < numLeaves; ++l) {PetscBTSet(neighbors, remotePoints[l].rank);}
417:   /* Calculate edges */
418:   PetscBTClear(neighbors, rank);
419:   for(proc = 0, numNeighbors = 0; proc < numProcs; ++proc) {if (PetscBTLookup(neighbors, proc)) ++numNeighbors;}
420:   PetscMalloc1(numNeighbors, &ranksNew);
421:   PetscMalloc1(numNeighbors, &localPointsNew);
422:   PetscMalloc1(numNeighbors, &remotePointsNew);
423:   for(proc = 0, n = 0; proc < numProcs; ++proc) {
424:     if (PetscBTLookup(neighbors, proc)) {
425:       ranksNew[n]              = proc;
426:       localPointsNew[n]        = proc;
427:       remotePointsNew[n].index = rank;
428:       remotePointsNew[n].rank  = proc;
429:       ++n;
430:     }
431:   }
432:   PetscBTDestroy(&neighbors);
433:   if (processRanks) {ISCreateGeneral(PetscObjectComm((PetscObject)dm), numNeighbors, ranksNew, PETSC_OWN_POINTER, processRanks);}
434:   else              {PetscFree(ranksNew);}
435:   if (sfProcess) {
436:     PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess);
437:     PetscObjectSetName((PetscObject) *sfProcess, "Two-Sided Process SF");
438:     PetscSFSetFromOptions(*sfProcess);
439:     PetscSFSetGraph(*sfProcess, numProcs, numNeighbors, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);
440:   }
441:   return(0);
442: }

444: /*@
445:   DMPlexDistributeOwnership - Compute owner information for shared points. This basically gets two-sided for an SF.

447:   Collective on DM

449:   Input Parameter:
450: . dm - The DM

452:   Output Parameters:
453: + rootSection - The number of leaves for a given root point
454: . rootrank    - The rank of each edge into the root point
455: . leafSection - The number of processes sharing a given leaf point
456: - leafrank    - The rank of each process sharing a leaf point

458:   Level: developer

460: .seealso: DMPlexCreateOverlap()
461: @*/
462: PetscErrorCode DMPlexDistributeOwnership(DM dm, PetscSection rootSection, IS *rootrank, PetscSection leafSection, IS *leafrank)
463: {
464:   MPI_Comm        comm;
465:   PetscSF         sfPoint;
466:   const PetscInt *rootdegree;
467:   PetscInt       *myrank, *remoterank;
468:   PetscInt        pStart, pEnd, p, nedges;
469:   PetscMPIInt     rank;
470:   PetscErrorCode  ierr;

473:   PetscObjectGetComm((PetscObject) dm, &comm);
474:   MPI_Comm_rank(comm, &rank);
475:   DMPlexGetChart(dm, &pStart, &pEnd);
476:   DMGetPointSF(dm, &sfPoint);
477:   /* Compute number of leaves for each root */
478:   PetscObjectSetName((PetscObject) rootSection, "Root Section");
479:   PetscSectionSetChart(rootSection, pStart, pEnd);
480:   PetscSFComputeDegreeBegin(sfPoint, &rootdegree);
481:   PetscSFComputeDegreeEnd(sfPoint, &rootdegree);
482:   for (p = pStart; p < pEnd; ++p) {PetscSectionSetDof(rootSection, p, rootdegree[p-pStart]);}
483:   PetscSectionSetUp(rootSection);
484:   /* Gather rank of each leaf to root */
485:   PetscSectionGetStorageSize(rootSection, &nedges);
486:   PetscMalloc1(pEnd-pStart, &myrank);
487:   PetscMalloc1(nedges,  &remoterank);
488:   for (p = 0; p < pEnd-pStart; ++p) myrank[p] = rank;
489:   PetscSFGatherBegin(sfPoint, MPIU_INT, myrank, remoterank);
490:   PetscSFGatherEnd(sfPoint, MPIU_INT, myrank, remoterank);
491:   PetscFree(myrank);
492:   ISCreateGeneral(comm, nedges, remoterank, PETSC_OWN_POINTER, rootrank);
493:   /* Distribute remote ranks to leaves */
494:   PetscObjectSetName((PetscObject) leafSection, "Leaf Section");
495:   DMPlexDistributeFieldIS(dm, sfPoint, rootSection, *rootrank, leafSection, leafrank);
496:   return(0);
497: }

499: /*@C
500:   DMPlexCreateOverlap - Compute owner information for shared points. This basically gets two-sided for an SF.

502:   Collective on DM

504:   Input Parameters:
505: + dm          - The DM
506: . levels      - Number of overlap levels
507: . rootSection - The number of leaves for a given root point
508: . rootrank    - The rank of each edge into the root point
509: . leafSection - The number of processes sharing a given leaf point
510: - leafrank    - The rank of each process sharing a leaf point

512:   Output Parameters:
513: + ovLabel     - DMLabel containing remote overlap contributions as point/rank pairings

515:   Level: developer

517: .seealso: DMPlexDistributeOwnership(), DMPlexDistribute()
518: @*/
519: PetscErrorCode DMPlexCreateOverlap(DM dm, PetscInt levels, PetscSection rootSection, IS rootrank, PetscSection leafSection, IS leafrank, DMLabel *ovLabel)
520: {
521:   MPI_Comm           comm;
522:   DMLabel            ovAdjByRank; /* A DMLabel containing all points adjacent to shared points, separated by rank (value in label) */
523:   PetscSF            sfPoint, sfProc;
524:   const PetscSFNode *remote;
525:   const PetscInt    *local;
526:   const PetscInt    *nrank, *rrank;
527:   PetscInt          *adj = NULL;
528:   PetscInt           pStart, pEnd, p, sStart, sEnd, nleaves, l;
529:   PetscMPIInt        rank, numProcs;
530:   PetscBool          useCone, useClosure, flg;
531:   PetscErrorCode     ierr;

534:   PetscObjectGetComm((PetscObject) dm, &comm);
535:   MPI_Comm_size(comm, &numProcs);
536:   MPI_Comm_rank(comm, &rank);
537:   DMGetPointSF(dm, &sfPoint);
538:   DMPlexGetChart(dm, &pStart, &pEnd);
539:   PetscSectionGetChart(leafSection, &sStart, &sEnd);
540:   PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);
541:   DMLabelCreate("Overlap adjacency", &ovAdjByRank);
542:   /* Handle leaves: shared with the root point */
543:   for (l = 0; l < nleaves; ++l) {
544:     PetscInt adjSize = PETSC_DETERMINE, a;

546:     DMPlexGetAdjacency(dm, local[l], &adjSize, &adj);
547:     for (a = 0; a < adjSize; ++a) {DMLabelSetValue(ovAdjByRank, adj[a], remote[l].rank);}
548:   }
549:   ISGetIndices(rootrank, &rrank);
550:   ISGetIndices(leafrank, &nrank);
551:   /* Handle roots */
552:   for (p = pStart; p < pEnd; ++p) {
553:     PetscInt adjSize = PETSC_DETERMINE, neighbors = 0, noff, n, a;

555:     if ((p >= sStart) && (p < sEnd)) {
556:       /* Some leaves share a root with other leaves on different processes */
557:       PetscSectionGetDof(leafSection, p, &neighbors);
558:       if (neighbors) {
559:         PetscSectionGetOffset(leafSection, p, &noff);
560:         DMPlexGetAdjacency(dm, p, &adjSize, &adj);
561:         for (n = 0; n < neighbors; ++n) {
562:           const PetscInt remoteRank = nrank[noff+n];

564:           if (remoteRank == rank) continue;
565:           for (a = 0; a < adjSize; ++a) {DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);}
566:         }
567:       }
568:     }
569:     /* Roots are shared with leaves */
570:     PetscSectionGetDof(rootSection, p, &neighbors);
571:     if (!neighbors) continue;
572:     PetscSectionGetOffset(rootSection, p, &noff);
573:     DMPlexGetAdjacency(dm, p, &adjSize, &adj);
574:     for (n = 0; n < neighbors; ++n) {
575:       const PetscInt remoteRank = rrank[noff+n];

577:       if (remoteRank == rank) continue;
578:       for (a = 0; a < adjSize; ++a) {DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);}
579:     }
580:   }
581:   PetscFree(adj);
582:   ISRestoreIndices(rootrank, &rrank);
583:   ISRestoreIndices(leafrank, &nrank);
584:   /* Add additional overlap levels */
585:   for (l = 1; l < levels; l++) {
586:     /* Propagate point donations over SF to capture remote connections */
587:     DMPlexPartitionLabelPropagate(dm, ovAdjByRank);
588:     /* Add next level of point donations to the label */
589:     DMPlexPartitionLabelAdjacency(dm, ovAdjByRank);
590:   }
591:   /* We require the closure in the overlap */
592:   DMPlexGetAdjacencyUseCone(dm, &useCone);
593:   DMPlexGetAdjacencyUseClosure(dm, &useClosure);
594:   if (useCone || !useClosure) {
595:     DMPlexPartitionLabelClosure(dm, ovAdjByRank);
596:   }
597:   PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-overlap_view", &flg);
598:   if (flg) {
599:     DMLabelView(ovAdjByRank, PETSC_VIEWER_STDOUT_WORLD);
600:   }
601:   /* Make global process SF and invert sender to receiver label */
602:   {
603:     /* Build a global process SF */
604:     PetscSFNode *remoteProc;
605:     PetscMalloc1(numProcs, &remoteProc);
606:     for (p = 0; p < numProcs; ++p) {
607:       remoteProc[p].rank  = p;
608:       remoteProc[p].index = rank;
609:     }
610:     PetscSFCreate(comm, &sfProc);
611:     PetscObjectSetName((PetscObject) sfProc, "Process SF");
612:     PetscSFSetGraph(sfProc, numProcs, numProcs, NULL, PETSC_OWN_POINTER, remoteProc, PETSC_OWN_POINTER);
613:   }
614:   DMLabelCreate("Overlap label", ovLabel);
615:   DMPlexPartitionLabelInvert(dm, ovAdjByRank, sfProc, *ovLabel);
616:   /* Add owned points, except for shared local points */
617:   for (p = pStart; p < pEnd; ++p) {DMLabelSetValue(*ovLabel, p, rank);}
618:   for (l = 0; l < nleaves; ++l) {
619:     DMLabelClearValue(*ovLabel, local[l], rank);
620:     DMLabelSetValue(*ovLabel, remote[l].index, remote[l].rank);
621:   }
622:   /* Clean up */
623:   DMLabelDestroy(&ovAdjByRank);
624:   PetscSFDestroy(&sfProc);
625:   return(0);
626: }

628: PetscErrorCode DMPlexCreateOverlapMigrationSF(DM dm, PetscSF overlapSF, PetscSF *migrationSF)
629: {
630:   MPI_Comm           comm;
631:   PetscMPIInt        rank, numProcs;
632:   PetscInt           d, dim, p, pStart, pEnd, nroots, nleaves, newLeaves, point, numSharedPoints;
633:   PetscInt          *pointDepths, *remoteDepths, *ilocal;
634:   PetscInt          *depthRecv, *depthShift, *depthIdx;
635:   PetscSFNode       *iremote;
636:   PetscSF            pointSF;
637:   const PetscInt    *sharedLocal;
638:   const PetscSFNode *overlapRemote, *sharedRemote;
639:   PetscErrorCode     ierr;

643:   PetscObjectGetComm((PetscObject)dm, &comm);
644:   MPI_Comm_rank(comm, &rank);
645:   MPI_Comm_size(comm, &numProcs);
646:   DMGetDimension(dm, &dim);

648:   /* Before building the migration SF we need to know the new stratum offsets */
649:   PetscSFGetGraph(overlapSF, &nroots, &nleaves, NULL, &overlapRemote);
650:   PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);
651:   for (d=0; d<dim+1; d++) {
652:     DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
653:     for (p=pStart; p<pEnd; p++) pointDepths[p] = d;
654:   }
655:   for (p=0; p<nleaves; p++) remoteDepths[p] = -1;
656:   PetscSFBcastBegin(overlapSF, MPIU_INT, pointDepths, remoteDepths);
657:   PetscSFBcastEnd(overlapSF, MPIU_INT, pointDepths, remoteDepths);

659:   /* Count recevied points in each stratum and compute the internal strata shift */
660:   PetscMalloc3(dim+1, &depthRecv, dim+1, &depthShift, dim+1, &depthIdx);
661:   for (d=0; d<dim+1; d++) depthRecv[d]=0;
662:   for (p=0; p<nleaves; p++) depthRecv[remoteDepths[p]]++;
663:   depthShift[dim] = 0;
664:   for (d=0; d<dim; d++) depthShift[d] = depthRecv[dim];
665:   for (d=1; d<dim; d++) depthShift[d] += depthRecv[0];
666:   for (d=dim-2; d>0; d--) depthShift[d] += depthRecv[d+1];
667:   for (d=0; d<dim+1; d++) {
668:     DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
669:     depthIdx[d] = pStart + depthShift[d];
670:   }

672:   /* Form the overlap SF build an SF that describes the full overlap migration SF */
673:   DMPlexGetChart(dm, &pStart, &pEnd);
674:   newLeaves = pEnd - pStart + nleaves;
675:   PetscMalloc1(newLeaves, &ilocal);
676:   PetscMalloc1(newLeaves, &iremote);
677:   /* First map local points to themselves */
678:   for (d=0; d<dim+1; d++) {
679:     DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
680:     for (p=pStart; p<pEnd; p++) {
681:       point = p + depthShift[d];
682:       ilocal[point] = point;
683:       iremote[point].index = p;
684:       iremote[point].rank = rank;
685:       depthIdx[d]++;
686:     }
687:   }

689:   /* Add in the remote roots for currently shared points */
690:   DMGetPointSF(dm, &pointSF);
691:   PetscSFGetGraph(pointSF, NULL, &numSharedPoints, &sharedLocal, &sharedRemote);
692:   for (d=0; d<dim+1; d++) {
693:     DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
694:     for (p=0; p<numSharedPoints; p++) {
695:       if (pStart <= sharedLocal[p] && sharedLocal[p] < pEnd) {
696:         point = sharedLocal[p] + depthShift[d];
697:         iremote[point].index = sharedRemote[p].index;
698:         iremote[point].rank = sharedRemote[p].rank;
699:       }
700:     }
701:   }

703:   /* Now add the incoming overlap points */
704:   for (p=0; p<nleaves; p++) {
705:     point = depthIdx[remoteDepths[p]];
706:     ilocal[point] = point;
707:     iremote[point].index = overlapRemote[p].index;
708:     iremote[point].rank = overlapRemote[p].rank;
709:     depthIdx[remoteDepths[p]]++;
710:   }
711:   PetscFree2(pointDepths,remoteDepths);

713:   PetscSFCreate(comm, migrationSF);
714:   PetscObjectSetName((PetscObject) *migrationSF, "Overlap Migration SF");
715:   PetscSFSetFromOptions(*migrationSF);
716:   DMPlexGetChart(dm, &pStart, &pEnd);
717:   PetscSFSetGraph(*migrationSF, pEnd-pStart, newLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);

719:   PetscFree3(depthRecv, depthShift, depthIdx);
720:   return(0);
721: }

723: /*@
724:   DMPlexStratifyMigrationSF - Rearrange the leaves of a migration sf for stratification.

726:   Input Parameter:
727: + dm          - The DM
728: - sf          - A star forest with non-ordered leaves, usually defining a DM point migration

730:   Output Parameter:
731: . migrationSF - A star forest with added leaf indirection that ensures the resulting DM is stratified

733:   Level: developer

735: .seealso: DMPlexPartitionLabelCreateSF(), DMPlexDistribute(), DMPlexDistributeOverlap()
736: @*/
737: PetscErrorCode DMPlexStratifyMigrationSF(DM dm, PetscSF sf, PetscSF *migrationSF)
738: {
739:   MPI_Comm           comm;
740:   PetscMPIInt        rank, numProcs;
741:   PetscInt           d, ldepth, depth, p, pStart, pEnd, nroots, nleaves;
742:   PetscInt          *pointDepths, *remoteDepths, *ilocal;
743:   PetscInt          *depthRecv, *depthShift, *depthIdx;
744:   PetscInt           hybEnd[4];
745:   const PetscSFNode *iremote;
746:   PetscErrorCode     ierr;

750:   PetscObjectGetComm((PetscObject) dm, &comm);
751:   MPI_Comm_rank(comm, &rank);
752:   MPI_Comm_size(comm, &numProcs);
753:   DMPlexGetDepth(dm, &ldepth);
754:   MPIU_Allreduce(&ldepth, &depth, 1, MPIU_INT, MPI_MAX, comm);
755:   if ((ldepth >= 0) && (depth != ldepth)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent Plex depth %d != %d", ldepth, depth);

757:   /* Before building the migration SF we need to know the new stratum offsets */
758:   PetscSFGetGraph(sf, &nroots, &nleaves, NULL, &iremote);
759:   PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);
760:   DMPlexGetHybridBounds(dm,&hybEnd[depth],&hybEnd[depth-1],&hybEnd[1],&hybEnd[0]);
761:   for (d = 0; d < depth+1; ++d) {
762:     DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);
763:     for (p = pStart; p < pEnd; ++p) {
764:       if (hybEnd[d] >= 0 && p >= hybEnd[d]) { /* put in a separate value for hybrid points */
765:         pointDepths[p] = 2 * d;
766:       } else {
767:         pointDepths[p] = 2 * d + 1;
768:       }
769:     }
770:   }
771:   for (p = 0; p < nleaves; ++p) remoteDepths[p] = -1;
772:   PetscSFBcastBegin(sf, MPIU_INT, pointDepths, remoteDepths);
773:   PetscSFBcastEnd(sf, MPIU_INT, pointDepths, remoteDepths);
774:   /* Count recevied points in each stratum and compute the internal strata shift */
775:   PetscMalloc3(2*(depth+1), &depthRecv, 2*(depth+1), &depthShift, 2*(depth+1), &depthIdx);
776:   for (d = 0; d < 2*(depth+1); ++d) depthRecv[d] = 0;
777:   for (p = 0; p < nleaves; ++p) depthRecv[remoteDepths[p]]++;
778:   depthShift[2*depth+1] = 0;
779:   for (d = 0; d < 2*depth+1; ++d) depthShift[d] = depthRecv[2 * depth + 1];
780:   for (d = 0; d < 2*depth; ++d) depthShift[d] += depthRecv[2 * depth];
781:   depthShift[0] += depthRecv[1];
782:   for (d = 2; d < 2*depth; ++d) depthShift[d] += depthRecv[1];
783:   for (d = 2; d < 2*depth; ++d) depthShift[d] += depthRecv[0];
784:   for (d = 2 * depth-1; d > 2; --d) {
785:     PetscInt e;

787:     for (e = d -1; e > 1; --e) depthShift[e] += depthRecv[d];
788:   }
789:   for (d = 0; d < 2*(depth+1); ++d) {depthIdx[d] = 0;}
790:   /* Derive a new local permutation based on stratified indices */
791:   PetscMalloc1(nleaves, &ilocal);
792:   for (p = 0; p < nleaves; ++p) {
793:     const PetscInt dep = remoteDepths[p];

795:     ilocal[p] = depthShift[dep] + depthIdx[dep];
796:     depthIdx[dep]++;
797:   }
798:   PetscSFCreate(comm, migrationSF);
799:   PetscObjectSetName((PetscObject) *migrationSF, "Migration SF");
800:   PetscSFSetGraph(*migrationSF, nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_COPY_VALUES);
801:   PetscFree2(pointDepths,remoteDepths);
802:   PetscFree3(depthRecv, depthShift, depthIdx);
803:   return(0);
804: }

806: /*@
807:   DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution

809:   Collective on DM

811:   Input Parameters:
812: + dm - The DMPlex object
813: . pointSF - The PetscSF describing the communication pattern
814: . originalSection - The PetscSection for existing data layout
815: - originalVec - The existing data

817:   Output Parameters:
818: + newSection - The PetscSF describing the new data layout
819: - newVec - The new data

821:   Level: developer

823: .seealso: DMPlexDistribute(), DMPlexDistributeFieldIS(), DMPlexDistributeData()
824: @*/
825: PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec)
826: {
827:   PetscSF        fieldSF;
828:   PetscInt      *remoteOffsets, fieldSize;
829:   PetscScalar   *originalValues, *newValues;

833:   PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);
834:   PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);

836:   PetscSectionGetStorageSize(newSection, &fieldSize);
837:   VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);
838:   VecSetType(newVec,dm->vectype);

840:   VecGetArray(originalVec, &originalValues);
841:   VecGetArray(newVec, &newValues);
842:   PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);
843:   PetscFree(remoteOffsets);
844:   PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);
845:   PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);
846:   PetscSFDestroy(&fieldSF);
847:   VecRestoreArray(newVec, &newValues);
848:   VecRestoreArray(originalVec, &originalValues);
849:   PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);
850:   return(0);
851: }

853: /*@
854:   DMPlexDistributeFieldIS - Distribute field data to match a given PetscSF, usually the SF from mesh distribution

856:   Collective on DM

858:   Input Parameters:
859: + dm - The DMPlex object
860: . pointSF - The PetscSF describing the communication pattern
861: . originalSection - The PetscSection for existing data layout
862: - originalIS - The existing data

864:   Output Parameters:
865: + newSection - The PetscSF describing the new data layout
866: - newIS - The new data

868:   Level: developer

870: .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexDistributeData()
871: @*/
872: PetscErrorCode DMPlexDistributeFieldIS(DM dm, PetscSF pointSF, PetscSection originalSection, IS originalIS, PetscSection newSection, IS *newIS)
873: {
874:   PetscSF         fieldSF;
875:   PetscInt       *newValues, *remoteOffsets, fieldSize;
876:   const PetscInt *originalValues;
877:   PetscErrorCode  ierr;

880:   PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);
881:   PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);

883:   PetscSectionGetStorageSize(newSection, &fieldSize);
884:   PetscMalloc1(fieldSize, &newValues);

886:   ISGetIndices(originalIS, &originalValues);
887:   PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);
888:   PetscFree(remoteOffsets);
889:   PetscSFBcastBegin(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);
890:   PetscSFBcastEnd(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);
891:   PetscSFDestroy(&fieldSF);
892:   ISRestoreIndices(originalIS, &originalValues);
893:   ISCreateGeneral(PetscObjectComm((PetscObject) pointSF), fieldSize, newValues, PETSC_OWN_POINTER, newIS);
894:   PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);
895:   return(0);
896: }

898: /*@
899:   DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution

901:   Collective on DM

903:   Input Parameters:
904: + dm - The DMPlex object
905: . pointSF - The PetscSF describing the communication pattern
906: . originalSection - The PetscSection for existing data layout
907: . datatype - The type of data
908: - originalData - The existing data

910:   Output Parameters:
911: + newSection - The PetscSection describing the new data layout
912: - newData - The new data

914:   Level: developer

916: .seealso: DMPlexDistribute(), DMPlexDistributeField()
917: @*/
918: PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData)
919: {
920:   PetscSF        fieldSF;
921:   PetscInt      *remoteOffsets, fieldSize;
922:   PetscMPIInt    dataSize;

926:   PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);
927:   PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);

929:   PetscSectionGetStorageSize(newSection, &fieldSize);
930:   MPI_Type_size(datatype, &dataSize);
931:   PetscMalloc(fieldSize * dataSize, newData);

933:   PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);
934:   PetscFree(remoteOffsets);
935:   PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);
936:   PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);
937:   PetscSFDestroy(&fieldSF);
938:   PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);
939:   return(0);
940: }

942: PetscErrorCode DMPlexDistributeCones(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping original, ISLocalToGlobalMapping renumbering, DM dmParallel)
943: {
944:   DM_Plex               *mesh  = (DM_Plex*) dm->data;
945:   DM_Plex               *pmesh = (DM_Plex*) (dmParallel)->data;
946:   MPI_Comm               comm;
947:   PetscSF                coneSF;
948:   PetscSection           originalConeSection, newConeSection;
949:   PetscInt              *remoteOffsets, *cones, *globCones, *newCones, newConesSize;
950:   PetscBool              flg;
951:   PetscErrorCode         ierr;

956:   PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);

958:   /* Distribute cone section */
959:   PetscObjectGetComm((PetscObject)dm, &comm);
960:   DMPlexGetConeSection(dm, &originalConeSection);
961:   DMPlexGetConeSection(dmParallel, &newConeSection);
962:   PetscSFDistributeSection(migrationSF, originalConeSection, &remoteOffsets, newConeSection);
963:   DMSetUp(dmParallel);
964:   {
965:     PetscInt pStart, pEnd, p;

967:     PetscSectionGetChart(newConeSection, &pStart, &pEnd);
968:     for (p = pStart; p < pEnd; ++p) {
969:       PetscInt coneSize;
970:       PetscSectionGetDof(newConeSection, p, &coneSize);
971:       pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize);
972:     }
973:   }
974:   /* Communicate and renumber cones */
975:   PetscSFCreateSectionSF(migrationSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);
976:   PetscFree(remoteOffsets);
977:   DMPlexGetCones(dm, &cones);
978:   if (original) {
979:     PetscInt numCones;

981:     PetscSectionGetStorageSize(originalConeSection,&numCones); PetscMalloc1(numCones,&globCones);
982:     ISLocalToGlobalMappingApplyBlock(original, numCones, cones, globCones);
983:   }
984:   else {
985:     globCones = cones;
986:   }
987:   DMPlexGetCones(dmParallel, &newCones);
988:   PetscSFBcastBegin(coneSF, MPIU_INT, globCones, newCones);
989:   PetscSFBcastEnd(coneSF, MPIU_INT, globCones, newCones);
990:   if (original) {
991:     PetscFree(globCones);
992:   }
993:   PetscSectionGetStorageSize(newConeSection, &newConesSize);
994:   ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);
995: #if PETSC_USE_DEBUG
996:   {
997:     PetscInt  p;
998:     PetscBool valid = PETSC_TRUE;
999:     for (p = 0; p < newConesSize; ++p) {
1000:       if (newCones[p] < 0) {valid = PETSC_FALSE; PetscPrintf(PETSC_COMM_SELF, "Point %d not in overlap SF\n", p);}
1001:     }
1002:     if (!valid) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map");
1003:   }
1004: #endif
1005:   PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-cones_view", &flg);
1006:   if (flg) {
1007:     PetscPrintf(comm, "Serial Cone Section:\n");
1008:     PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);
1009:     PetscPrintf(comm, "Parallel Cone Section:\n");
1010:     PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);
1011:     PetscSFView(coneSF, NULL);
1012:   }
1013:   DMPlexGetConeOrientations(dm, &cones);
1014:   DMPlexGetConeOrientations(dmParallel, &newCones);
1015:   PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);
1016:   PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);
1017:   PetscSFDestroy(&coneSF);
1018:   PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);
1019:   /* Create supports and stratify sieve */
1020:   {
1021:     PetscInt pStart, pEnd;

1023:     PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);
1024:     PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);
1025:   }
1026:   DMPlexSymmetrize(dmParallel);
1027:   DMPlexStratify(dmParallel);
1028:   pmesh->useCone    = mesh->useCone;
1029:   pmesh->useClosure = mesh->useClosure;
1030:   return(0);
1031: }

1033: PetscErrorCode DMPlexDistributeCoordinates(DM dm, PetscSF migrationSF, DM dmParallel)
1034: {
1035:   MPI_Comm         comm;
1036:   PetscSection     originalCoordSection, newCoordSection;
1037:   Vec              originalCoordinates, newCoordinates;
1038:   PetscInt         bs;
1039:   const char      *name;
1040:   const PetscReal *maxCell, *L;
1041:   const DMBoundaryType *bd;
1042:   PetscErrorCode   ierr;


1048:   PetscObjectGetComm((PetscObject)dm, &comm);
1049:   DMGetCoordinateSection(dm, &originalCoordSection);
1050:   DMGetCoordinateSection(dmParallel, &newCoordSection);
1051:   DMGetCoordinatesLocal(dm, &originalCoordinates);
1052:   if (originalCoordinates) {
1053:     VecCreate(PETSC_COMM_SELF, &newCoordinates);
1054:     PetscObjectGetName((PetscObject) originalCoordinates, &name);
1055:     PetscObjectSetName((PetscObject) newCoordinates, name);

1057:     DMPlexDistributeField(dm, migrationSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);
1058:     DMSetCoordinatesLocal(dmParallel, newCoordinates);
1059:     VecGetBlockSize(originalCoordinates, &bs);
1060:     VecSetBlockSize(newCoordinates, bs);
1061:     VecDestroy(&newCoordinates);
1062:   }
1063:   DMGetPeriodicity(dm, &maxCell, &L, &bd);
1064:   if (L) {DMSetPeriodicity(dmParallel, maxCell, L, bd);}
1065:   return(0);
1066: }

1068: /* Here we are assuming that process 0 always has everything */
1069: PetscErrorCode DMPlexDistributeLabels(DM dm, PetscSF migrationSF, DM dmParallel)
1070: {
1071:   DM_Plex         *mesh = (DM_Plex*) dm->data;
1072:   MPI_Comm         comm;
1073:   DMLabel          depthLabel;
1074:   PetscMPIInt      rank;
1075:   PetscInt         depth, d, numLabels, numLocalLabels, l;
1076:   PetscBool        hasLabels = PETSC_FALSE, lsendDepth, sendDepth;
1077:   PetscObjectState depthState = -1;
1078:   PetscErrorCode   ierr;

1083:   PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);
1084:   PetscObjectGetComm((PetscObject)dm, &comm);
1085:   MPI_Comm_rank(comm, &rank);

1087:   /* If the user has changed the depth label, communicate it instead */
1088:   DMPlexGetDepth(dm, &depth);
1089:   DMPlexGetDepthLabel(dm, &depthLabel);
1090:   if (depthLabel) {DMLabelGetState(depthLabel, &depthState);}
1091:   lsendDepth = mesh->depthState != depthState ? PETSC_TRUE : PETSC_FALSE;
1092:   MPIU_Allreduce(&lsendDepth, &sendDepth, 1, MPIU_BOOL, MPI_LOR, comm);
1093:   if (sendDepth) {
1094:     DMRemoveLabel(dmParallel, "depth", &depthLabel);
1095:     DMLabelDestroy(&depthLabel);
1096:   }
1097:   /* Everyone must have either the same number of labels, or none */
1098:   DMGetNumLabels(dm, &numLocalLabels);
1099:   numLabels = numLocalLabels;
1100:   MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);
1101:   if (numLabels == numLocalLabels) hasLabels = PETSC_TRUE;
1102:   for (l = numLabels-1; l >= 0; --l) {
1103:     DMLabel     label = NULL, labelNew = NULL;
1104:     PetscBool   isdepth;

1106:     if (hasLabels) {
1107:       DMGetLabelByNum(dm, l, &label);
1108:       /* Skip "depth" because it is recreated */
1109:       PetscStrcmp(label->name, "depth", &isdepth);
1110:     }
1111:     MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);
1112:     if (isdepth && !sendDepth) continue;
1113:     DMLabelDistribute(label, migrationSF, &labelNew);
1114:     if (isdepth) {
1115:       /* Put in any missing strata which can occur if users are managing the depth label themselves */
1116:       PetscInt gdepth;

1118:       MPIU_Allreduce(&depth, &gdepth, 1, MPIU_INT, MPI_MAX, comm);
1119:       if ((depth >= 0) && (gdepth != depth)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent Plex depth %d != %d", depth, gdepth);
1120:       for (d = 0; d <= gdepth; ++d) {
1121:         PetscBool has;

1123:         DMLabelHasStratum(labelNew, d, &has);
1124:         if (!has) DMLabelAddStratum(labelNew, d);
1125:       }
1126:     }
1127:     DMAddLabel(dmParallel, labelNew);
1128:   }
1129:   PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);
1130:   return(0);
1131: }

1133: PetscErrorCode DMPlexDistributeSetupHybrid(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
1134: {
1135:   DM_Plex        *mesh  = (DM_Plex*) dm->data;
1136:   DM_Plex        *pmesh = (DM_Plex*) (dmParallel)->data;
1137:   PetscBool      *isHybrid, *isHybridParallel;
1138:   PetscInt        dim, depth, d;
1139:   PetscInt        pStart, pEnd, pStartP, pEndP;
1140:   PetscErrorCode  ierr;


1146:   DMGetDimension(dm, &dim);
1147:   DMPlexGetDepth(dm, &depth);
1148:   DMPlexGetChart(dm,&pStart,&pEnd);
1149:   DMPlexGetChart(dmParallel,&pStartP,&pEndP);
1150:   PetscCalloc2(pEnd-pStart,&isHybrid,pEndP-pStartP,&isHybridParallel);
1151:   for (d = 0; d <= depth; d++) {
1152:     PetscInt hybridMax = (depth == 1 && d == 1) ? mesh->hybridPointMax[dim] : mesh->hybridPointMax[d];

1154:     if (hybridMax >= 0) {
1155:       PetscInt sStart, sEnd, p;

1157:       DMPlexGetDepthStratum(dm,d,&sStart,&sEnd);
1158:       for (p = hybridMax; p < sEnd; p++) isHybrid[p-pStart] = PETSC_TRUE;
1159:     }
1160:   }
1161:   PetscSFBcastBegin(migrationSF,MPIU_BOOL,isHybrid,isHybridParallel);
1162:   PetscSFBcastEnd(migrationSF,MPIU_BOOL,isHybrid,isHybridParallel);
1163:   for (d = 0; d <= dim; d++) pmesh->hybridPointMax[d] = -1;
1164:   for (d = 0; d <= depth; d++) {
1165:     PetscInt sStart, sEnd, p, dd;

1167:     DMPlexGetDepthStratum(dmParallel,d,&sStart,&sEnd);
1168:     dd = (depth == 1 && d == 1) ? dim : d;
1169:     for (p = sStart; p < sEnd; p++) {
1170:       if (isHybridParallel[p-pStartP]) {
1171:         pmesh->hybridPointMax[dd] = p;
1172:         break;
1173:       }
1174:     }
1175:   }
1176:   PetscFree2(isHybrid,isHybridParallel);
1177:   return(0);
1178: }

1180: PetscErrorCode DMPlexDistributeSetupTree(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping original, ISLocalToGlobalMapping renumbering, DM dmParallel)
1181: {
1182:   DM_Plex        *mesh  = (DM_Plex*) dm->data;
1183:   DM_Plex        *pmesh = (DM_Plex*) (dmParallel)->data;
1184:   MPI_Comm        comm;
1185:   DM              refTree;
1186:   PetscSection    origParentSection, newParentSection;
1187:   PetscInt        *origParents, *origChildIDs;
1188:   PetscBool       flg;
1189:   PetscErrorCode  ierr;

1194:   PetscObjectGetComm((PetscObject)dm, &comm);

1196:   /* Set up tree */
1197:   DMPlexGetReferenceTree(dm,&refTree);
1198:   DMPlexSetReferenceTree(dmParallel,refTree);
1199:   DMPlexGetTree(dm,&origParentSection,&origParents,&origChildIDs,NULL,NULL);
1200:   if (origParentSection) {
1201:     PetscInt        pStart, pEnd;
1202:     PetscInt        *newParents, *newChildIDs, *globParents;
1203:     PetscInt        *remoteOffsetsParents, newParentSize;
1204:     PetscSF         parentSF;

1206:     DMPlexGetChart(dmParallel, &pStart, &pEnd);
1207:     PetscSectionCreate(PetscObjectComm((PetscObject)dmParallel),&newParentSection);
1208:     PetscSectionSetChart(newParentSection,pStart,pEnd);
1209:     PetscSFDistributeSection(migrationSF, origParentSection, &remoteOffsetsParents, newParentSection);
1210:     PetscSFCreateSectionSF(migrationSF, origParentSection, remoteOffsetsParents, newParentSection, &parentSF);
1211:     PetscFree(remoteOffsetsParents);
1212:     PetscSectionGetStorageSize(newParentSection,&newParentSize);
1213:     PetscMalloc2(newParentSize,&newParents,newParentSize,&newChildIDs);
1214:     if (original) {
1215:       PetscInt numParents;

1217:       PetscSectionGetStorageSize(origParentSection,&numParents);
1218:       PetscMalloc1(numParents,&globParents);
1219:       ISLocalToGlobalMappingApplyBlock(original, numParents, origParents, globParents);
1220:     }
1221:     else {
1222:       globParents = origParents;
1223:     }
1224:     PetscSFBcastBegin(parentSF, MPIU_INT, globParents, newParents);
1225:     PetscSFBcastEnd(parentSF, MPIU_INT, globParents, newParents);
1226:     if (original) {
1227:       PetscFree(globParents);
1228:     }
1229:     PetscSFBcastBegin(parentSF, MPIU_INT, origChildIDs, newChildIDs);
1230:     PetscSFBcastEnd(parentSF, MPIU_INT, origChildIDs, newChildIDs);
1231:     ISGlobalToLocalMappingApplyBlock(renumbering,IS_GTOLM_MASK, newParentSize, newParents, NULL, newParents);
1232: #if PETSC_USE_DEBUG
1233:     {
1234:       PetscInt  p;
1235:       PetscBool valid = PETSC_TRUE;
1236:       for (p = 0; p < newParentSize; ++p) {
1237:         if (newParents[p] < 0) {valid = PETSC_FALSE; PetscPrintf(PETSC_COMM_SELF, "Point %d not in overlap SF\n", p);}
1238:       }
1239:       if (!valid) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map");
1240:     }
1241: #endif
1242:     PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-parents_view", &flg);
1243:     if (flg) {
1244:       PetscPrintf(comm, "Serial Parent Section: \n");
1245:       PetscSectionView(origParentSection, PETSC_VIEWER_STDOUT_WORLD);
1246:       PetscPrintf(comm, "Parallel Parent Section: \n");
1247:       PetscSectionView(newParentSection, PETSC_VIEWER_STDOUT_WORLD);
1248:       PetscSFView(parentSF, NULL);
1249:     }
1250:     DMPlexSetTree(dmParallel,newParentSection,newParents,newChildIDs);
1251:     PetscSectionDestroy(&newParentSection);
1252:     PetscFree2(newParents,newChildIDs);
1253:     PetscSFDestroy(&parentSF);
1254:   }
1255:   pmesh->useAnchors = mesh->useAnchors;
1256:   return(0);
1257: }

1259: PetscErrorCode DMPlexDistributeSF(DM dm, PetscSF migrationSF, DM dmParallel)
1260: {
1261:   DM_Plex               *mesh  = (DM_Plex*) dm->data;
1262:   DM_Plex               *pmesh = (DM_Plex*) (dmParallel)->data;
1263:   PetscMPIInt            rank, numProcs;
1264:   MPI_Comm               comm;
1265:   PetscErrorCode         ierr;


1271:   /* Create point SF for parallel mesh */
1272:   PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);
1273:   PetscObjectGetComm((PetscObject)dm, &comm);
1274:   MPI_Comm_rank(comm, &rank);
1275:   MPI_Comm_size(comm, &numProcs);
1276:   {
1277:     const PetscInt *leaves;
1278:     PetscSFNode    *remotePoints, *rowners, *lowners;
1279:     PetscInt        numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints;
1280:     PetscInt        pStart, pEnd;

1282:     DMPlexGetChart(dmParallel, &pStart, &pEnd);
1283:     PetscSFGetGraph(migrationSF, &numRoots, &numLeaves, &leaves, NULL);
1284:     PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);
1285:     for (p=0; p<numRoots; p++) {
1286:       rowners[p].rank  = -1;
1287:       rowners[p].index = -1;
1288:     }
1289:     PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);
1290:     PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);
1291:     for (p = 0; p < numLeaves; ++p) {
1292:       if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */
1293:         lowners[p].rank  = rank;
1294:         lowners[p].index = leaves ? leaves[p] : p;
1295:       } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */
1296:         lowners[p].rank  = -2;
1297:         lowners[p].index = -2;
1298:       }
1299:     }
1300:     for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */
1301:       rowners[p].rank  = -3;
1302:       rowners[p].index = -3;
1303:     }
1304:     PetscSFReduceBegin(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);
1305:     PetscSFReduceEnd(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);
1306:     PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);
1307:     PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);
1308:     for (p = 0; p < numLeaves; ++p) {
1309:       if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed");
1310:       if (lowners[p].rank != rank) ++numGhostPoints;
1311:     }
1312:     PetscMalloc1(numGhostPoints, &ghostPoints);
1313:     PetscMalloc1(numGhostPoints, &remotePoints);
1314:     for (p = 0, gp = 0; p < numLeaves; ++p) {
1315:       if (lowners[p].rank != rank) {
1316:         ghostPoints[gp]        = leaves ? leaves[p] : p;
1317:         remotePoints[gp].rank  = lowners[p].rank;
1318:         remotePoints[gp].index = lowners[p].index;
1319:         ++gp;
1320:       }
1321:     }
1322:     PetscFree2(rowners,lowners);
1323:     PetscSFSetGraph((dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);
1324:     PetscSFSetFromOptions((dmParallel)->sf);
1325:   }
1326:   pmesh->useCone    = mesh->useCone;
1327:   pmesh->useClosure = mesh->useClosure;
1328:   PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);
1329:   return(0);
1330: }

1332: /*@C
1333:   DMPlexDerivePointSF - Build a point SF from an SF describing a point migration

1335:   Input Parameter:
1336: + dm          - The source DMPlex object
1337: . migrationSF - The star forest that describes the parallel point remapping
1338: . ownership   - Flag causing a vote to determine point ownership

1340:   Output Parameter:
1341: - pointSF     - The star forest describing the point overlap in the remapped DM

1343:   Level: developer

1345: .seealso: DMPlexDistribute(), DMPlexDistributeOverlap()
1346: @*/
1347: PetscErrorCode DMPlexCreatePointSF(DM dm, PetscSF migrationSF, PetscBool ownership, PetscSF *pointSF)
1348: {
1349:   PetscMPIInt        rank;
1350:   PetscInt           p, nroots, nleaves, idx, npointLeaves;
1351:   PetscInt          *pointLocal;
1352:   const PetscInt    *leaves;
1353:   const PetscSFNode *roots;
1354:   PetscSFNode       *rootNodes, *leafNodes, *pointRemote;
1355:   PetscErrorCode     ierr;

1359:   MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);

1361:   PetscSFGetGraph(migrationSF, &nroots, &nleaves, &leaves, &roots);
1362:   PetscMalloc2(nroots, &rootNodes, nleaves, &leafNodes);
1363:   if (ownership) {
1364:     /* Point ownership vote: Process with highest rank ownes shared points */
1365:     for (p = 0; p < nleaves; ++p) {
1366:       /* Either put in a bid or we know we own it */
1367:       leafNodes[p].rank  = rank;
1368:       leafNodes[p].index = p;
1369:     }
1370:     for (p = 0; p < nroots; p++) {
1371:       /* Root must not participate in the reduction, flag so that MAXLOC does not use */
1372:       rootNodes[p].rank  = -3;
1373:       rootNodes[p].index = -3;
1374:     }
1375:     PetscSFReduceBegin(migrationSF, MPIU_2INT, leafNodes, rootNodes, MPI_MAXLOC);
1376:     PetscSFReduceEnd(migrationSF, MPIU_2INT, leafNodes, rootNodes, MPI_MAXLOC);
1377:   } else {
1378:     for (p = 0; p < nroots; p++) {
1379:       rootNodes[p].index = -1;
1380:       rootNodes[p].rank = rank;
1381:     };
1382:     for (p = 0; p < nleaves; p++) {
1383:       /* Write new local id into old location */
1384:       if (roots[p].rank == rank) {
1385:         rootNodes[roots[p].index].index = leaves ? leaves[p] : p;
1386:       }
1387:     }
1388:   }
1389:   PetscSFBcastBegin(migrationSF, MPIU_2INT, rootNodes, leafNodes);
1390:   PetscSFBcastEnd(migrationSF, MPIU_2INT, rootNodes, leafNodes);

1392:   for (npointLeaves = 0, p = 0; p < nleaves; p++) {if (leafNodes[p].rank != rank) npointLeaves++;}
1393:   PetscMalloc1(npointLeaves, &pointLocal);
1394:   PetscMalloc1(npointLeaves, &pointRemote);
1395:   for (idx = 0, p = 0; p < nleaves; p++) {
1396:     if (leafNodes[p].rank != rank) {
1397:       pointLocal[idx] = p;
1398:       pointRemote[idx] = leafNodes[p];
1399:       idx++;
1400:     }
1401:   }
1402:   PetscSFCreate(PetscObjectComm((PetscObject) dm), pointSF);
1403:   PetscSFSetFromOptions(*pointSF);
1404:   PetscSFSetGraph(*pointSF, nleaves, npointLeaves, pointLocal, PETSC_OWN_POINTER, pointRemote, PETSC_OWN_POINTER);
1405:   PetscFree2(rootNodes, leafNodes);
1406:   return(0);
1407: }

1409: /*@C
1410:   DMPlexMigrate  - Migrates internal DM data over the supplied star forest

1412:   Input Parameter:
1413: + dm       - The source DMPlex object
1414: . sf       - The star forest communication context describing the migration pattern

1416:   Output Parameter:
1417: - targetDM - The target DMPlex object

1419:   Level: intermediate

1421: .seealso: DMPlexDistribute(), DMPlexDistributeOverlap()
1422: @*/
1423: PetscErrorCode DMPlexMigrate(DM dm, PetscSF sf, DM targetDM)
1424: {
1425:   MPI_Comm               comm;
1426:   PetscInt               dim, nroots;
1427:   PetscSF                sfPoint;
1428:   ISLocalToGlobalMapping ltogMigration;
1429:   ISLocalToGlobalMapping ltogOriginal = NULL;
1430:   PetscBool              flg;
1431:   PetscErrorCode         ierr;

1435:   PetscLogEventBegin(DMPLEX_Migrate, dm, 0, 0, 0);
1436:   PetscObjectGetComm((PetscObject) dm, &comm);
1437:   DMGetDimension(dm, &dim);
1438:   DMSetDimension(targetDM, dim);

1440:   /* Check for a one-to-all distribution pattern */
1441:   DMGetPointSF(dm, &sfPoint);
1442:   PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL);
1443:   if (nroots >= 0) {
1444:     IS                     isOriginal;
1445:     PetscInt               n, size, nleaves;
1446:     PetscInt              *numbering_orig, *numbering_new;
1447:     /* Get the original point numbering */
1448:     DMPlexCreatePointNumbering(dm, &isOriginal);
1449:     ISLocalToGlobalMappingCreateIS(isOriginal, &ltogOriginal);
1450:     ISLocalToGlobalMappingGetSize(ltogOriginal, &size);
1451:     ISLocalToGlobalMappingGetBlockIndices(ltogOriginal, (const PetscInt**)&numbering_orig);
1452:     /* Convert to positive global numbers */
1453:     for (n=0; n<size; n++) {if (numbering_orig[n] < 0) numbering_orig[n] = -(numbering_orig[n]+1);}
1454:     /* Derive the new local-to-global mapping from the old one */
1455:     PetscSFGetGraph(sf, NULL, &nleaves, NULL, NULL);
1456:     PetscMalloc1(nleaves, &numbering_new);
1457:     PetscSFBcastBegin(sf, MPIU_INT, (PetscInt *) numbering_orig, numbering_new);
1458:     PetscSFBcastEnd(sf, MPIU_INT, (PetscInt *) numbering_orig, numbering_new);
1459:     ISLocalToGlobalMappingCreate(comm, 1, nleaves, (const PetscInt*) numbering_new, PETSC_OWN_POINTER, &ltogMigration);
1460:     ISLocalToGlobalMappingRestoreIndices(ltogOriginal, (const PetscInt**)&numbering_orig);
1461:     ISDestroy(&isOriginal);
1462:   } else {
1463:     /* One-to-all distribution pattern: We can derive LToG from SF */
1464:     ISLocalToGlobalMappingCreateSF(sf, 0, &ltogMigration);
1465:   }
1466:   PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-partition_view", &flg);
1467:   if (flg) {
1468:     PetscPrintf(comm, "Point renumbering for DM migration:\n");
1469:     ISLocalToGlobalMappingView(ltogMigration, NULL);
1470:   }
1471:   /* Migrate DM data to target DM */
1472:   DMPlexDistributeCones(dm, sf, ltogOriginal, ltogMigration, targetDM);
1473:   DMPlexDistributeLabels(dm, sf, targetDM);
1474:   DMPlexDistributeCoordinates(dm, sf, targetDM);
1475:   DMPlexDistributeSetupHybrid(dm, sf, ltogMigration, targetDM);
1476:   DMPlexDistributeSetupTree(dm, sf, ltogOriginal, ltogMigration, targetDM);
1477:   ISLocalToGlobalMappingDestroy(&ltogOriginal);
1478:   ISLocalToGlobalMappingDestroy(&ltogMigration);
1479:   PetscLogEventEnd(DMPLEX_Migrate, dm, 0, 0, 0);
1480:   return(0);
1481: }

1483: /*@C
1484:   DMPlexDistribute - Distributes the mesh and any associated sections.

1486:   Not Collective

1488:   Input Parameter:
1489: + dm  - The original DMPlex object
1490: - overlap - The overlap of partitions, 0 is the default

1492:   Output Parameter:
1493: + sf - The PetscSF used for point distribution
1494: - parallelMesh - The distributed DMPlex object, or NULL

1496:   Note: If the mesh was not distributed, the return value is NULL.

1498:   The user can control the definition of adjacency for the mesh using DMPlexSetAdjacencyUseCone() and
1499:   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
1500:   representation on the mesh.

1502:   Level: intermediate

1504: .keywords: mesh, elements
1505: .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
1506: @*/
1507: PetscErrorCode DMPlexDistribute(DM dm, PetscInt overlap, PetscSF *sf, DM *dmParallel)
1508: {
1509:   MPI_Comm               comm;
1510:   PetscPartitioner       partitioner;
1511:   IS                     cellPart;
1512:   PetscSection           cellPartSection;
1513:   DM                     dmCoord;
1514:   DMLabel                lblPartition, lblMigration;
1515:   PetscSF                sfProcess, sfMigration, sfStratified, sfPoint;
1516:   PetscBool              flg;
1517:   PetscMPIInt            rank, numProcs, p;
1518:   PetscErrorCode         ierr;


1525:   PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);
1526:   PetscObjectGetComm((PetscObject)dm,&comm);
1527:   MPI_Comm_rank(comm, &rank);
1528:   MPI_Comm_size(comm, &numProcs);

1530:   *dmParallel = NULL;
1531:   if (numProcs == 1) {
1532:     PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);
1533:     return(0);
1534:   }

1536:   /* Create cell partition */
1537:   PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);
1538:   PetscSectionCreate(comm, &cellPartSection);
1539:   DMPlexGetPartitioner(dm, &partitioner);
1540:   PetscPartitionerPartition(partitioner, dm, cellPartSection, &cellPart);
1541:   {
1542:     /* Convert partition to DMLabel */
1543:     PetscInt proc, pStart, pEnd, npoints, poffset;
1544:     const PetscInt *points;
1545:     DMLabelCreate("Point Partition", &lblPartition);
1546:     ISGetIndices(cellPart, &points);
1547:     PetscSectionGetChart(cellPartSection, &pStart, &pEnd);
1548:     for (proc = pStart; proc < pEnd; proc++) {
1549:       PetscSectionGetDof(cellPartSection, proc, &npoints);
1550:       PetscSectionGetOffset(cellPartSection, proc, &poffset);
1551:       for (p = poffset; p < poffset+npoints; p++) {
1552:         DMLabelSetValue(lblPartition, points[p], proc);
1553:       }
1554:     }
1555:     ISRestoreIndices(cellPart, &points);
1556:   }
1557:   DMPlexPartitionLabelClosure(dm, lblPartition);
1558:   {
1559:     /* Build a global process SF */
1560:     PetscSFNode *remoteProc;
1561:     PetscMalloc1(numProcs, &remoteProc);
1562:     for (p = 0; p < numProcs; ++p) {
1563:       remoteProc[p].rank  = p;
1564:       remoteProc[p].index = rank;
1565:     }
1566:     PetscSFCreate(comm, &sfProcess);
1567:     PetscObjectSetName((PetscObject) sfProcess, "Process SF");
1568:     PetscSFSetGraph(sfProcess, numProcs, numProcs, NULL, PETSC_OWN_POINTER, remoteProc, PETSC_OWN_POINTER);
1569:   }
1570:   DMLabelCreate("Point migration", &lblMigration);
1571:   DMPlexPartitionLabelInvert(dm, lblPartition, sfProcess, lblMigration);
1572:   DMPlexPartitionLabelCreateSF(dm, lblMigration, &sfMigration);
1573:   /* Stratify the SF in case we are migrating an already parallel plex */
1574:   DMPlexStratifyMigrationSF(dm, sfMigration, &sfStratified);
1575:   PetscSFDestroy(&sfMigration);
1576:   sfMigration = sfStratified;
1577:   PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);
1578:   PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-partition_view", &flg);
1579:   if (flg) {
1580:     DMLabelView(lblPartition, PETSC_VIEWER_STDOUT_WORLD);
1581:     PetscSFView(sfMigration, PETSC_VIEWER_STDOUT_WORLD);
1582:   }

1584:   /* Create non-overlapping parallel DM and migrate internal data */
1585:   DMPlexCreate(comm, dmParallel);
1586:   PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");
1587:   DMPlexMigrate(dm, sfMigration, *dmParallel);

1589:   /* Build the point SF without overlap */
1590:   DMPlexCreatePointSF(*dmParallel, sfMigration, PETSC_TRUE, &sfPoint);
1591:   DMSetPointSF(*dmParallel, sfPoint);
1592:   DMGetCoordinateDM(*dmParallel, &dmCoord);
1593:   if (dmCoord) {DMSetPointSF(dmCoord, sfPoint);}
1594:   if (flg) {PetscSFView(sfPoint, PETSC_VIEWER_STDOUT_WORLD);}

1596:   if (overlap > 0) {
1597:     DM                 dmOverlap;
1598:     PetscInt           nroots, nleaves;
1599:     PetscSFNode       *newRemote;
1600:     const PetscSFNode *oldRemote;
1601:     PetscSF            sfOverlap, sfOverlapPoint;
1602:     /* Add the partition overlap to the distributed DM */
1603:     DMPlexDistributeOverlap(*dmParallel, overlap, &sfOverlap, &dmOverlap);
1604:     DMDestroy(dmParallel);
1605:     *dmParallel = dmOverlap;
1606:     if (flg) {
1607:       PetscPrintf(comm, "Overlap Migration SF:\n");
1608:       PetscSFView(sfOverlap, NULL);
1609:     }

1611:     /* Re-map the migration SF to establish the full migration pattern */
1612:     PetscSFGetGraph(sfMigration, &nroots, NULL, NULL, &oldRemote);
1613:     PetscSFGetGraph(sfOverlap, NULL, &nleaves, NULL, NULL);
1614:     PetscMalloc1(nleaves, &newRemote);
1615:     PetscSFBcastBegin(sfOverlap, MPIU_2INT, oldRemote, newRemote);
1616:     PetscSFBcastEnd(sfOverlap, MPIU_2INT, oldRemote, newRemote);
1617:     PetscSFCreate(comm, &sfOverlapPoint);
1618:     PetscSFSetGraph(sfOverlapPoint, nroots, nleaves, NULL, PETSC_OWN_POINTER, newRemote, PETSC_OWN_POINTER);
1619:     PetscSFDestroy(&sfOverlap);
1620:     PetscSFDestroy(&sfMigration);
1621:     sfMigration = sfOverlapPoint;
1622:   }
1623:   /* Cleanup Partition */
1624:   PetscSFDestroy(&sfProcess);
1625:   DMLabelDestroy(&lblPartition);
1626:   DMLabelDestroy(&lblMigration);
1627:   PetscSectionDestroy(&cellPartSection);
1628:   ISDestroy(&cellPart);
1629:   /* Copy BC */
1630:   DMCopyBoundary(dm, *dmParallel);
1631:   /* Create sfNatural */
1632:   if (dm->useNatural) {
1633:     PetscSection section;

1635:     DMGetDefaultSection(dm, &section);
1636:     DMPlexCreateGlobalToNaturalSF(*dmParallel, section, sfMigration, &(*dmParallel)->sfNatural);
1637:   }
1638:   /* Cleanup */
1639:   if (sf) {*sf = sfMigration;}
1640:   else    {PetscSFDestroy(&sfMigration);}
1641:   PetscSFDestroy(&sfPoint);
1642:   PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);
1643:   return(0);
1644: }

1646: /*@C
1647:   DMPlexDistribute - Add partition overlap to a distributed non-overlapping DM.

1649:   Not Collective

1651:   Input Parameter:
1652: + dm  - The non-overlapping distrbuted DMPlex object
1653: - overlap - The overlap of partitions, 0 is the default

1655:   Output Parameter:
1656: + sf - The PetscSF used for point distribution
1657: - dmOverlap - The overlapping distributed DMPlex object, or NULL

1659:   Note: If the mesh was not distributed, the return value is NULL.

1661:   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
1662:   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
1663:   representation on the mesh.

1665:   Level: intermediate

1667: .keywords: mesh, elements
1668: .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
1669: @*/
1670: PetscErrorCode DMPlexDistributeOverlap(DM dm, PetscInt overlap, PetscSF *sf, DM *dmOverlap)
1671: {
1672:   MPI_Comm               comm;
1673:   PetscMPIInt            rank;
1674:   PetscSection           rootSection, leafSection;
1675:   IS                     rootrank, leafrank;
1676:   DM                     dmCoord;
1677:   DMLabel                lblOverlap;
1678:   PetscSF                sfOverlap, sfStratified, sfPoint;
1679:   PetscErrorCode         ierr;


1686:   PetscLogEventBegin(DMPLEX_DistributeOverlap, dm, 0, 0, 0);
1687:   PetscObjectGetComm((PetscObject)dm,&comm);
1688:   MPI_Comm_rank(comm, &rank);

1690:   /* Compute point overlap with neighbouring processes on the distributed DM */
1691:   PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);
1692:   PetscSectionCreate(comm, &rootSection);
1693:   PetscSectionCreate(comm, &leafSection);
1694:   DMPlexDistributeOwnership(dm, rootSection, &rootrank, leafSection, &leafrank);
1695:   DMPlexCreateOverlap(dm, overlap, rootSection, rootrank, leafSection, leafrank, &lblOverlap);
1696:   /* Convert overlap label to stratified migration SF */
1697:   DMPlexPartitionLabelCreateSF(dm, lblOverlap, &sfOverlap);
1698:   DMPlexStratifyMigrationSF(dm, sfOverlap, &sfStratified);
1699:   PetscSFDestroy(&sfOverlap);
1700:   sfOverlap = sfStratified;
1701:   PetscObjectSetName((PetscObject) sfOverlap, "Overlap SF");
1702:   PetscSFSetFromOptions(sfOverlap);

1704:   PetscSectionDestroy(&rootSection);
1705:   PetscSectionDestroy(&leafSection);
1706:   ISDestroy(&rootrank);
1707:   ISDestroy(&leafrank);
1708:   PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);

1710:   /* Build the overlapping DM */
1711:   DMPlexCreate(comm, dmOverlap);
1712:   PetscObjectSetName((PetscObject) *dmOverlap, "Parallel Mesh");
1713:   DMPlexMigrate(dm, sfOverlap, *dmOverlap);
1714:   /* Build the new point SF */
1715:   DMPlexCreatePointSF(*dmOverlap, sfOverlap, PETSC_FALSE, &sfPoint);
1716:   DMSetPointSF(*dmOverlap, sfPoint);
1717:   DMGetCoordinateDM(*dmOverlap, &dmCoord);
1718:   if (dmCoord) {DMSetPointSF(dmCoord, sfPoint);}
1719:   PetscSFDestroy(&sfPoint);
1720:   /* Cleanup overlap partition */
1721:   DMLabelDestroy(&lblOverlap);
1722:   if (sf) *sf = sfOverlap;
1723:   else    {PetscSFDestroy(&sfOverlap);}
1724:   PetscLogEventEnd(DMPLEX_DistributeOverlap, dm, 0, 0, 0);
1725:   return(0);
1726: }

1728: /*@C
1729:   DMPlexGetGatherDM - Get a copy of the DMPlex that gathers all points on the
1730:   root process of the original's communicator.

1732:   Input Parameters:
1733: . dm - the original DMPlex object

1735:   Output Parameters:
1736: . gatherMesh - the gathered DM object, or NULL

1738:   Level: intermediate

1740: .keywords: mesh
1741: .seealso: DMPlexDistribute(), DMPlexGetRedundantDM()
1742: @*/
1743: PetscErrorCode DMPlexGetGatherDM(DM dm, DM * gatherMesh)
1744: {
1745:   MPI_Comm       comm;
1746:   PetscMPIInt    size;
1747:   PetscPartitioner oldPart, gatherPart;

1752:   comm = PetscObjectComm((PetscObject)dm);
1753:   MPI_Comm_size(comm,&size);
1754:   *gatherMesh = NULL;
1755:   if (size == 1) return(0);
1756:   DMPlexGetPartitioner(dm,&oldPart);
1757:   PetscObjectReference((PetscObject)oldPart);
1758:   PetscPartitionerCreate(comm,&gatherPart);
1759:   PetscPartitionerSetType(gatherPart,PETSCPARTITIONERGATHER);
1760:   DMPlexSetPartitioner(dm,gatherPart);
1761:   DMPlexDistribute(dm,0,NULL,gatherMesh);
1762:   DMPlexSetPartitioner(dm,oldPart);
1763:   PetscPartitionerDestroy(&gatherPart);
1764:   PetscPartitionerDestroy(&oldPart);
1765:   return(0);
1766: }

1768: /*@C
1769:   DMPlexGetRedundantDM - Get a copy of the DMPlex that is completely copied on each process.

1771:   Input Parameters:
1772: . dm - the original DMPlex object

1774:   Output Parameters:
1775: . redundantMesh - the redundant DM object, or NULL

1777:   Level: intermediate

1779: .keywords: mesh
1780: .seealso: DMPlexDistribute(), DMPlexGetGatherDM()
1781: @*/
1782: PetscErrorCode DMPlexGetRedundantDM(DM dm, DM * redundantMesh)
1783: {
1784:   MPI_Comm       comm;
1785:   PetscMPIInt    size, rank;
1786:   PetscInt       pStart, pEnd, p;
1787:   PetscInt       numPoints = -1;
1788:   PetscSF        migrationSF, sfPoint;
1789:   DM             gatherDM, dmCoord;
1790:   PetscSFNode    *points;

1795:   comm = PetscObjectComm((PetscObject)dm);
1796:   MPI_Comm_size(comm,&size);
1797:   *redundantMesh = NULL;
1798:   if (size == 1) {
1799:     PetscObjectReference((PetscObject) dm);
1800:     *redundantMesh = dm;
1801:     return(0);
1802:   }
1803:   DMPlexGetGatherDM(dm,&gatherDM);
1804:   if (!gatherDM) return(0);
1805:   MPI_Comm_rank(comm,&rank);
1806:   DMPlexGetChart(gatherDM,&pStart,&pEnd);
1807:   numPoints = pEnd - pStart;
1808:   MPI_Bcast(&numPoints,1,MPIU_INT,0,comm);
1809:   PetscMalloc1(numPoints,&points);
1810:   PetscSFCreate(comm,&migrationSF);
1811:   for (p = 0; p < numPoints; p++) {
1812:     points[p].index = p;
1813:     points[p].rank  = 0;
1814:   }
1815:   PetscSFSetGraph(migrationSF,pEnd-pStart,numPoints,NULL,PETSC_OWN_POINTER,points,PETSC_OWN_POINTER);
1816:   DMPlexCreate(comm, redundantMesh);
1817:   PetscObjectSetName((PetscObject) *redundantMesh, "Redundant Mesh");
1818:   DMPlexMigrate(gatherDM, migrationSF, *redundantMesh);
1819:   DMPlexCreatePointSF(*redundantMesh, migrationSF, PETSC_FALSE, &sfPoint);
1820:   DMSetPointSF(*redundantMesh, sfPoint);
1821:   DMGetCoordinateDM(*redundantMesh, &dmCoord);
1822:   if (dmCoord) {DMSetPointSF(dmCoord, sfPoint);}
1823:   PetscSFDestroy(&sfPoint);
1824:   PetscSFDestroy(&migrationSF);
1825:   DMDestroy(&gatherDM);
1826:   return(0);
1827: }