xref: /petsc/src/dm/impls/plex/plexdistribute.c (revision 5de52c6d2b8d6de382140bd9fae367e1f02f2f13)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>    /*I      "petscdmplex.h"   I*/
2c58f1c22SToby Isaac #include <petsc/private/dmlabelimpl.h>   /*I      "petscdmlabel.h"  I*/
370034214SMatthew G. Knepley 
43c1f0c11SLawrence Mitchell /*@C
53c1f0c11SLawrence Mitchell   DMPlexSetAdjacencyUser - Define adjacency in the mesh using a user-provided callback
63c1f0c11SLawrence Mitchell 
73c1f0c11SLawrence Mitchell   Input Parameters:
83c1f0c11SLawrence Mitchell + dm      - The DM object
93c1f0c11SLawrence Mitchell . user    - The user callback, may be NULL (to clear the callback)
103c1f0c11SLawrence Mitchell - ctx     - context for callback evaluation, may be NULL
113c1f0c11SLawrence Mitchell 
123c1f0c11SLawrence Mitchell   Level: advanced
133c1f0c11SLawrence Mitchell 
143c1f0c11SLawrence Mitchell   Notes:
153c1f0c11SLawrence Mitchell      The caller of DMPlexGetAdjacency may need to arrange that a large enough array is available for the adjacency.
163c1f0c11SLawrence Mitchell 
173c1f0c11SLawrence Mitchell      Any setting here overrides other configuration of DMPlex adjacency determination.
183c1f0c11SLawrence Mitchell 
19db781477SPatrick Sanan .seealso: `DMSetAdjacency()`, `DMPlexDistribute()`, `DMPlexPreallocateOperator()`, `DMPlexGetAdjacency()`, `DMPlexGetAdjacencyUser()`
203c1f0c11SLawrence Mitchell @*/
213c1f0c11SLawrence Mitchell PetscErrorCode DMPlexSetAdjacencyUser(DM dm,PetscErrorCode (*user)(DM,PetscInt,PetscInt*,PetscInt[],void*),void *ctx)
223c1f0c11SLawrence Mitchell {
233c1f0c11SLawrence Mitchell   DM_Plex *mesh = (DM_Plex *)dm->data;
243c1f0c11SLawrence Mitchell 
253c1f0c11SLawrence Mitchell   PetscFunctionBegin;
263c1f0c11SLawrence Mitchell   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
273c1f0c11SLawrence Mitchell   mesh->useradjacency = user;
283c1f0c11SLawrence Mitchell   mesh->useradjacencyctx = ctx;
293c1f0c11SLawrence Mitchell   PetscFunctionReturn(0);
303c1f0c11SLawrence Mitchell }
313c1f0c11SLawrence Mitchell 
323c1f0c11SLawrence Mitchell /*@C
333c1f0c11SLawrence Mitchell   DMPlexGetAdjacencyUser - get the user-defined adjacency callback
343c1f0c11SLawrence Mitchell 
353c1f0c11SLawrence Mitchell   Input Parameter:
363c1f0c11SLawrence Mitchell . dm      - The DM object
373c1f0c11SLawrence Mitchell 
383c1f0c11SLawrence Mitchell   Output Parameters:
396b867d5aSJose E. Roman + user    - The user callback
403c1f0c11SLawrence Mitchell - ctx     - context for callback evaluation
413c1f0c11SLawrence Mitchell 
423c1f0c11SLawrence Mitchell   Level: advanced
433c1f0c11SLawrence Mitchell 
44db781477SPatrick Sanan .seealso: `DMSetAdjacency()`, `DMPlexDistribute()`, `DMPlexPreallocateOperator()`, `DMPlexGetAdjacency()`, `DMPlexSetAdjacencyUser()`
453c1f0c11SLawrence Mitchell @*/
463c1f0c11SLawrence Mitchell PetscErrorCode DMPlexGetAdjacencyUser(DM dm, PetscErrorCode (**user)(DM,PetscInt,PetscInt*,PetscInt[],void*), void **ctx)
473c1f0c11SLawrence Mitchell {
483c1f0c11SLawrence Mitchell   DM_Plex *mesh = (DM_Plex *)dm->data;
493c1f0c11SLawrence Mitchell 
503c1f0c11SLawrence Mitchell   PetscFunctionBegin;
513c1f0c11SLawrence Mitchell   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
523c1f0c11SLawrence Mitchell   if (user) *user = mesh->useradjacency;
533c1f0c11SLawrence Mitchell   if (ctx) *ctx = mesh->useradjacencyctx;
543c1f0c11SLawrence Mitchell   PetscFunctionReturn(0);
553c1f0c11SLawrence Mitchell }
563c1f0c11SLawrence Mitchell 
5770034214SMatthew G. Knepley /*@
58a17985deSToby Isaac   DMPlexSetAdjacencyUseAnchors - Define adjacency in the mesh using the point-to-point constraints.
598b0b4c70SToby Isaac 
608b0b4c70SToby Isaac   Input Parameters:
618b0b4c70SToby Isaac + dm      - The DM object
625b317d89SToby Isaac - useAnchors - Flag to use the constraints.  If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place.
638b0b4c70SToby Isaac 
648b0b4c70SToby Isaac   Level: intermediate
658b0b4c70SToby Isaac 
66db781477SPatrick Sanan .seealso: `DMGetAdjacency()`, `DMSetAdjacency()`, `DMPlexDistribute()`, `DMPlexPreallocateOperator()`, `DMPlexSetAnchors()`
678b0b4c70SToby Isaac @*/
685b317d89SToby Isaac PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM dm, PetscBool useAnchors)
698b0b4c70SToby Isaac {
708b0b4c70SToby Isaac   DM_Plex *mesh = (DM_Plex *) dm->data;
718b0b4c70SToby Isaac 
728b0b4c70SToby Isaac   PetscFunctionBegin;
738b0b4c70SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
745b317d89SToby Isaac   mesh->useAnchors = useAnchors;
758b0b4c70SToby Isaac   PetscFunctionReturn(0);
768b0b4c70SToby Isaac }
778b0b4c70SToby Isaac 
788b0b4c70SToby Isaac /*@
79a17985deSToby Isaac   DMPlexGetAdjacencyUseAnchors - Query whether adjacency in the mesh uses the point-to-point constraints.
808b0b4c70SToby Isaac 
818b0b4c70SToby Isaac   Input Parameter:
828b0b4c70SToby Isaac . dm      - The DM object
838b0b4c70SToby Isaac 
848b0b4c70SToby Isaac   Output Parameter:
855b317d89SToby Isaac . useAnchors - Flag to use the closure.  If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place.
868b0b4c70SToby Isaac 
878b0b4c70SToby Isaac   Level: intermediate
888b0b4c70SToby Isaac 
89db781477SPatrick Sanan .seealso: `DMPlexSetAdjacencyUseAnchors()`, `DMSetAdjacency()`, `DMGetAdjacency()`, `DMPlexDistribute()`, `DMPlexPreallocateOperator()`, `DMPlexSetAnchors()`
908b0b4c70SToby Isaac @*/
915b317d89SToby Isaac PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM dm, PetscBool *useAnchors)
928b0b4c70SToby Isaac {
938b0b4c70SToby Isaac   DM_Plex *mesh = (DM_Plex *) dm->data;
948b0b4c70SToby Isaac 
958b0b4c70SToby Isaac   PetscFunctionBegin;
968b0b4c70SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
97064a246eSJacob Faibussowitsch   PetscValidBoolPointer(useAnchors, 2);
985b317d89SToby Isaac   *useAnchors = mesh->useAnchors;
998b0b4c70SToby Isaac   PetscFunctionReturn(0);
1008b0b4c70SToby Isaac }
1018b0b4c70SToby Isaac 
10270034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Cone_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
10370034214SMatthew G. Knepley {
10470034214SMatthew G. Knepley   const PetscInt *cone = NULL;
10570034214SMatthew G. Knepley   PetscInt        numAdj = 0, maxAdjSize = *adjSize, coneSize, c;
10670034214SMatthew G. Knepley 
10770034214SMatthew G. Knepley   PetscFunctionBeginHot;
1089566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSize(dm, p, &coneSize));
1099566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCone(dm, p, &cone));
1104b6b44bdSMatthew G. Knepley   for (c = 0; c <= coneSize; ++c) {
1114b6b44bdSMatthew G. Knepley     const PetscInt  point   = !c ? p : cone[c-1];
11270034214SMatthew G. Knepley     const PetscInt *support = NULL;
11370034214SMatthew G. Knepley     PetscInt        supportSize, s, q;
11470034214SMatthew G. Knepley 
1159566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupportSize(dm, point, &supportSize));
1169566063dSJacob Faibussowitsch     PetscCall(DMPlexGetSupport(dm, point, &support));
11770034214SMatthew G. Knepley     for (s = 0; s < supportSize; ++s) {
118527e7439SSatish Balay       for (q = 0; q < numAdj || ((void)(adj[numAdj++] = support[s]),0); ++q) {
11970034214SMatthew G. Knepley         if (support[s] == adj[q]) break;
12070034214SMatthew G. Knepley       }
12163a3b9bcSJacob Faibussowitsch       PetscCheck(numAdj <= maxAdjSize,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%" PetscInt_FMT ")", maxAdjSize);
12270034214SMatthew G. Knepley     }
12370034214SMatthew G. Knepley   }
12470034214SMatthew G. Knepley   *adjSize = numAdj;
12570034214SMatthew G. Knepley   PetscFunctionReturn(0);
12670034214SMatthew G. Knepley }
12770034214SMatthew G. Knepley 
12870034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Support_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
12970034214SMatthew G. Knepley {
13070034214SMatthew G. Knepley   const PetscInt *support = NULL;
13170034214SMatthew G. Knepley   PetscInt        numAdj   = 0, maxAdjSize = *adjSize, supportSize, s;
13270034214SMatthew G. Knepley 
13370034214SMatthew G. Knepley   PetscFunctionBeginHot;
1349566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSupportSize(dm, p, &supportSize));
1359566063dSJacob Faibussowitsch   PetscCall(DMPlexGetSupport(dm, p, &support));
1364b6b44bdSMatthew G. Knepley   for (s = 0; s <= supportSize; ++s) {
1374b6b44bdSMatthew G. Knepley     const PetscInt  point = !s ? p : support[s-1];
13870034214SMatthew G. Knepley     const PetscInt *cone  = NULL;
13970034214SMatthew G. Knepley     PetscInt        coneSize, c, q;
14070034214SMatthew G. Knepley 
1419566063dSJacob Faibussowitsch     PetscCall(DMPlexGetConeSize(dm, point, &coneSize));
1429566063dSJacob Faibussowitsch     PetscCall(DMPlexGetCone(dm, point, &cone));
14370034214SMatthew G. Knepley     for (c = 0; c < coneSize; ++c) {
144527e7439SSatish Balay       for (q = 0; q < numAdj || ((void)(adj[numAdj++] = cone[c]),0); ++q) {
14570034214SMatthew G. Knepley         if (cone[c] == adj[q]) break;
14670034214SMatthew G. Knepley       }
14763a3b9bcSJacob Faibussowitsch       PetscCheck(numAdj <= maxAdjSize,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%" PetscInt_FMT ")", maxAdjSize);
14870034214SMatthew G. Knepley     }
14970034214SMatthew G. Knepley   }
15070034214SMatthew G. Knepley   *adjSize = numAdj;
15170034214SMatthew G. Knepley   PetscFunctionReturn(0);
15270034214SMatthew G. Knepley }
15370034214SMatthew G. Knepley 
15470034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Transitive_Internal(DM dm, PetscInt p, PetscBool useClosure, PetscInt *adjSize, PetscInt adj[])
15570034214SMatthew G. Knepley {
15670034214SMatthew G. Knepley   PetscInt      *star = NULL;
15770034214SMatthew G. Knepley   PetscInt       numAdj = 0, maxAdjSize = *adjSize, starSize, s;
15870034214SMatthew G. Knepley 
15970034214SMatthew G. Knepley   PetscFunctionBeginHot;
1609566063dSJacob Faibussowitsch   PetscCall(DMPlexGetTransitiveClosure(dm, p, useClosure, &starSize, &star));
16170034214SMatthew G. Knepley   for (s = 0; s < starSize*2; s += 2) {
16270034214SMatthew G. Knepley     const PetscInt *closure = NULL;
16370034214SMatthew G. Knepley     PetscInt        closureSize, c, q;
16470034214SMatthew G. Knepley 
1659566063dSJacob Faibussowitsch     PetscCall(DMPlexGetTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure));
16670034214SMatthew G. Knepley     for (c = 0; c < closureSize*2; c += 2) {
167527e7439SSatish Balay       for (q = 0; q < numAdj || ((void)(adj[numAdj++] = closure[c]),0); ++q) {
16870034214SMatthew G. Knepley         if (closure[c] == adj[q]) break;
16970034214SMatthew G. Knepley       }
17063a3b9bcSJacob Faibussowitsch       PetscCheck(numAdj <= maxAdjSize,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%" PetscInt_FMT ")", maxAdjSize);
17170034214SMatthew G. Knepley     }
1729566063dSJacob Faibussowitsch     PetscCall(DMPlexRestoreTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure));
17370034214SMatthew G. Knepley   }
1749566063dSJacob Faibussowitsch   PetscCall(DMPlexRestoreTransitiveClosure(dm, p, useClosure, &starSize, &star));
17570034214SMatthew G. Knepley   *adjSize = numAdj;
17670034214SMatthew G. Knepley   PetscFunctionReturn(0);
17770034214SMatthew G. Knepley }
17870034214SMatthew G. Knepley 
1795b317d89SToby Isaac PetscErrorCode DMPlexGetAdjacency_Internal(DM dm, PetscInt p, PetscBool useCone, PetscBool useTransitiveClosure, PetscBool useAnchors, PetscInt *adjSize, PetscInt *adj[])
18070034214SMatthew G. Knepley {
18179bad331SMatthew G. Knepley   static PetscInt asiz = 0;
1828b0b4c70SToby Isaac   PetscInt maxAnchors = 1;
1838b0b4c70SToby Isaac   PetscInt aStart = -1, aEnd = -1;
1848b0b4c70SToby Isaac   PetscInt maxAdjSize;
1858b0b4c70SToby Isaac   PetscSection aSec = NULL;
1868b0b4c70SToby Isaac   IS aIS = NULL;
1878b0b4c70SToby Isaac   const PetscInt *anchors;
1883c1f0c11SLawrence Mitchell   DM_Plex *mesh = (DM_Plex *)dm->data;
18970034214SMatthew G. Knepley 
19070034214SMatthew G. Knepley   PetscFunctionBeginHot;
1915b317d89SToby Isaac   if (useAnchors) {
1929566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAnchors(dm,&aSec,&aIS));
1938b0b4c70SToby Isaac     if (aSec) {
1949566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetMaxDof(aSec,&maxAnchors));
19524c766afSToby Isaac       maxAnchors = PetscMax(1,maxAnchors);
1969566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetChart(aSec,&aStart,&aEnd));
1979566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(aIS,&anchors));
1988b0b4c70SToby Isaac     }
1998b0b4c70SToby Isaac   }
20079bad331SMatthew G. Knepley   if (!*adj) {
20124c766afSToby Isaac     PetscInt depth, coneSeries, supportSeries, maxC, maxS, pStart, pEnd;
20279bad331SMatthew G. Knepley 
2039566063dSJacob Faibussowitsch     PetscCall(DMPlexGetChart(dm, &pStart,&pEnd));
2049566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepth(dm, &depth));
205412e9a14SMatthew G. Knepley     depth = PetscMax(depth, -depth);
2069566063dSJacob Faibussowitsch     PetscCall(DMPlexGetMaxSizes(dm, &maxC, &maxS));
20724c766afSToby Isaac     coneSeries    = (maxC > 1) ? ((PetscPowInt(maxC,depth+1)-1)/(maxC-1)) : depth+1;
20824c766afSToby Isaac     supportSeries = (maxS > 1) ? ((PetscPowInt(maxS,depth+1)-1)/(maxS-1)) : depth+1;
20924c766afSToby Isaac     asiz  = PetscMax(PetscPowInt(maxS,depth)*coneSeries,PetscPowInt(maxC,depth)*supportSeries);
2108b0b4c70SToby Isaac     asiz *= maxAnchors;
21124c766afSToby Isaac     asiz  = PetscMin(asiz,pEnd-pStart);
2129566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(asiz,adj));
21379bad331SMatthew G. Knepley   }
21479bad331SMatthew G. Knepley   if (*adjSize < 0) *adjSize = asiz;
2158b0b4c70SToby Isaac   maxAdjSize = *adjSize;
2163c1f0c11SLawrence Mitchell   if (mesh->useradjacency) {
2179566063dSJacob Faibussowitsch     PetscCall((*mesh->useradjacency)(dm, p, adjSize, *adj, mesh->useradjacencyctx));
2183c1f0c11SLawrence Mitchell   } else if (useTransitiveClosure) {
2199566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Transitive_Internal(dm, p, useCone, adjSize, *adj));
22070034214SMatthew G. Knepley   } else if (useCone) {
2219566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Cone_Internal(dm, p, adjSize, *adj));
22270034214SMatthew G. Knepley   } else {
2239566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency_Support_Internal(dm, p, adjSize, *adj));
22470034214SMatthew G. Knepley   }
2255b317d89SToby Isaac   if (useAnchors && aSec) {
2268b0b4c70SToby Isaac     PetscInt origSize = *adjSize;
2278b0b4c70SToby Isaac     PetscInt numAdj = origSize;
2288b0b4c70SToby Isaac     PetscInt i = 0, j;
2298b0b4c70SToby Isaac     PetscInt *orig = *adj;
2308b0b4c70SToby Isaac 
2318b0b4c70SToby Isaac     while (i < origSize) {
2328b0b4c70SToby Isaac       PetscInt p = orig[i];
2338b0b4c70SToby Isaac       PetscInt aDof = 0;
2348b0b4c70SToby Isaac 
2358b0b4c70SToby Isaac       if (p >= aStart && p < aEnd) {
2369566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetDof(aSec,p,&aDof));
2378b0b4c70SToby Isaac       }
2388b0b4c70SToby Isaac       if (aDof) {
2398b0b4c70SToby Isaac         PetscInt aOff;
2408b0b4c70SToby Isaac         PetscInt s, q;
2418b0b4c70SToby Isaac 
2428b0b4c70SToby Isaac         for (j = i + 1; j < numAdj; j++) {
2438b0b4c70SToby Isaac           orig[j - 1] = orig[j];
2448b0b4c70SToby Isaac         }
2458b0b4c70SToby Isaac         origSize--;
2468b0b4c70SToby Isaac         numAdj--;
2479566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(aSec,p,&aOff));
2488b0b4c70SToby Isaac         for (s = 0; s < aDof; ++s) {
249527e7439SSatish Balay           for (q = 0; q < numAdj || ((void)(orig[numAdj++] = anchors[aOff+s]),0); ++q) {
2508b0b4c70SToby Isaac             if (anchors[aOff+s] == orig[q]) break;
2518b0b4c70SToby Isaac           }
25263a3b9bcSJacob Faibussowitsch           PetscCheck(numAdj <= maxAdjSize,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%" PetscInt_FMT ")", maxAdjSize);
2538b0b4c70SToby Isaac         }
2548b0b4c70SToby Isaac       }
2558b0b4c70SToby Isaac       else {
2568b0b4c70SToby Isaac         i++;
2578b0b4c70SToby Isaac       }
2588b0b4c70SToby Isaac     }
2598b0b4c70SToby Isaac     *adjSize = numAdj;
2609566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(aIS,&anchors));
2618b0b4c70SToby Isaac   }
26270034214SMatthew G. Knepley   PetscFunctionReturn(0);
26370034214SMatthew G. Knepley }
26470034214SMatthew G. Knepley 
26570034214SMatthew G. Knepley /*@
26670034214SMatthew G. Knepley   DMPlexGetAdjacency - Return all points adjacent to the given point
26770034214SMatthew G. Knepley 
26870034214SMatthew G. Knepley   Input Parameters:
26970034214SMatthew G. Knepley + dm - The DM object
2706b867d5aSJose E. Roman - p  - The point
27170034214SMatthew G. Knepley 
2726b867d5aSJose E. Roman   Input/Output Parameters:
2736b867d5aSJose E. Roman + adjSize - The maximum size of adj if it is non-NULL, or PETSC_DETERMINE;
2746b867d5aSJose E. Roman             on output the number of adjacent points
2756b867d5aSJose E. Roman - adj - Either NULL so that the array is allocated, or an existing array with size adjSize;
2766b867d5aSJose E. Roman         on output contains the adjacent points
27770034214SMatthew G. Knepley 
27870034214SMatthew G. Knepley   Level: advanced
27970034214SMatthew G. Knepley 
28095452b02SPatrick Sanan   Notes:
28195452b02SPatrick Sanan     The user must PetscFree the adj array if it was not passed in.
28270034214SMatthew G. Knepley 
283db781477SPatrick Sanan .seealso: `DMSetAdjacency()`, `DMPlexDistribute()`, `DMCreateMatrix()`, `DMPlexPreallocateOperator()`
28470034214SMatthew G. Knepley @*/
28570034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacency(DM dm, PetscInt p, PetscInt *adjSize, PetscInt *adj[])
28670034214SMatthew G. Knepley {
2871cf84007SMatthew G. Knepley   PetscBool      useCone, useClosure, useAnchors;
28870034214SMatthew G. Knepley 
28970034214SMatthew G. Knepley   PetscFunctionBeginHot;
29070034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
291dadcf809SJacob Faibussowitsch   PetscValidIntPointer(adjSize,3);
29270034214SMatthew G. Knepley   PetscValidPointer(adj,4);
2939566063dSJacob Faibussowitsch   PetscCall(DMGetBasicAdjacency(dm, &useCone, &useClosure));
2949566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAdjacencyUseAnchors(dm, &useAnchors));
2959566063dSJacob Faibussowitsch   PetscCall(DMPlexGetAdjacency_Internal(dm, p, useCone, useClosure, useAnchors, adjSize, adj));
29670034214SMatthew G. Knepley   PetscFunctionReturn(0);
29770034214SMatthew G. Knepley }
29808633170SToby Isaac 
299b0a623aaSMatthew G. Knepley /*@
300b0a623aaSMatthew G. Knepley   DMPlexCreateTwoSidedProcessSF - Create an SF which just has process connectivity
301b0a623aaSMatthew G. Knepley 
302d083f849SBarry Smith   Collective on dm
303b0a623aaSMatthew G. Knepley 
304b0a623aaSMatthew G. Knepley   Input Parameters:
305b0a623aaSMatthew G. Knepley + dm      - The DM
3066b867d5aSJose E. Roman . sfPoint - The PetscSF which encodes point connectivity
3076b867d5aSJose E. Roman . rootRankSection -
3086b867d5aSJose E. Roman . rootRanks -
3096b867d5aSJose E. Roman . leftRankSection -
3106b867d5aSJose E. Roman - leafRanks -
311b0a623aaSMatthew G. Knepley 
312b0a623aaSMatthew G. Knepley   Output Parameters:
313b0a623aaSMatthew G. Knepley + processRanks - A list of process neighbors, or NULL
314b0a623aaSMatthew G. Knepley - sfProcess    - An SF encoding the two-sided process connectivity, or NULL
315b0a623aaSMatthew G. Knepley 
316b0a623aaSMatthew G. Knepley   Level: developer
317b0a623aaSMatthew G. Knepley 
318db781477SPatrick Sanan .seealso: `PetscSFCreate()`, `DMPlexCreateProcessSF()`
319b0a623aaSMatthew G. Knepley @*/
320b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM dm, PetscSF sfPoint, PetscSection rootRankSection, IS rootRanks, PetscSection leafRankSection, IS leafRanks, IS *processRanks, PetscSF *sfProcess)
321b0a623aaSMatthew G. Knepley {
322b0a623aaSMatthew G. Knepley   const PetscSFNode *remotePoints;
323b0a623aaSMatthew G. Knepley   PetscInt          *localPointsNew;
324b0a623aaSMatthew G. Knepley   PetscSFNode       *remotePointsNew;
325b0a623aaSMatthew G. Knepley   const PetscInt    *nranks;
326b0a623aaSMatthew G. Knepley   PetscInt          *ranksNew;
327b0a623aaSMatthew G. Knepley   PetscBT            neighbors;
328b0a623aaSMatthew G. Knepley   PetscInt           pStart, pEnd, p, numLeaves, l, numNeighbors, n;
3299852e123SBarry Smith   PetscMPIInt        size, proc, rank;
330b0a623aaSMatthew G. Knepley 
331b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
332b0a623aaSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
333b0a623aaSMatthew G. Knepley   PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2);
334064a246eSJacob Faibussowitsch   if (processRanks) {PetscValidPointer(processRanks, 7);}
335064a246eSJacob Faibussowitsch   if (sfProcess)    {PetscValidPointer(sfProcess, 8);}
3369566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size));
3379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank));
3389566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfPoint, NULL, &numLeaves, NULL, &remotePoints));
3399566063dSJacob Faibussowitsch   PetscCall(PetscBTCreate(size, &neighbors));
3409566063dSJacob Faibussowitsch   PetscCall(PetscBTMemzero(size, neighbors));
341b0a623aaSMatthew G. Knepley   /* Compute root-to-leaf process connectivity */
3429566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(rootRankSection, &pStart, &pEnd));
3439566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(rootRanks, &nranks));
344b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
345b0a623aaSMatthew G. Knepley     PetscInt ndof, noff, n;
346b0a623aaSMatthew G. Knepley 
3479566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(rootRankSection, p, &ndof));
3489566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(rootRankSection, p, &noff));
3499566063dSJacob Faibussowitsch     for (n = 0; n < ndof; ++n) PetscCall(PetscBTSet(neighbors, nranks[noff+n]));
350b0a623aaSMatthew G. Knepley   }
3519566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(rootRanks, &nranks));
352b0a623aaSMatthew G. Knepley   /* Compute leaf-to-neighbor process connectivity */
3539566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(leafRankSection, &pStart, &pEnd));
3549566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(leafRanks, &nranks));
355b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
356b0a623aaSMatthew G. Knepley     PetscInt ndof, noff, n;
357b0a623aaSMatthew G. Knepley 
3589566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(leafRankSection, p, &ndof));
3599566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(leafRankSection, p, &noff));
3609566063dSJacob Faibussowitsch     for (n = 0; n < ndof; ++n) PetscCall(PetscBTSet(neighbors, nranks[noff+n]));
361b0a623aaSMatthew G. Knepley   }
3629566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(leafRanks, &nranks));
363b0a623aaSMatthew G. Knepley   /* Compute leaf-to-root process connectivity */
364b0a623aaSMatthew G. Knepley   for (l = 0; l < numLeaves; ++l) {PetscBTSet(neighbors, remotePoints[l].rank);}
365b0a623aaSMatthew G. Knepley   /* Calculate edges */
366b0a623aaSMatthew G. Knepley   PetscBTClear(neighbors, rank);
3679852e123SBarry Smith   for (proc = 0, numNeighbors = 0; proc < size; ++proc) {if (PetscBTLookup(neighbors, proc)) ++numNeighbors;}
3689566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numNeighbors, &ranksNew));
3699566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numNeighbors, &localPointsNew));
3709566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numNeighbors, &remotePointsNew));
3719852e123SBarry Smith   for (proc = 0, n = 0; proc < size; ++proc) {
372b0a623aaSMatthew G. Knepley     if (PetscBTLookup(neighbors, proc)) {
373b0a623aaSMatthew G. Knepley       ranksNew[n]              = proc;
374b0a623aaSMatthew G. Knepley       localPointsNew[n]        = proc;
375b0a623aaSMatthew G. Knepley       remotePointsNew[n].index = rank;
376b0a623aaSMatthew G. Knepley       remotePointsNew[n].rank  = proc;
377b0a623aaSMatthew G. Knepley       ++n;
378b0a623aaSMatthew G. Knepley     }
379b0a623aaSMatthew G. Knepley   }
3809566063dSJacob Faibussowitsch   PetscCall(PetscBTDestroy(&neighbors));
3819566063dSJacob Faibussowitsch   if (processRanks) PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)dm), numNeighbors, ranksNew, PETSC_OWN_POINTER, processRanks));
3829566063dSJacob Faibussowitsch   else              PetscCall(PetscFree(ranksNew));
383b0a623aaSMatthew G. Knepley   if (sfProcess) {
3849566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess));
3859566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) *sfProcess, "Two-Sided Process SF"));
3869566063dSJacob Faibussowitsch     PetscCall(PetscSFSetFromOptions(*sfProcess));
3879566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(*sfProcess, size, numNeighbors, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER));
388b0a623aaSMatthew G. Knepley   }
389b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
390b0a623aaSMatthew G. Knepley }
391b0a623aaSMatthew G. Knepley 
392b0a623aaSMatthew G. Knepley /*@
393b0a623aaSMatthew G. Knepley   DMPlexDistributeOwnership - Compute owner information for shared points. This basically gets two-sided for an SF.
394b0a623aaSMatthew G. Knepley 
395d083f849SBarry Smith   Collective on dm
396b0a623aaSMatthew G. Knepley 
397b0a623aaSMatthew G. Knepley   Input Parameter:
398b0a623aaSMatthew G. Knepley . dm - The DM
399b0a623aaSMatthew G. Knepley 
400b0a623aaSMatthew G. Knepley   Output Parameters:
401b0a623aaSMatthew G. Knepley + rootSection - The number of leaves for a given root point
402b0a623aaSMatthew G. Knepley . rootrank    - The rank of each edge into the root point
403b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point
404b0a623aaSMatthew G. Knepley - leafrank    - The rank of each process sharing a leaf point
405b0a623aaSMatthew G. Knepley 
406b0a623aaSMatthew G. Knepley   Level: developer
407b0a623aaSMatthew G. Knepley 
408db781477SPatrick Sanan .seealso: `DMPlexCreateOverlapLabel()`, `DMPlexDistribute()`, `DMPlexDistributeOverlap()`
409b0a623aaSMatthew G. Knepley @*/
410b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexDistributeOwnership(DM dm, PetscSection rootSection, IS *rootrank, PetscSection leafSection, IS *leafrank)
411b0a623aaSMatthew G. Knepley {
412b0a623aaSMatthew G. Knepley   MPI_Comm        comm;
413b0a623aaSMatthew G. Knepley   PetscSF         sfPoint;
414b0a623aaSMatthew G. Knepley   const PetscInt *rootdegree;
415b0a623aaSMatthew G. Knepley   PetscInt       *myrank, *remoterank;
416b0a623aaSMatthew G. Knepley   PetscInt        pStart, pEnd, p, nedges;
417b0a623aaSMatthew G. Knepley   PetscMPIInt     rank;
418b0a623aaSMatthew G. Knepley 
419b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
4209566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
4219566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
4229566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
4239566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &sfPoint));
424b0a623aaSMatthew G. Knepley   /* Compute number of leaves for each root */
4259566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) rootSection, "Root Section"));
4269566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetChart(rootSection, pStart, pEnd));
4279566063dSJacob Faibussowitsch   PetscCall(PetscSFComputeDegreeBegin(sfPoint, &rootdegree));
4289566063dSJacob Faibussowitsch   PetscCall(PetscSFComputeDegreeEnd(sfPoint, &rootdegree));
4299566063dSJacob Faibussowitsch   for (p = pStart; p < pEnd; ++p) PetscCall(PetscSectionSetDof(rootSection, p, rootdegree[p-pStart]));
4309566063dSJacob Faibussowitsch   PetscCall(PetscSectionSetUp(rootSection));
431b0a623aaSMatthew G. Knepley   /* Gather rank of each leaf to root */
4329566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(rootSection, &nedges));
4339566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(pEnd-pStart, &myrank));
4349566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nedges,  &remoterank));
435b0a623aaSMatthew G. Knepley   for (p = 0; p < pEnd-pStart; ++p) myrank[p] = rank;
4369566063dSJacob Faibussowitsch   PetscCall(PetscSFGatherBegin(sfPoint, MPIU_INT, myrank, remoterank));
4379566063dSJacob Faibussowitsch   PetscCall(PetscSFGatherEnd(sfPoint, MPIU_INT, myrank, remoterank));
4389566063dSJacob Faibussowitsch   PetscCall(PetscFree(myrank));
4399566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(comm, nedges, remoterank, PETSC_OWN_POINTER, rootrank));
440b0a623aaSMatthew G. Knepley   /* Distribute remote ranks to leaves */
4419566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) leafSection, "Leaf Section"));
4429566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeFieldIS(dm, sfPoint, rootSection, *rootrank, leafSection, leafrank));
443b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
444b0a623aaSMatthew G. Knepley }
445b0a623aaSMatthew G. Knepley 
446278397a0SMatthew G. Knepley /*@C
447dccdeccaSVaclav Hapla   DMPlexCreateOverlapLabel - Compute owner information for shared points. This basically gets two-sided for an SF.
448b0a623aaSMatthew G. Knepley 
449d083f849SBarry Smith   Collective on dm
450b0a623aaSMatthew G. Knepley 
451b0a623aaSMatthew G. Knepley   Input Parameters:
452b0a623aaSMatthew G. Knepley + dm          - The DM
45324d039d7SMichael Lange . levels      - Number of overlap levels
454b0a623aaSMatthew G. Knepley . rootSection - The number of leaves for a given root point
455b0a623aaSMatthew G. Knepley . rootrank    - The rank of each edge into the root point
456b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point
457b0a623aaSMatthew G. Knepley - leafrank    - The rank of each process sharing a leaf point
458b0a623aaSMatthew G. Knepley 
459064ec1f3Sprj-   Output Parameter:
460b4ec6ac8SBarry Smith . ovLabel     - DMLabel containing remote overlap contributions as point/rank pairings
461b0a623aaSMatthew G. Knepley 
462b0a623aaSMatthew G. Knepley   Level: developer
463b0a623aaSMatthew G. Knepley 
464db781477SPatrick Sanan .seealso: `DMPlexDistributeOwnership()`, `DMPlexDistribute()`
465b0a623aaSMatthew G. Knepley @*/
466dccdeccaSVaclav Hapla PetscErrorCode DMPlexCreateOverlapLabel(DM dm, PetscInt levels, PetscSection rootSection, IS rootrank, PetscSection leafSection, IS leafrank, DMLabel *ovLabel)
467b0a623aaSMatthew G. Knepley {
468e540f424SMichael Lange   MPI_Comm           comm;
469b0a623aaSMatthew G. Knepley   DMLabel            ovAdjByRank; /* A DMLabel containing all points adjacent to shared points, separated by rank (value in label) */
470874ddda9SLisandro Dalcin   PetscSF            sfPoint;
471b0a623aaSMatthew G. Knepley   const PetscSFNode *remote;
472b0a623aaSMatthew G. Knepley   const PetscInt    *local;
4731fd9873aSMichael Lange   const PetscInt    *nrank, *rrank;
474b0a623aaSMatthew G. Knepley   PetscInt          *adj = NULL;
4751fd9873aSMichael Lange   PetscInt           pStart, pEnd, p, sStart, sEnd, nleaves, l;
4769852e123SBarry Smith   PetscMPIInt        rank, size;
47731bc6364SLisandro Dalcin   PetscBool          flg;
478b0a623aaSMatthew G. Knepley 
479b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
4806ba1a4c7SVaclav Hapla   *ovLabel = NULL;
4819566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
4829566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
4839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
4846ba1a4c7SVaclav Hapla   if (size == 1) PetscFunctionReturn(0);
4859566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &sfPoint));
4869566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
487d1674c6dSMatthew Knepley   if (!levels) {
488d1674c6dSMatthew Knepley     /* Add owned points */
4899566063dSJacob Faibussowitsch     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Overlap label", ovLabel));
4909566063dSJacob Faibussowitsch     for (p = pStart; p < pEnd; ++p) PetscCall(DMLabelSetValue(*ovLabel, p, rank));
491d1674c6dSMatthew Knepley     PetscFunctionReturn(0);
492d1674c6dSMatthew Knepley   }
4939566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetChart(leafSection, &sStart, &sEnd));
4949566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote));
4959566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Overlap adjacency", &ovAdjByRank));
496b0a623aaSMatthew G. Knepley   /* Handle leaves: shared with the root point */
497b0a623aaSMatthew G. Knepley   for (l = 0; l < nleaves; ++l) {
498b0a623aaSMatthew G. Knepley     PetscInt adjSize = PETSC_DETERMINE, a;
499b0a623aaSMatthew G. Knepley 
5009566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency(dm, local ? local[l] : l, &adjSize, &adj));
5019566063dSJacob Faibussowitsch     for (a = 0; a < adjSize; ++a) PetscCall(DMLabelSetValue(ovAdjByRank, adj[a], remote[l].rank));
502b0a623aaSMatthew G. Knepley   }
5039566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(rootrank, &rrank));
5049566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(leafrank, &nrank));
505b0a623aaSMatthew G. Knepley   /* Handle roots */
506b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
507b0a623aaSMatthew G. Knepley     PetscInt adjSize = PETSC_DETERMINE, neighbors = 0, noff, n, a;
508b0a623aaSMatthew G. Knepley 
509b0a623aaSMatthew G. Knepley     if ((p >= sStart) && (p < sEnd)) {
510b0a623aaSMatthew G. Knepley       /* Some leaves share a root with other leaves on different processes */
5119566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(leafSection, p, &neighbors));
512b0a623aaSMatthew G. Knepley       if (neighbors) {
5139566063dSJacob Faibussowitsch         PetscCall(PetscSectionGetOffset(leafSection, p, &noff));
5149566063dSJacob Faibussowitsch         PetscCall(DMPlexGetAdjacency(dm, p, &adjSize, &adj));
515b0a623aaSMatthew G. Knepley         for (n = 0; n < neighbors; ++n) {
516b0a623aaSMatthew G. Knepley           const PetscInt remoteRank = nrank[noff+n];
517b0a623aaSMatthew G. Knepley 
518b0a623aaSMatthew G. Knepley           if (remoteRank == rank) continue;
5199566063dSJacob Faibussowitsch           for (a = 0; a < adjSize; ++a) PetscCall(DMLabelSetValue(ovAdjByRank, adj[a], remoteRank));
520b0a623aaSMatthew G. Knepley         }
521b0a623aaSMatthew G. Knepley       }
522b0a623aaSMatthew G. Knepley     }
523b0a623aaSMatthew G. Knepley     /* Roots are shared with leaves */
5249566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetDof(rootSection, p, &neighbors));
525b0a623aaSMatthew G. Knepley     if (!neighbors) continue;
5269566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetOffset(rootSection, p, &noff));
5279566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacency(dm, p, &adjSize, &adj));
528b0a623aaSMatthew G. Knepley     for (n = 0; n < neighbors; ++n) {
529b0a623aaSMatthew G. Knepley       const PetscInt remoteRank = rrank[noff+n];
530b0a623aaSMatthew G. Knepley 
531b0a623aaSMatthew G. Knepley       if (remoteRank == rank) continue;
5329566063dSJacob Faibussowitsch       for (a = 0; a < adjSize; ++a) PetscCall(DMLabelSetValue(ovAdjByRank, adj[a], remoteRank));
533b0a623aaSMatthew G. Knepley     }
534b0a623aaSMatthew G. Knepley   }
5359566063dSJacob Faibussowitsch   PetscCall(PetscFree(adj));
5369566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(rootrank, &rrank));
5379566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(leafrank, &nrank));
53824d039d7SMichael Lange   /* Add additional overlap levels */
539be200f8dSMichael Lange   for (l = 1; l < levels; l++) {
540be200f8dSMichael Lange     /* Propagate point donations over SF to capture remote connections */
5419566063dSJacob Faibussowitsch     PetscCall(DMPlexPartitionLabelPropagate(dm, ovAdjByRank));
542be200f8dSMichael Lange     /* Add next level of point donations to the label */
5439566063dSJacob Faibussowitsch     PetscCall(DMPlexPartitionLabelAdjacency(dm, ovAdjByRank));
544be200f8dSMichael Lange   }
54526a7d390SMatthew G. Knepley   /* We require the closure in the overlap */
5469566063dSJacob Faibussowitsch   PetscCall(DMPlexPartitionLabelClosure(dm, ovAdjByRank));
5479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-overlap_view", &flg));
548e540f424SMichael Lange   if (flg) {
549825f8a23SLisandro Dalcin     PetscViewer viewer;
5509566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)dm), &viewer));
5519566063dSJacob Faibussowitsch     PetscCall(DMLabelView(ovAdjByRank, viewer));
552b0a623aaSMatthew G. Knepley   }
553874ddda9SLisandro Dalcin   /* Invert sender to receiver label */
5549566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Overlap label", ovLabel));
5559566063dSJacob Faibussowitsch   PetscCall(DMPlexPartitionLabelInvert(dm, ovAdjByRank, NULL, *ovLabel));
556a9f1d5b2SMichael Lange   /* Add owned points, except for shared local points */
5579566063dSJacob Faibussowitsch   for (p = pStart; p < pEnd; ++p) PetscCall(DMLabelSetValue(*ovLabel, p, rank));
558e540f424SMichael Lange   for (l = 0; l < nleaves; ++l) {
5599566063dSJacob Faibussowitsch     PetscCall(DMLabelClearValue(*ovLabel, local[l], rank));
5609566063dSJacob Faibussowitsch     PetscCall(DMLabelSetValue(*ovLabel, remote[l].index, remote[l].rank));
561e540f424SMichael Lange   }
562e540f424SMichael Lange   /* Clean up */
5639566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&ovAdjByRank));
564b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
565b0a623aaSMatthew G. Knepley }
56670034214SMatthew G. Knepley 
56724cc2ca5SMatthew G. Knepley /*@C
56824cc2ca5SMatthew G. Knepley   DMPlexCreateOverlapMigrationSF - Create an SF describing the new mesh distribution to make the overlap described by the input SF
56924cc2ca5SMatthew G. Knepley 
570d083f849SBarry Smith   Collective on dm
57124cc2ca5SMatthew G. Knepley 
57224cc2ca5SMatthew G. Knepley   Input Parameters:
57324cc2ca5SMatthew G. Knepley + dm          - The DM
57424cc2ca5SMatthew G. Knepley - overlapSF   - The SF mapping ghost points in overlap to owner points on other processes
57524cc2ca5SMatthew G. Knepley 
576064ec1f3Sprj-   Output Parameter:
577a2b725a8SWilliam Gropp . migrationSF - An SF that maps original points in old locations to points in new locations
57824cc2ca5SMatthew G. Knepley 
57924cc2ca5SMatthew G. Knepley   Level: developer
58024cc2ca5SMatthew G. Knepley 
581db781477SPatrick Sanan .seealso: `DMPlexCreateOverlapLabel()`, `DMPlexDistributeOverlap()`, `DMPlexDistribute()`
58224cc2ca5SMatthew G. Knepley @*/
58346f9b1c3SMichael Lange PetscErrorCode DMPlexCreateOverlapMigrationSF(DM dm, PetscSF overlapSF, PetscSF *migrationSF)
58446f9b1c3SMichael Lange {
58546f9b1c3SMichael Lange   MPI_Comm           comm;
5869852e123SBarry Smith   PetscMPIInt        rank, size;
58746f9b1c3SMichael Lange   PetscInt           d, dim, p, pStart, pEnd, nroots, nleaves, newLeaves, point, numSharedPoints;
58846f9b1c3SMichael Lange   PetscInt          *pointDepths, *remoteDepths, *ilocal;
58946f9b1c3SMichael Lange   PetscInt          *depthRecv, *depthShift, *depthIdx;
59046f9b1c3SMichael Lange   PetscSFNode       *iremote;
59146f9b1c3SMichael Lange   PetscSF            pointSF;
59246f9b1c3SMichael Lange   const PetscInt    *sharedLocal;
59346f9b1c3SMichael Lange   const PetscSFNode *overlapRemote, *sharedRemote;
59446f9b1c3SMichael Lange 
59546f9b1c3SMichael Lange   PetscFunctionBegin;
59646f9b1c3SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5979566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
5989566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
5999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
6009566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
60146f9b1c3SMichael Lange 
60246f9b1c3SMichael Lange   /* Before building the migration SF we need to know the new stratum offsets */
6039566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(overlapSF, &nroots, &nleaves, NULL, &overlapRemote));
6049566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths));
60546f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
6069566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
60746f9b1c3SMichael Lange     for (p=pStart; p<pEnd; p++) pointDepths[p] = d;
60846f9b1c3SMichael Lange   }
60946f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) remoteDepths[p] = -1;
6109566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(overlapSF, MPIU_INT, pointDepths, remoteDepths,MPI_REPLACE));
6119566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(overlapSF, MPIU_INT, pointDepths, remoteDepths,MPI_REPLACE));
61246f9b1c3SMichael Lange 
6132d4ee042Sprj-   /* Count received points in each stratum and compute the internal strata shift */
6149566063dSJacob Faibussowitsch   PetscCall(PetscMalloc3(dim+1, &depthRecv, dim+1, &depthShift, dim+1, &depthIdx));
61546f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) depthRecv[d]=0;
61646f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) depthRecv[remoteDepths[p]]++;
61746f9b1c3SMichael Lange   depthShift[dim] = 0;
61846f9b1c3SMichael Lange   for (d=0; d<dim; d++) depthShift[d] = depthRecv[dim];
61946f9b1c3SMichael Lange   for (d=1; d<dim; d++) depthShift[d] += depthRecv[0];
62046f9b1c3SMichael Lange   for (d=dim-2; d>0; d--) depthShift[d] += depthRecv[d+1];
62146f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
6229566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
62346f9b1c3SMichael Lange     depthIdx[d] = pStart + depthShift[d];
62446f9b1c3SMichael Lange   }
62546f9b1c3SMichael Lange 
62646f9b1c3SMichael Lange   /* Form the overlap SF build an SF that describes the full overlap migration SF */
6279566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
62846f9b1c3SMichael Lange   newLeaves = pEnd - pStart + nleaves;
6299566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(newLeaves, &ilocal));
6309566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(newLeaves, &iremote));
63146f9b1c3SMichael Lange   /* First map local points to themselves */
63246f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
6339566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
63446f9b1c3SMichael Lange     for (p=pStart; p<pEnd; p++) {
63546f9b1c3SMichael Lange       point = p + depthShift[d];
63646f9b1c3SMichael Lange       ilocal[point] = point;
63746f9b1c3SMichael Lange       iremote[point].index = p;
63846f9b1c3SMichael Lange       iremote[point].rank = rank;
63946f9b1c3SMichael Lange       depthIdx[d]++;
64046f9b1c3SMichael Lange     }
64146f9b1c3SMichael Lange   }
64246f9b1c3SMichael Lange 
64346f9b1c3SMichael Lange   /* Add in the remote roots for currently shared points */
6449566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &pointSF));
6459566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(pointSF, NULL, &numSharedPoints, &sharedLocal, &sharedRemote));
64646f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
6479566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
64846f9b1c3SMichael Lange     for (p=0; p<numSharedPoints; p++) {
64946f9b1c3SMichael Lange       if (pStart <= sharedLocal[p] && sharedLocal[p] < pEnd) {
65046f9b1c3SMichael Lange         point = sharedLocal[p] + depthShift[d];
65146f9b1c3SMichael Lange         iremote[point].index = sharedRemote[p].index;
65246f9b1c3SMichael Lange         iremote[point].rank = sharedRemote[p].rank;
65346f9b1c3SMichael Lange       }
65446f9b1c3SMichael Lange     }
65546f9b1c3SMichael Lange   }
65646f9b1c3SMichael Lange 
65746f9b1c3SMichael Lange   /* Now add the incoming overlap points */
65846f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) {
65946f9b1c3SMichael Lange     point = depthIdx[remoteDepths[p]];
66046f9b1c3SMichael Lange     ilocal[point] = point;
66146f9b1c3SMichael Lange     iremote[point].index = overlapRemote[p].index;
66246f9b1c3SMichael Lange     iremote[point].rank = overlapRemote[p].rank;
66346f9b1c3SMichael Lange     depthIdx[remoteDepths[p]]++;
66446f9b1c3SMichael Lange   }
6659566063dSJacob Faibussowitsch   PetscCall(PetscFree2(pointDepths,remoteDepths));
66646f9b1c3SMichael Lange 
6679566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm, migrationSF));
6689566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *migrationSF, "Overlap Migration SF"));
6699566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(*migrationSF));
6709566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
6719566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*migrationSF, pEnd-pStart, newLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER));
67246f9b1c3SMichael Lange 
6739566063dSJacob Faibussowitsch   PetscCall(PetscFree3(depthRecv, depthShift, depthIdx));
67470034214SMatthew G. Knepley   PetscFunctionReturn(0);
67570034214SMatthew G. Knepley }
67670034214SMatthew G. Knepley 
677a9f1d5b2SMichael Lange /*@
678f0e73a3dSToby Isaac   DMPlexStratifyMigrationSF - Rearrange the leaves of a migration sf for stratification.
679a9f1d5b2SMichael Lange 
680064ec1f3Sprj-   Input Parameters:
681a9f1d5b2SMichael Lange + dm          - The DM
682a9f1d5b2SMichael Lange - sf          - A star forest with non-ordered leaves, usually defining a DM point migration
683a9f1d5b2SMichael Lange 
684a9f1d5b2SMichael Lange   Output Parameter:
685a9f1d5b2SMichael Lange . migrationSF - A star forest with added leaf indirection that ensures the resulting DM is stratified
686a9f1d5b2SMichael Lange 
687412e9a14SMatthew G. Knepley   Note:
688412e9a14SMatthew G. Knepley   This lexicographically sorts by (depth, cellType)
689412e9a14SMatthew G. Knepley 
690a9f1d5b2SMichael Lange   Level: developer
691a9f1d5b2SMichael Lange 
692db781477SPatrick Sanan .seealso: `DMPlexPartitionLabelCreateSF()`, `DMPlexDistribute()`, `DMPlexDistributeOverlap()`
693a9f1d5b2SMichael Lange @*/
694a9f1d5b2SMichael Lange PetscErrorCode DMPlexStratifyMigrationSF(DM dm, PetscSF sf, PetscSF *migrationSF)
695a9f1d5b2SMichael Lange {
696a9f1d5b2SMichael Lange   MPI_Comm           comm;
6979852e123SBarry Smith   PetscMPIInt        rank, size;
698412e9a14SMatthew G. Knepley   PetscInt           d, ldepth, depth, dim, p, pStart, pEnd, nroots, nleaves;
699412e9a14SMatthew G. Knepley   PetscSFNode       *pointDepths, *remoteDepths;
700412e9a14SMatthew G. Knepley   PetscInt          *ilocal;
701a9f1d5b2SMichael Lange   PetscInt          *depthRecv, *depthShift, *depthIdx;
702412e9a14SMatthew G. Knepley   PetscInt          *ctRecv,    *ctShift,    *ctIdx;
703a9f1d5b2SMichael Lange   const PetscSFNode *iremote;
704a9f1d5b2SMichael Lange 
705a9f1d5b2SMichael Lange   PetscFunctionBegin;
706a9f1d5b2SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
7079566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
7089566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
7099566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
7109566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &ldepth));
7119566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
7121c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&ldepth, &depth, 1, MPIU_INT, MPI_MAX, comm));
71363a3b9bcSJacob Faibussowitsch   PetscCheck(!(ldepth >= 0) || !(depth != ldepth),PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent Plex depth %" PetscInt_FMT " != %" PetscInt_FMT, ldepth, depth);
7149566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_PartStratSF,dm,0,0,0));
715a9f1d5b2SMichael Lange 
716a9f1d5b2SMichael Lange   /* Before building the migration SF we need to know the new stratum offsets */
7179566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sf, &nroots, &nleaves, NULL, &iremote));
7189566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths));
7197fab53ddSMatthew G. Knepley   for (d = 0; d < depth+1; ++d) {
7209566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthStratum(dm, d, &pStart, &pEnd));
721f0e73a3dSToby Isaac     for (p = pStart; p < pEnd; ++p) {
722412e9a14SMatthew G. Knepley       DMPolytopeType ct;
723f0e73a3dSToby Isaac 
7249566063dSJacob Faibussowitsch       PetscCall(DMPlexGetCellType(dm, p, &ct));
725412e9a14SMatthew G. Knepley       pointDepths[p].index = d;
726412e9a14SMatthew G. Knepley       pointDepths[p].rank  = ct;
727f0e73a3dSToby Isaac     }
728412e9a14SMatthew G. Knepley   }
729412e9a14SMatthew G. Knepley   for (p = 0; p < nleaves; ++p) {remoteDepths[p].index = -1; remoteDepths[p].rank = -1;}
7309566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(sf, MPIU_2INT, pointDepths, remoteDepths,MPI_REPLACE));
7319566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(sf, MPIU_2INT, pointDepths, remoteDepths,MPI_REPLACE));
732412e9a14SMatthew G. Knepley   /* Count received points in each stratum and compute the internal strata shift */
7339566063dSJacob Faibussowitsch   PetscCall(PetscCalloc6(depth+1, &depthRecv, depth+1, &depthShift, depth+1, &depthIdx, DM_NUM_POLYTOPES, &ctRecv, DM_NUM_POLYTOPES, &ctShift, DM_NUM_POLYTOPES, &ctIdx));
734412e9a14SMatthew G. Knepley   for (p = 0; p < nleaves; ++p) {
735412e9a14SMatthew G. Knepley     if (remoteDepths[p].rank < 0) {
736412e9a14SMatthew G. Knepley       ++depthRecv[remoteDepths[p].index];
737412e9a14SMatthew G. Knepley     } else {
738412e9a14SMatthew G. Knepley       ++ctRecv[remoteDepths[p].rank];
739412e9a14SMatthew G. Knepley     }
740412e9a14SMatthew G. Knepley   }
741412e9a14SMatthew G. Knepley   {
742412e9a14SMatthew G. Knepley     PetscInt depths[4], dims[4], shift = 0, i, c;
743412e9a14SMatthew G. Knepley 
7448238f61eSMatthew G. Knepley     /* Cells (depth), Vertices (0), Faces (depth-1), Edges (1)
7458238f61eSMatthew G. Knepley          Consider DM_POLYTOPE_FV_GHOST and DM_POLYTOPE_INTERIOR_GHOST as cells
7468238f61eSMatthew G. Knepley      */
747412e9a14SMatthew G. Knepley     depths[0] = depth; depths[1] = 0; depths[2] = depth-1; depths[3] = 1;
748412e9a14SMatthew G. Knepley     dims[0]   = dim;   dims[1]   = 0; dims[2]   = dim-1;   dims[3]   = 1;
749412e9a14SMatthew G. Knepley     for (i = 0; i <= depth; ++i) {
750412e9a14SMatthew G. Knepley       const PetscInt dep = depths[i];
751412e9a14SMatthew G. Knepley       const PetscInt dim = dims[i];
752412e9a14SMatthew G. Knepley 
753412e9a14SMatthew G. Knepley       for (c = 0; c < DM_NUM_POLYTOPES; ++c) {
7548238f61eSMatthew G. Knepley         if (DMPolytopeTypeGetDim((DMPolytopeType) c) != dim && !(i == 0 && (c == DM_POLYTOPE_FV_GHOST || c == DM_POLYTOPE_INTERIOR_GHOST))) continue;
755412e9a14SMatthew G. Knepley         ctShift[c] = shift;
756412e9a14SMatthew G. Knepley         shift     += ctRecv[c];
757412e9a14SMatthew G. Knepley       }
758412e9a14SMatthew G. Knepley       depthShift[dep] = shift;
759412e9a14SMatthew G. Knepley       shift          += depthRecv[dep];
760412e9a14SMatthew G. Knepley     }
761412e9a14SMatthew G. Knepley     for (c = 0; c < DM_NUM_POLYTOPES; ++c) {
762412e9a14SMatthew G. Knepley       const PetscInt ctDim = DMPolytopeTypeGetDim((DMPolytopeType) c);
763412e9a14SMatthew G. Knepley 
7648238f61eSMatthew G. Knepley       if ((ctDim < 0 || ctDim > dim) && (c != DM_POLYTOPE_FV_GHOST && c != DM_POLYTOPE_INTERIOR_GHOST)) {
765412e9a14SMatthew G. Knepley         ctShift[c] = shift;
766412e9a14SMatthew G. Knepley         shift     += ctRecv[c];
767412e9a14SMatthew G. Knepley       }
768412e9a14SMatthew G. Knepley     }
769412e9a14SMatthew G. Knepley   }
770a9f1d5b2SMichael Lange   /* Derive a new local permutation based on stratified indices */
7719566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(nleaves, &ilocal));
7727fab53ddSMatthew G. Knepley   for (p = 0; p < nleaves; ++p) {
773412e9a14SMatthew G. Knepley     const PetscInt       dep = remoteDepths[p].index;
774412e9a14SMatthew G. Knepley     const DMPolytopeType ct  = (DMPolytopeType) remoteDepths[p].rank;
7757fab53ddSMatthew G. Knepley 
776412e9a14SMatthew G. Knepley     if ((PetscInt) ct < 0) {
7777fab53ddSMatthew G. Knepley       ilocal[p] = depthShift[dep] + depthIdx[dep];
778412e9a14SMatthew G. Knepley       ++depthIdx[dep];
779412e9a14SMatthew G. Knepley     } else {
780412e9a14SMatthew G. Knepley       ilocal[p] = ctShift[ct] + ctIdx[ct];
781412e9a14SMatthew G. Knepley       ++ctIdx[ct];
782412e9a14SMatthew G. Knepley     }
783a9f1d5b2SMichael Lange   }
7849566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm, migrationSF));
7859566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *migrationSF, "Migration SF"));
7869566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*migrationSF, nroots, nleaves, ilocal, PETSC_OWN_POINTER, (PetscSFNode*)iremote, PETSC_COPY_VALUES));
7879566063dSJacob Faibussowitsch   PetscCall(PetscFree2(pointDepths,remoteDepths));
7889566063dSJacob Faibussowitsch   PetscCall(PetscFree6(depthRecv, depthShift, depthIdx, ctRecv, ctShift, ctIdx));
7899566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_PartStratSF,dm,0,0,0));
790a9f1d5b2SMichael Lange   PetscFunctionReturn(0);
791a9f1d5b2SMichael Lange }
792a9f1d5b2SMichael Lange 
79370034214SMatthew G. Knepley /*@
79470034214SMatthew G. Knepley   DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
79570034214SMatthew G. Knepley 
796d083f849SBarry Smith   Collective on dm
79770034214SMatthew G. Knepley 
79870034214SMatthew G. Knepley   Input Parameters:
79970034214SMatthew G. Knepley + dm - The DMPlex object
80070034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
80170034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
802cb15cd0eSMatthew G. Knepley - originalVec - The existing data in a local vector
80370034214SMatthew G. Knepley 
80470034214SMatthew G. Knepley   Output Parameters:
80570034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout
806cb15cd0eSMatthew G. Knepley - newVec - The new data in a local vector
80770034214SMatthew G. Knepley 
80870034214SMatthew G. Knepley   Level: developer
80970034214SMatthew G. Knepley 
810db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeFieldIS()`, `DMPlexDistributeData()`
81170034214SMatthew G. Knepley @*/
81270034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec)
81370034214SMatthew G. Knepley {
81470034214SMatthew G. Knepley   PetscSF        fieldSF;
81570034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
81670034214SMatthew G. Knepley   PetscScalar   *originalValues, *newValues;
81770034214SMatthew G. Knepley 
81870034214SMatthew G. Knepley   PetscFunctionBegin;
8199566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0));
8209566063dSJacob Faibussowitsch   PetscCall(PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection));
82170034214SMatthew G. Knepley 
8229566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(newSection, &fieldSize));
8239566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(newVec, fieldSize, PETSC_DETERMINE));
8249566063dSJacob Faibussowitsch   PetscCall(VecSetType(newVec,dm->vectype));
82570034214SMatthew G. Knepley 
8269566063dSJacob Faibussowitsch   PetscCall(VecGetArray(originalVec, &originalValues));
8279566063dSJacob Faibussowitsch   PetscCall(VecGetArray(newVec, &newValues));
8289566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF));
8299566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
8309566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues,MPI_REPLACE));
8319566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues,MPI_REPLACE));
8329566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&fieldSF));
8339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(newVec, &newValues));
8349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(originalVec, &originalValues));
8359566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0));
83670034214SMatthew G. Knepley   PetscFunctionReturn(0);
83770034214SMatthew G. Knepley }
83870034214SMatthew G. Knepley 
83970034214SMatthew G. Knepley /*@
84070034214SMatthew G. Knepley   DMPlexDistributeFieldIS - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
84170034214SMatthew G. Knepley 
842d083f849SBarry Smith   Collective on dm
84370034214SMatthew G. Knepley 
84470034214SMatthew G. Knepley   Input Parameters:
84570034214SMatthew G. Knepley + dm - The DMPlex object
84670034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
84770034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
84870034214SMatthew G. Knepley - originalIS - The existing data
84970034214SMatthew G. Knepley 
85070034214SMatthew G. Knepley   Output Parameters:
85170034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout
85270034214SMatthew G. Knepley - newIS - The new data
85370034214SMatthew G. Knepley 
85470034214SMatthew G. Knepley   Level: developer
85570034214SMatthew G. Knepley 
856db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`, `DMPlexDistributeData()`
85770034214SMatthew G. Knepley @*/
85870034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeFieldIS(DM dm, PetscSF pointSF, PetscSection originalSection, IS originalIS, PetscSection newSection, IS *newIS)
85970034214SMatthew G. Knepley {
86070034214SMatthew G. Knepley   PetscSF         fieldSF;
86170034214SMatthew G. Knepley   PetscInt       *newValues, *remoteOffsets, fieldSize;
86270034214SMatthew G. Knepley   const PetscInt *originalValues;
86370034214SMatthew G. Knepley 
86470034214SMatthew G. Knepley   PetscFunctionBegin;
8659566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0));
8669566063dSJacob Faibussowitsch   PetscCall(PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection));
86770034214SMatthew G. Knepley 
8689566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(newSection, &fieldSize));
8699566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(fieldSize, &newValues));
87070034214SMatthew G. Knepley 
8719566063dSJacob Faibussowitsch   PetscCall(ISGetIndices(originalIS, &originalValues));
8729566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF));
8739566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
8749566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues,MPI_REPLACE));
8759566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues,MPI_REPLACE));
8769566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&fieldSF));
8779566063dSJacob Faibussowitsch   PetscCall(ISRestoreIndices(originalIS, &originalValues));
8789566063dSJacob Faibussowitsch   PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject) pointSF), fieldSize, newValues, PETSC_OWN_POINTER, newIS));
8799566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0));
88070034214SMatthew G. Knepley   PetscFunctionReturn(0);
88170034214SMatthew G. Knepley }
88270034214SMatthew G. Knepley 
88370034214SMatthew G. Knepley /*@
88470034214SMatthew G. Knepley   DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
88570034214SMatthew G. Knepley 
886d083f849SBarry Smith   Collective on dm
88770034214SMatthew G. Knepley 
88870034214SMatthew G. Knepley   Input Parameters:
88970034214SMatthew G. Knepley + dm - The DMPlex object
89070034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
89170034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
89270034214SMatthew G. Knepley . datatype - The type of data
89370034214SMatthew G. Knepley - originalData - The existing data
89470034214SMatthew G. Knepley 
89570034214SMatthew G. Knepley   Output Parameters:
89670034214SMatthew G. Knepley + newSection - The PetscSection describing the new data layout
89770034214SMatthew G. Knepley - newData - The new data
89870034214SMatthew G. Knepley 
89970034214SMatthew G. Knepley   Level: developer
90070034214SMatthew G. Knepley 
901db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeField()`
90270034214SMatthew G. Knepley @*/
90370034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData)
90470034214SMatthew G. Knepley {
90570034214SMatthew G. Knepley   PetscSF        fieldSF;
90670034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
90770034214SMatthew G. Knepley   PetscMPIInt    dataSize;
90870034214SMatthew G. Knepley 
90970034214SMatthew G. Knepley   PetscFunctionBegin;
9109566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0));
9119566063dSJacob Faibussowitsch   PetscCall(PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection));
91270034214SMatthew G. Knepley 
9139566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(newSection, &fieldSize));
9149566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Type_size(datatype, &dataSize));
9159566063dSJacob Faibussowitsch   PetscCall(PetscMalloc(fieldSize * dataSize, newData));
91670034214SMatthew G. Knepley 
9179566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF));
9189566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
9199566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(fieldSF, datatype, originalData, *newData,MPI_REPLACE));
9209566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(fieldSF, datatype, originalData, *newData,MPI_REPLACE));
9219566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&fieldSF));
9229566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0));
92370034214SMatthew G. Knepley   PetscFunctionReturn(0);
92470034214SMatthew G. Knepley }
92570034214SMatthew G. Knepley 
92624cc2ca5SMatthew G. Knepley static PetscErrorCode DMPlexDistributeCones(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping original, ISLocalToGlobalMapping renumbering, DM dmParallel)
927cc71bff1SMichael Lange {
928cc71bff1SMichael Lange   DM_Plex               *pmesh = (DM_Plex*) (dmParallel)->data;
929cc71bff1SMichael Lange   MPI_Comm               comm;
930cc71bff1SMichael Lange   PetscSF                coneSF;
931cc71bff1SMichael Lange   PetscSection           originalConeSection, newConeSection;
932ac04eaf7SToby Isaac   PetscInt              *remoteOffsets, *cones, *globCones, *newCones, newConesSize;
933cc71bff1SMichael Lange   PetscBool              flg;
934cc71bff1SMichael Lange 
935cc71bff1SMichael Lange   PetscFunctionBegin;
936cc71bff1SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9370c86c063SLisandro Dalcin   PetscValidHeaderSpecific(dmParallel, DM_CLASSID, 5);
9389566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0));
939cc71bff1SMichael Lange   /* Distribute cone section */
9409566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
9419566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSection(dm, &originalConeSection));
9429566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeSection(dmParallel, &newConeSection));
9439566063dSJacob Faibussowitsch   PetscCall(PetscSFDistributeSection(migrationSF, originalConeSection, &remoteOffsets, newConeSection));
9449566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dmParallel));
945cc71bff1SMichael Lange   /* Communicate and renumber cones */
9469566063dSJacob Faibussowitsch   PetscCall(PetscSFCreateSectionSF(migrationSF, originalConeSection, remoteOffsets, newConeSection, &coneSF));
9479566063dSJacob Faibussowitsch   PetscCall(PetscFree(remoteOffsets));
9489566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCones(dm, &cones));
949ac04eaf7SToby Isaac   if (original) {
950ac04eaf7SToby Isaac     PetscInt numCones;
951ac04eaf7SToby Isaac 
9529566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetStorageSize(originalConeSection,&numCones));
9539566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numCones,&globCones));
9549566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingApplyBlock(original, numCones, cones, globCones));
955367003a6SStefano Zampini   } else {
956ac04eaf7SToby Isaac     globCones = cones;
957ac04eaf7SToby Isaac   }
9589566063dSJacob Faibussowitsch   PetscCall(DMPlexGetCones(dmParallel, &newCones));
9599566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(coneSF, MPIU_INT, globCones, newCones,MPI_REPLACE));
9609566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(coneSF, MPIU_INT, globCones, newCones,MPI_REPLACE));
961ac04eaf7SToby Isaac   if (original) {
9629566063dSJacob Faibussowitsch     PetscCall(PetscFree(globCones));
963ac04eaf7SToby Isaac   }
9649566063dSJacob Faibussowitsch   PetscCall(PetscSectionGetStorageSize(newConeSection, &newConesSize));
9659566063dSJacob Faibussowitsch   PetscCall(ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones));
96676bd3646SJed Brown   if (PetscDefined(USE_DEBUG)) {
9673533c52bSMatthew G. Knepley     PetscInt  p;
9683533c52bSMatthew G. Knepley     PetscBool valid = PETSC_TRUE;
9693533c52bSMatthew G. Knepley     for (p = 0; p < newConesSize; ++p) {
97063a3b9bcSJacob Faibussowitsch       if (newCones[p] < 0) {valid = PETSC_FALSE; PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] Point %" PetscInt_FMT " not in overlap SF\n", PetscGlobalRank,p));}
9713533c52bSMatthew G. Knepley     }
97228b400f6SJacob Faibussowitsch     PetscCheck(valid,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map");
9733533c52bSMatthew G. Knepley   }
9749566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-cones_view", &flg));
975cc71bff1SMichael Lange   if (flg) {
9769566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Serial Cone Section:\n"));
9779566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_(comm)));
9789566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Parallel Cone Section:\n"));
9799566063dSJacob Faibussowitsch     PetscCall(PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_(comm)));
9809566063dSJacob Faibussowitsch     PetscCall(PetscSFView(coneSF, NULL));
981cc71bff1SMichael Lange   }
9829566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeOrientations(dm, &cones));
9839566063dSJacob Faibussowitsch   PetscCall(DMPlexGetConeOrientations(dmParallel, &newCones));
9849566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones,MPI_REPLACE));
9859566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones,MPI_REPLACE));
9869566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&coneSF));
9879566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0));
988eaf898f9SPatrick Sanan   /* Create supports and stratify DMPlex */
989cc71bff1SMichael Lange   {
990cc71bff1SMichael Lange     PetscInt pStart, pEnd;
991cc71bff1SMichael Lange 
9929566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd));
9939566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetChart(pmesh->supportSection, pStart, pEnd));
994cc71bff1SMichael Lange   }
9959566063dSJacob Faibussowitsch   PetscCall(DMPlexSymmetrize(dmParallel));
9969566063dSJacob Faibussowitsch   PetscCall(DMPlexStratify(dmParallel));
9971cf84007SMatthew G. Knepley   {
9981cf84007SMatthew G. Knepley     PetscBool useCone, useClosure, useAnchors;
9991cf84007SMatthew G. Knepley 
10009566063dSJacob Faibussowitsch     PetscCall(DMGetBasicAdjacency(dm, &useCone, &useClosure));
10019566063dSJacob Faibussowitsch     PetscCall(DMSetBasicAdjacency(dmParallel, useCone, useClosure));
10029566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacencyUseAnchors(dm, &useAnchors));
10039566063dSJacob Faibussowitsch     PetscCall(DMPlexSetAdjacencyUseAnchors(dmParallel, useAnchors));
10041cf84007SMatthew G. Knepley   }
1005cc71bff1SMichael Lange   PetscFunctionReturn(0);
1006cc71bff1SMichael Lange }
1007cc71bff1SMichael Lange 
100824cc2ca5SMatthew G. Knepley static PetscErrorCode DMPlexDistributeCoordinates(DM dm, PetscSF migrationSF, DM dmParallel)
10090df0e737SMichael Lange {
10100df0e737SMichael Lange   MPI_Comm         comm;
10119318fe57SMatthew G. Knepley   DM               cdm, cdmParallel;
10120df0e737SMichael Lange   PetscSection     originalCoordSection, newCoordSection;
10130df0e737SMichael Lange   Vec              originalCoordinates, newCoordinates;
10140df0e737SMichael Lange   PetscInt         bs;
101590b157c4SStefano Zampini   PetscBool        isper;
10160df0e737SMichael Lange   const char      *name;
10170df0e737SMichael Lange   const PetscReal *maxCell, *L;
10185dc8c3f7SMatthew G. Knepley   const DMBoundaryType *bd;
10190df0e737SMichael Lange 
10200df0e737SMichael Lange   PetscFunctionBegin;
10210df0e737SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10220c86c063SLisandro Dalcin   PetscValidHeaderSpecific(dmParallel, DM_CLASSID, 3);
10230df0e737SMichael Lange 
10249566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
10259566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dm, &originalCoordSection));
10269566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateSection(dmParallel, &newCoordSection));
10279566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(dm, &originalCoordinates));
10280df0e737SMichael Lange   if (originalCoordinates) {
10299566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &newCoordinates));
10309566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject) originalCoordinates, &name));
10319566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) newCoordinates, name));
10320df0e737SMichael Lange 
10339566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeField(dm, migrationSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates));
10349566063dSJacob Faibussowitsch     PetscCall(DMSetCoordinatesLocal(dmParallel, newCoordinates));
10359566063dSJacob Faibussowitsch     PetscCall(VecGetBlockSize(originalCoordinates, &bs));
10369566063dSJacob Faibussowitsch     PetscCall(VecSetBlockSize(newCoordinates, bs));
10379566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&newCoordinates));
10380df0e737SMichael Lange   }
10399566063dSJacob Faibussowitsch   PetscCall(DMGetPeriodicity(dm, &isper, &maxCell, &L, &bd));
10409566063dSJacob Faibussowitsch   PetscCall(DMSetPeriodicity(dmParallel, isper, maxCell, L, bd));
10419566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dm, &cdm));
10429566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(dmParallel, &cdmParallel));
10439566063dSJacob Faibussowitsch   PetscCall(DMCopyDisc(cdm, cdmParallel));
10440df0e737SMichael Lange   PetscFunctionReturn(0);
10450df0e737SMichael Lange }
10460df0e737SMichael Lange 
104724cc2ca5SMatthew G. Knepley static PetscErrorCode DMPlexDistributeLabels(DM dm, PetscSF migrationSF, DM dmParallel)
10480df0e737SMichael Lange {
1049df0420ecSMatthew G. Knepley   DM_Plex         *mesh = (DM_Plex*) dm->data;
10500df0e737SMichael Lange   MPI_Comm         comm;
10517980c9d4SMatthew G. Knepley   DMLabel          depthLabel;
10520df0e737SMichael Lange   PetscMPIInt      rank;
10537980c9d4SMatthew G. Knepley   PetscInt         depth, d, numLabels, numLocalLabels, l;
1054df0420ecSMatthew G. Knepley   PetscBool        hasLabels = PETSC_FALSE, lsendDepth, sendDepth;
1055df0420ecSMatthew G. Knepley   PetscObjectState depthState = -1;
10560df0e737SMichael Lange 
10570df0e737SMichael Lange   PetscFunctionBegin;
10580df0e737SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10590c86c063SLisandro Dalcin   PetscValidHeaderSpecific(dmParallel, DM_CLASSID, 3);
10600c86c063SLisandro Dalcin 
10619566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0));
10629566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
10639566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
10640df0e737SMichael Lange 
1065df0420ecSMatthew G. Knepley   /* If the user has changed the depth label, communicate it instead */
10669566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepth(dm, &depth));
10679566063dSJacob Faibussowitsch   PetscCall(DMPlexGetDepthLabel(dm, &depthLabel));
10689566063dSJacob Faibussowitsch   if (depthLabel) PetscCall(PetscObjectStateGet((PetscObject) depthLabel, &depthState));
1069df0420ecSMatthew G. Knepley   lsendDepth = mesh->depthState != depthState ? PETSC_TRUE : PETSC_FALSE;
10701c2dc1cbSBarry Smith   PetscCall(MPIU_Allreduce(&lsendDepth, &sendDepth, 1, MPIU_BOOL, MPI_LOR, comm));
1071df0420ecSMatthew G. Knepley   if (sendDepth) {
10729566063dSJacob Faibussowitsch     PetscCall(DMPlexGetDepthLabel(dmParallel, &dmParallel->depthLabel));
10739566063dSJacob Faibussowitsch     PetscCall(DMRemoveLabelBySelf(dmParallel, &dmParallel->depthLabel, PETSC_FALSE));
1074df0420ecSMatthew G. Knepley   }
1075d995df53SMatthew G. Knepley   /* Everyone must have either the same number of labels, or none */
10769566063dSJacob Faibussowitsch   PetscCall(DMGetNumLabels(dm, &numLocalLabels));
1077d995df53SMatthew G. Knepley   numLabels = numLocalLabels;
10789566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm));
1079627847f0SMatthew G. Knepley   if (numLabels == numLocalLabels) hasLabels = PETSC_TRUE;
10805d80c0bfSVaclav Hapla   for (l = 0; l < numLabels; ++l) {
1081bdd2d751SMichael Lange     DMLabel     label = NULL, labelNew = NULL;
108283e10cb3SLisandro Dalcin     PetscBool   isDepth, lisOutput = PETSC_TRUE, isOutput;
1083d67d17b1SMatthew G. Knepley     const char *name = NULL;
10840df0e737SMichael Lange 
1085d67d17b1SMatthew G. Knepley     if (hasLabels) {
10869566063dSJacob Faibussowitsch       PetscCall(DMGetLabelByNum(dm, l, &label));
10870df0e737SMichael Lange       /* Skip "depth" because it is recreated */
10889566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject) label, &name));
10899566063dSJacob Faibussowitsch       PetscCall(PetscStrcmp(name, "depth", &isDepth));
1090d67d17b1SMatthew G. Knepley     }
10919566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&isDepth, 1, MPIU_BOOL, 0, comm));
109283e10cb3SLisandro Dalcin     if (isDepth && !sendDepth) continue;
10939566063dSJacob Faibussowitsch     PetscCall(DMLabelDistribute(label, migrationSF, &labelNew));
109483e10cb3SLisandro Dalcin     if (isDepth) {
10957980c9d4SMatthew G. Knepley       /* Put in any missing strata which can occur if users are managing the depth label themselves */
10967980c9d4SMatthew G. Knepley       PetscInt gdepth;
10977980c9d4SMatthew G. Knepley 
10981c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(&depth, &gdepth, 1, MPIU_INT, MPI_MAX, comm));
109963a3b9bcSJacob Faibussowitsch       PetscCheck(!(depth >= 0) || !(gdepth != depth),PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent Plex depth %" PetscInt_FMT " != %" PetscInt_FMT, depth, gdepth);
11007980c9d4SMatthew G. Knepley       for (d = 0; d <= gdepth; ++d) {
11017980c9d4SMatthew G. Knepley         PetscBool has;
11027980c9d4SMatthew G. Knepley 
11039566063dSJacob Faibussowitsch         PetscCall(DMLabelHasStratum(labelNew, d, &has));
11049566063dSJacob Faibussowitsch         if (!has) PetscCall(DMLabelAddStratum(labelNew, d));
11057980c9d4SMatthew G. Knepley       }
11067980c9d4SMatthew G. Knepley     }
11079566063dSJacob Faibussowitsch     PetscCall(DMAddLabel(dmParallel, labelNew));
110883e10cb3SLisandro Dalcin     /* Put the output flag in the new label */
11099566063dSJacob Faibussowitsch     if (hasLabels) PetscCall(DMGetLabelOutput(dm, name, &lisOutput));
11101c2dc1cbSBarry Smith     PetscCall(MPIU_Allreduce(&lisOutput, &isOutput, 1, MPIU_BOOL, MPI_LAND, comm));
11119566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName((PetscObject) labelNew, &name));
11129566063dSJacob Faibussowitsch     PetscCall(DMSetLabelOutput(dmParallel, name, isOutput));
11139566063dSJacob Faibussowitsch     PetscCall(DMLabelDestroy(&labelNew));
11140df0e737SMichael Lange   }
11159566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0));
11160df0e737SMichael Lange   PetscFunctionReturn(0);
11170df0e737SMichael Lange }
11180df0e737SMichael Lange 
111924cc2ca5SMatthew G. Knepley static PetscErrorCode DMPlexDistributeSetupTree(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping original, ISLocalToGlobalMapping renumbering, DM dmParallel)
1120a6f36705SMichael Lange {
112115078cd4SMichael Lange   DM_Plex        *mesh  = (DM_Plex*) dm->data;
112215078cd4SMichael Lange   DM_Plex        *pmesh = (DM_Plex*) (dmParallel)->data;
1123a6f36705SMichael Lange   MPI_Comm        comm;
1124a6f36705SMichael Lange   DM              refTree;
1125a6f36705SMichael Lange   PetscSection    origParentSection, newParentSection;
1126a6f36705SMichael Lange   PetscInt        *origParents, *origChildIDs;
1127a6f36705SMichael Lange   PetscBool       flg;
1128a6f36705SMichael Lange 
1129a6f36705SMichael Lange   PetscFunctionBegin;
1130a6f36705SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11310c86c063SLisandro Dalcin   PetscValidHeaderSpecific(dmParallel, DM_CLASSID, 5);
11329566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
1133a6f36705SMichael Lange 
1134a6f36705SMichael Lange   /* Set up tree */
11359566063dSJacob Faibussowitsch   PetscCall(DMPlexGetReferenceTree(dm,&refTree));
11369566063dSJacob Faibussowitsch   PetscCall(DMPlexSetReferenceTree(dmParallel,refTree));
11379566063dSJacob Faibussowitsch   PetscCall(DMPlexGetTree(dm,&origParentSection,&origParents,&origChildIDs,NULL,NULL));
1138a6f36705SMichael Lange   if (origParentSection) {
1139a6f36705SMichael Lange     PetscInt        pStart, pEnd;
114008633170SToby Isaac     PetscInt        *newParents, *newChildIDs, *globParents;
1141a6f36705SMichael Lange     PetscInt        *remoteOffsetsParents, newParentSize;
1142a6f36705SMichael Lange     PetscSF         parentSF;
1143a6f36705SMichael Lange 
11449566063dSJacob Faibussowitsch     PetscCall(DMPlexGetChart(dmParallel, &pStart, &pEnd));
11459566063dSJacob Faibussowitsch     PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)dmParallel),&newParentSection));
11469566063dSJacob Faibussowitsch     PetscCall(PetscSectionSetChart(newParentSection,pStart,pEnd));
11479566063dSJacob Faibussowitsch     PetscCall(PetscSFDistributeSection(migrationSF, origParentSection, &remoteOffsetsParents, newParentSection));
11489566063dSJacob Faibussowitsch     PetscCall(PetscSFCreateSectionSF(migrationSF, origParentSection, remoteOffsetsParents, newParentSection, &parentSF));
11499566063dSJacob Faibussowitsch     PetscCall(PetscFree(remoteOffsetsParents));
11509566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetStorageSize(newParentSection,&newParentSize));
11519566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(newParentSize,&newParents,newParentSize,&newChildIDs));
115208633170SToby Isaac     if (original) {
115308633170SToby Isaac       PetscInt numParents;
115408633170SToby Isaac 
11559566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetStorageSize(origParentSection,&numParents));
11569566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(numParents,&globParents));
11579566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingApplyBlock(original, numParents, origParents, globParents));
115808633170SToby Isaac     }
115908633170SToby Isaac     else {
116008633170SToby Isaac       globParents = origParents;
116108633170SToby Isaac     }
11629566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(parentSF, MPIU_INT, globParents, newParents,MPI_REPLACE));
11639566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(parentSF, MPIU_INT, globParents, newParents,MPI_REPLACE));
116408633170SToby Isaac     if (original) {
11659566063dSJacob Faibussowitsch       PetscCall(PetscFree(globParents));
116608633170SToby Isaac     }
11679566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(parentSF, MPIU_INT, origChildIDs, newChildIDs,MPI_REPLACE));
11689566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(parentSF, MPIU_INT, origChildIDs, newChildIDs,MPI_REPLACE));
11699566063dSJacob Faibussowitsch     PetscCall(ISGlobalToLocalMappingApplyBlock(renumbering,IS_GTOLM_MASK, newParentSize, newParents, NULL, newParents));
117076bd3646SJed Brown     if (PetscDefined(USE_DEBUG)) {
11714a54e071SToby Isaac       PetscInt  p;
11724a54e071SToby Isaac       PetscBool valid = PETSC_TRUE;
11734a54e071SToby Isaac       for (p = 0; p < newParentSize; ++p) {
117463a3b9bcSJacob Faibussowitsch         if (newParents[p] < 0) {valid = PETSC_FALSE; PetscCall(PetscPrintf(PETSC_COMM_SELF, "Point %" PetscInt_FMT " not in overlap SF\n", p));}
11754a54e071SToby Isaac       }
117628b400f6SJacob Faibussowitsch       PetscCheck(valid,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map");
11774a54e071SToby Isaac     }
11789566063dSJacob Faibussowitsch     PetscCall(PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-parents_view", &flg));
1179a6f36705SMichael Lange     if (flg) {
11809566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Serial Parent Section: \n"));
11819566063dSJacob Faibussowitsch       PetscCall(PetscSectionView(origParentSection, PETSC_VIEWER_STDOUT_(comm)));
11829566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Parallel Parent Section: \n"));
11839566063dSJacob Faibussowitsch       PetscCall(PetscSectionView(newParentSection, PETSC_VIEWER_STDOUT_(comm)));
11849566063dSJacob Faibussowitsch       PetscCall(PetscSFView(parentSF, NULL));
1185a6f36705SMichael Lange     }
11869566063dSJacob Faibussowitsch     PetscCall(DMPlexSetTree(dmParallel,newParentSection,newParents,newChildIDs));
11879566063dSJacob Faibussowitsch     PetscCall(PetscSectionDestroy(&newParentSection));
11889566063dSJacob Faibussowitsch     PetscCall(PetscFree2(newParents,newChildIDs));
11899566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&parentSF));
1190a6f36705SMichael Lange   }
119115078cd4SMichael Lange   pmesh->useAnchors = mesh->useAnchors;
1192a6f36705SMichael Lange   PetscFunctionReturn(0);
1193a6f36705SMichael Lange }
11940df0e737SMichael Lange 
119524cc2ca5SMatthew G. Knepley PETSC_UNUSED static PetscErrorCode DMPlexDistributeSF(DM dm, PetscSF migrationSF, DM dmParallel)
11964eca1733SMichael Lange {
11979852e123SBarry Smith   PetscMPIInt            rank, size;
11984eca1733SMichael Lange   MPI_Comm               comm;
11994eca1733SMichael Lange 
12004eca1733SMichael Lange   PetscFunctionBegin;
12014eca1733SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12020c86c063SLisandro Dalcin   PetscValidHeaderSpecific(dmParallel, DM_CLASSID, 3);
12034eca1733SMichael Lange 
12044eca1733SMichael Lange   /* Create point SF for parallel mesh */
12059566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0));
12069566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
12079566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
12089566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
12094eca1733SMichael Lange   {
12104eca1733SMichael Lange     const PetscInt *leaves;
12114eca1733SMichael Lange     PetscSFNode    *remotePoints, *rowners, *lowners;
12124eca1733SMichael Lange     PetscInt        numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints;
12134eca1733SMichael Lange     PetscInt        pStart, pEnd;
12144eca1733SMichael Lange 
12159566063dSJacob Faibussowitsch     PetscCall(DMPlexGetChart(dmParallel, &pStart, &pEnd));
12169566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(migrationSF, &numRoots, &numLeaves, &leaves, NULL));
12179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(numRoots,&rowners,numLeaves,&lowners));
12184eca1733SMichael Lange     for (p=0; p<numRoots; p++) {
12194eca1733SMichael Lange       rowners[p].rank  = -1;
12204eca1733SMichael Lange       rowners[p].index = -1;
12214eca1733SMichael Lange     }
12229566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners,MPI_REPLACE));
12239566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners,MPI_REPLACE));
12244eca1733SMichael Lange     for (p = 0; p < numLeaves; ++p) {
12254eca1733SMichael Lange       if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */
12264eca1733SMichael Lange         lowners[p].rank  = rank;
12274eca1733SMichael Lange         lowners[p].index = leaves ? leaves[p] : p;
12284eca1733SMichael Lange       } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */
12294eca1733SMichael Lange         lowners[p].rank  = -2;
12304eca1733SMichael Lange         lowners[p].index = -2;
12314eca1733SMichael Lange       }
12324eca1733SMichael Lange     }
12334eca1733SMichael Lange     for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */
12344eca1733SMichael Lange       rowners[p].rank  = -3;
12354eca1733SMichael Lange       rowners[p].index = -3;
12364eca1733SMichael Lange     }
12379566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC));
12389566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC));
12399566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners,MPI_REPLACE));
12409566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners,MPI_REPLACE));
12414eca1733SMichael Lange     for (p = 0; p < numLeaves; ++p) {
12421dca8a05SBarry Smith       PetscCheck(lowners[p].rank >= 0 && lowners[p].index >= 0,PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed");
12434eca1733SMichael Lange       if (lowners[p].rank != rank) ++numGhostPoints;
12444eca1733SMichael Lange     }
12459566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numGhostPoints, &ghostPoints));
12469566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(numGhostPoints, &remotePoints));
12474eca1733SMichael Lange     for (p = 0, gp = 0; p < numLeaves; ++p) {
12484eca1733SMichael Lange       if (lowners[p].rank != rank) {
12494eca1733SMichael Lange         ghostPoints[gp]        = leaves ? leaves[p] : p;
12504eca1733SMichael Lange         remotePoints[gp].rank  = lowners[p].rank;
12514eca1733SMichael Lange         remotePoints[gp].index = lowners[p].index;
12524eca1733SMichael Lange         ++gp;
12534eca1733SMichael Lange       }
12544eca1733SMichael Lange     }
12559566063dSJacob Faibussowitsch     PetscCall(PetscFree2(rowners,lowners));
12569566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph((dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER));
12579566063dSJacob Faibussowitsch     PetscCall(PetscSFSetFromOptions((dmParallel)->sf));
12584eca1733SMichael Lange   }
12591cf84007SMatthew G. Knepley   {
12601cf84007SMatthew G. Knepley     PetscBool useCone, useClosure, useAnchors;
12611cf84007SMatthew G. Knepley 
12629566063dSJacob Faibussowitsch     PetscCall(DMGetBasicAdjacency(dm, &useCone, &useClosure));
12639566063dSJacob Faibussowitsch     PetscCall(DMSetBasicAdjacency(dmParallel, useCone, useClosure));
12649566063dSJacob Faibussowitsch     PetscCall(DMPlexGetAdjacencyUseAnchors(dm, &useAnchors));
12659566063dSJacob Faibussowitsch     PetscCall(DMPlexSetAdjacencyUseAnchors(dmParallel, useAnchors));
12661cf84007SMatthew G. Knepley   }
12679566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0));
12684eca1733SMichael Lange   PetscFunctionReturn(0);
12694eca1733SMichael Lange }
12704eca1733SMichael Lange 
127198ba2d7fSLawrence Mitchell /*@
127298ba2d7fSLawrence Mitchell   DMPlexSetPartitionBalance - Should distribution of the DM attempt to balance the shared point partition?
127398ba2d7fSLawrence Mitchell 
127498ba2d7fSLawrence Mitchell   Input Parameters:
127598ba2d7fSLawrence Mitchell + dm - The DMPlex object
127698ba2d7fSLawrence Mitchell - flg - Balance the partition?
127798ba2d7fSLawrence Mitchell 
127898ba2d7fSLawrence Mitchell   Level: intermediate
127998ba2d7fSLawrence Mitchell 
1280db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexGetPartitionBalance()`
128198ba2d7fSLawrence Mitchell @*/
128298ba2d7fSLawrence Mitchell PetscErrorCode DMPlexSetPartitionBalance(DM dm, PetscBool flg)
128398ba2d7fSLawrence Mitchell {
128498ba2d7fSLawrence Mitchell   DM_Plex *mesh = (DM_Plex *)dm->data;
128598ba2d7fSLawrence Mitchell 
128698ba2d7fSLawrence Mitchell   PetscFunctionBegin;
128798ba2d7fSLawrence Mitchell   mesh->partitionBalance = flg;
128898ba2d7fSLawrence Mitchell   PetscFunctionReturn(0);
128998ba2d7fSLawrence Mitchell }
129098ba2d7fSLawrence Mitchell 
129198ba2d7fSLawrence Mitchell /*@
129298ba2d7fSLawrence Mitchell   DMPlexGetPartitionBalance - Does distribution of the DM attempt to balance the shared point partition?
129398ba2d7fSLawrence Mitchell 
129498ba2d7fSLawrence Mitchell   Input Parameter:
1295a2b725a8SWilliam Gropp . dm - The DMPlex object
129698ba2d7fSLawrence Mitchell 
129798ba2d7fSLawrence Mitchell   Output Parameter:
1298a2b725a8SWilliam Gropp . flg - Balance the partition?
129998ba2d7fSLawrence Mitchell 
130098ba2d7fSLawrence Mitchell   Level: intermediate
130198ba2d7fSLawrence Mitchell 
1302db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexSetPartitionBalance()`
130398ba2d7fSLawrence Mitchell @*/
130498ba2d7fSLawrence Mitchell PetscErrorCode DMPlexGetPartitionBalance(DM dm, PetscBool *flg)
130598ba2d7fSLawrence Mitchell {
130698ba2d7fSLawrence Mitchell   DM_Plex *mesh = (DM_Plex *)dm->data;
130798ba2d7fSLawrence Mitchell 
130898ba2d7fSLawrence Mitchell   PetscFunctionBegin;
130998ba2d7fSLawrence Mitchell   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1310534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg, 2);
131198ba2d7fSLawrence Mitchell   *flg = mesh->partitionBalance;
131298ba2d7fSLawrence Mitchell   PetscFunctionReturn(0);
131398ba2d7fSLawrence Mitchell }
131498ba2d7fSLawrence Mitchell 
1315fc02256fSLawrence Mitchell typedef struct {
1316fc02256fSLawrence Mitchell   PetscInt vote, rank, index;
1317fc02256fSLawrence Mitchell } Petsc3Int;
1318fc02256fSLawrence Mitchell 
1319fc02256fSLawrence Mitchell /* MaxLoc, but carry a third piece of information around */
13202eb0eadbSSatish Balay static void MPIAPI MaxLocCarry(void *in_, void *inout_, PetscMPIInt *len_, MPI_Datatype *dtype)
1321fc02256fSLawrence Mitchell {
1322fc02256fSLawrence Mitchell   Petsc3Int *a = (Petsc3Int *)inout_;
1323fc02256fSLawrence Mitchell   Petsc3Int *b = (Petsc3Int *)in_;
1324fc02256fSLawrence Mitchell   PetscInt i, len = *len_;
1325fc02256fSLawrence Mitchell   for (i = 0; i < len; i++) {
1326fc02256fSLawrence Mitchell     if (a[i].vote < b[i].vote) {
1327fc02256fSLawrence Mitchell       a[i].vote = b[i].vote;
1328fc02256fSLawrence Mitchell       a[i].rank = b[i].rank;
1329fc02256fSLawrence Mitchell       a[i].index = b[i].index;
1330fc02256fSLawrence Mitchell     } else if (a[i].vote <= b[i].vote) {
1331fc02256fSLawrence Mitchell       if (a[i].rank >= b[i].rank) {
1332fc02256fSLawrence Mitchell         a[i].rank = b[i].rank;
1333fc02256fSLawrence Mitchell         a[i].index = b[i].index;
1334fc02256fSLawrence Mitchell       }
1335fc02256fSLawrence Mitchell     }
1336fc02256fSLawrence Mitchell   }
1337fc02256fSLawrence Mitchell }
1338fc02256fSLawrence Mitchell 
1339f5bf2dbfSMichael Lange /*@C
1340a8c5bd36SStefano Zampini   DMPlexCreatePointSF - Build a point SF from an SF describing a point migration
1341f5bf2dbfSMichael Lange 
1342064ec1f3Sprj-   Input Parameters:
1343f5bf2dbfSMichael Lange + dm          - The source DMPlex object
1344f5bf2dbfSMichael Lange . migrationSF - The star forest that describes the parallel point remapping
1345d8d19677SJose E. Roman - ownership   - Flag causing a vote to determine point ownership
1346f5bf2dbfSMichael Lange 
1347f5bf2dbfSMichael Lange   Output Parameter:
1348d8d19677SJose E. Roman . pointSF     - The star forest describing the point overlap in the remapped DM
1349f5bf2dbfSMichael Lange 
13503618482eSVaclav Hapla   Notes:
13513618482eSVaclav Hapla   Output pointSF is guaranteed to have the array of local indices (leaves) sorted.
13523618482eSVaclav Hapla 
1353f5bf2dbfSMichael Lange   Level: developer
1354f5bf2dbfSMichael Lange 
1355db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeOverlap()`
1356f5bf2dbfSMichael Lange @*/
1357f5bf2dbfSMichael Lange PetscErrorCode DMPlexCreatePointSF(DM dm, PetscSF migrationSF, PetscBool ownership, PetscSF *pointSF)
1358f5bf2dbfSMichael Lange {
135923193802SMatthew G. Knepley   PetscMPIInt        rank, size;
13601627f6ccSMichael Lange   PetscInt           p, nroots, nleaves, idx, npointLeaves;
1361f5bf2dbfSMichael Lange   PetscInt          *pointLocal;
1362f5bf2dbfSMichael Lange   const PetscInt    *leaves;
1363f5bf2dbfSMichael Lange   const PetscSFNode *roots;
1364f5bf2dbfSMichael Lange   PetscSFNode       *rootNodes, *leafNodes, *pointRemote;
136523193802SMatthew G. Knepley   Vec                shifts;
1366cae3e4f3SLawrence Mitchell   const PetscInt     numShifts = 13759;
136723193802SMatthew G. Knepley   const PetscScalar *shift = NULL;
136823193802SMatthew G. Knepley   const PetscBool    shiftDebug = PETSC_FALSE;
136998ba2d7fSLawrence Mitchell   PetscBool          balance;
1370f5bf2dbfSMichael Lange 
1371f5bf2dbfSMichael Lange   PetscFunctionBegin;
1372f5bf2dbfSMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
13739566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank));
13749566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) dm), &size));
13759566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_CreatePointSF,dm,0,0,0));
1376f5bf2dbfSMichael Lange 
13779566063dSJacob Faibussowitsch   PetscCall(DMPlexGetPartitionBalance(dm, &balance));
13789566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(migrationSF, &nroots, &nleaves, &leaves, &roots));
13799566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(nroots, &rootNodes, nleaves, &leafNodes));
1380f5bf2dbfSMichael Lange   if (ownership) {
1381fc02256fSLawrence Mitchell     MPI_Op       op;
1382fc02256fSLawrence Mitchell     MPI_Datatype datatype;
1383fc02256fSLawrence Mitchell     Petsc3Int   *rootVote = NULL, *leafVote = NULL;
138423193802SMatthew G. Knepley     /* If balancing, we compute a random cyclic shift of the rank for each remote point. That way, the max will evenly distribute among ranks. */
138598ba2d7fSLawrence Mitchell     if (balance) {
138623193802SMatthew G. Knepley       PetscRandom r;
138723193802SMatthew G. Knepley 
13889566063dSJacob Faibussowitsch       PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &r));
13899566063dSJacob Faibussowitsch       PetscCall(PetscRandomSetInterval(r, 0, 2467*size));
13909566063dSJacob Faibussowitsch       PetscCall(VecCreate(PETSC_COMM_SELF, &shifts));
13919566063dSJacob Faibussowitsch       PetscCall(VecSetSizes(shifts, numShifts, numShifts));
13929566063dSJacob Faibussowitsch       PetscCall(VecSetType(shifts, VECSTANDARD));
13939566063dSJacob Faibussowitsch       PetscCall(VecSetRandom(shifts, r));
13949566063dSJacob Faibussowitsch       PetscCall(PetscRandomDestroy(&r));
13959566063dSJacob Faibussowitsch       PetscCall(VecGetArrayRead(shifts, &shift));
139623193802SMatthew G. Knepley     }
139723193802SMatthew G. Knepley 
13989566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nroots, &rootVote));
13999566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &leafVote));
140023193802SMatthew G. Knepley     /* Point ownership vote: Process with highest rank owns shared points */
1401f5bf2dbfSMichael Lange     for (p = 0; p < nleaves; ++p) {
140223193802SMatthew G. Knepley       if (shiftDebug) {
140363a3b9bcSJacob Faibussowitsch         PetscCall(PetscSynchronizedPrintf(PetscObjectComm((PetscObject) dm), "[%d] Point %" PetscInt_FMT " RemotePoint %" PetscInt_FMT " Shift %" PetscInt_FMT " MyRank %" PetscInt_FMT "\n", rank, leaves ? leaves[p] : p, roots[p].index, (PetscInt) PetscRealPart(shift[roots[p].index%numShifts]), (rank + (shift ? (PetscInt) PetscRealPart(shift[roots[p].index%numShifts]) : 0))%size));
140423193802SMatthew G. Knepley       }
1405f5bf2dbfSMichael Lange       /* Either put in a bid or we know we own it */
1406fc02256fSLawrence Mitchell       leafVote[p].vote  = (rank + (shift ? (PetscInt) PetscRealPart(shift[roots[p].index%numShifts]) : 0))%size;
1407fc02256fSLawrence Mitchell       leafVote[p].rank = rank;
1408fc02256fSLawrence Mitchell       leafVote[p].index = p;
1409f5bf2dbfSMichael Lange     }
1410f5bf2dbfSMichael Lange     for (p = 0; p < nroots; p++) {
14111627f6ccSMichael Lange       /* Root must not participate in the reduction, flag so that MAXLOC does not use */
1412fc02256fSLawrence Mitchell       rootVote[p].vote  = -3;
1413fc02256fSLawrence Mitchell       rootVote[p].rank  = -3;
1414fc02256fSLawrence Mitchell       rootVote[p].index = -3;
1415f5bf2dbfSMichael Lange     }
14169566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_contiguous(3, MPIU_INT, &datatype));
14179566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_commit(&datatype));
14189566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Op_create(&MaxLocCarry, 1, &op));
14199566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceBegin(migrationSF, datatype, leafVote, rootVote, op));
14209566063dSJacob Faibussowitsch     PetscCall(PetscSFReduceEnd(migrationSF, datatype, leafVote, rootVote, op));
14219566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Op_free(&op));
14229566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Type_free(&datatype));
1423c091126eSLawrence Mitchell     for (p = 0; p < nroots; p++) {
1424fc02256fSLawrence Mitchell       rootNodes[p].rank = rootVote[p].rank;
1425fc02256fSLawrence Mitchell       rootNodes[p].index = rootVote[p].index;
1426c091126eSLawrence Mitchell     }
14279566063dSJacob Faibussowitsch     PetscCall(PetscFree(leafVote));
14289566063dSJacob Faibussowitsch     PetscCall(PetscFree(rootVote));
1429f5bf2dbfSMichael Lange   } else {
1430f5bf2dbfSMichael Lange     for (p = 0; p < nroots; p++) {
1431f5bf2dbfSMichael Lange       rootNodes[p].index = -1;
1432f5bf2dbfSMichael Lange       rootNodes[p].rank = rank;
1433fc02256fSLawrence Mitchell     }
1434f5bf2dbfSMichael Lange     for (p = 0; p < nleaves; p++) {
1435f5bf2dbfSMichael Lange       /* Write new local id into old location */
1436f5bf2dbfSMichael Lange       if (roots[p].rank == rank) {
14376462276dSToby Isaac         rootNodes[roots[p].index].index = leaves ? leaves[p] : p;
1438f5bf2dbfSMichael Lange       }
1439f5bf2dbfSMichael Lange     }
1440f5bf2dbfSMichael Lange   }
14419566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastBegin(migrationSF, MPIU_2INT, rootNodes, leafNodes,MPI_REPLACE));
14429566063dSJacob Faibussowitsch   PetscCall(PetscSFBcastEnd(migrationSF, MPIU_2INT, rootNodes, leafNodes,MPI_REPLACE));
1443f5bf2dbfSMichael Lange 
144423193802SMatthew G. Knepley   for (npointLeaves = 0, p = 0; p < nleaves; p++) {
1445b0e1264bSMatthew G. Knepley     if (leafNodes[p].rank != rank) npointLeaves++;
144623193802SMatthew G. Knepley   }
14479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npointLeaves, &pointLocal));
14489566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(npointLeaves, &pointRemote));
1449f5bf2dbfSMichael Lange   for (idx = 0, p = 0; p < nleaves; p++) {
1450b0e1264bSMatthew G. Knepley     if (leafNodes[p].rank != rank) {
14513618482eSVaclav Hapla       /* Note that pointLocal is automatically sorted as it is sublist of 0, ..., nleaves-1 */
1452f5bf2dbfSMichael Lange       pointLocal[idx] = p;
1453f5bf2dbfSMichael Lange       pointRemote[idx] = leafNodes[p];
1454f5bf2dbfSMichael Lange       idx++;
1455f5bf2dbfSMichael Lange     }
1456f5bf2dbfSMichael Lange   }
145723193802SMatthew G. Knepley   if (shift) {
14589566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(shifts, &shift));
14599566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&shifts));
146023193802SMatthew G. Knepley   }
14619566063dSJacob Faibussowitsch   if (shiftDebug) PetscCall(PetscSynchronizedFlush(PetscObjectComm((PetscObject) dm), PETSC_STDOUT));
14629566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(PetscObjectComm((PetscObject) dm), pointSF));
14639566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(*pointSF));
14649566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(*pointSF, nleaves, npointLeaves, pointLocal, PETSC_OWN_POINTER, pointRemote, PETSC_OWN_POINTER));
14659566063dSJacob Faibussowitsch   PetscCall(PetscFree2(rootNodes, leafNodes));
14669566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_CreatePointSF,dm,0,0,0));
1467f5bf2dbfSMichael Lange   PetscFunctionReturn(0);
1468f5bf2dbfSMichael Lange }
1469f5bf2dbfSMichael Lange 
147015078cd4SMichael Lange /*@C
147115078cd4SMichael Lange   DMPlexMigrate  - Migrates internal DM data over the supplied star forest
147215078cd4SMichael Lange 
1473d083f849SBarry Smith   Collective on dm
147483655b49SVáclav Hapla 
1475064ec1f3Sprj-   Input Parameters:
147615078cd4SMichael Lange + dm       - The source DMPlex object
1477d8d19677SJose E. Roman - sf       - The star forest communication context describing the migration pattern
147815078cd4SMichael Lange 
147915078cd4SMichael Lange   Output Parameter:
1480d8d19677SJose E. Roman . targetDM - The target DMPlex object
148115078cd4SMichael Lange 
1482b9f40539SMichael Lange   Level: intermediate
148315078cd4SMichael Lange 
1484db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeOverlap()`
148515078cd4SMichael Lange @*/
1486b9f40539SMichael Lange PetscErrorCode DMPlexMigrate(DM dm, PetscSF sf, DM targetDM)
148715078cd4SMichael Lange {
1488b9f40539SMichael Lange   MPI_Comm               comm;
1489cc1750acSStefano Zampini   PetscInt               dim, cdim, nroots;
1490b9f40539SMichael Lange   PetscSF                sfPoint;
149115078cd4SMichael Lange   ISLocalToGlobalMapping ltogMigration;
1492ac04eaf7SToby Isaac   ISLocalToGlobalMapping ltogOriginal = NULL;
1493b9f40539SMichael Lange   PetscBool              flg;
149415078cd4SMichael Lange 
149515078cd4SMichael Lange   PetscFunctionBegin;
149615078cd4SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14979566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_Migrate, dm, 0, 0, 0));
14989566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject) dm, &comm));
14999566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(dm, &dim));
15009566063dSJacob Faibussowitsch   PetscCall(DMSetDimension(targetDM, dim));
15019566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDim(dm, &cdim));
15029566063dSJacob Faibussowitsch   PetscCall(DMSetCoordinateDim(targetDM, cdim));
150315078cd4SMichael Lange 
1504bfb0467fSMichael Lange   /* Check for a one-to-all distribution pattern */
15059566063dSJacob Faibussowitsch   PetscCall(DMGetPointSF(dm, &sfPoint));
15069566063dSJacob Faibussowitsch   PetscCall(PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL));
1507bfb0467fSMichael Lange   if (nroots >= 0) {
1508b9f40539SMichael Lange     IS        isOriginal;
1509ac04eaf7SToby Isaac     PetscInt  n, size, nleaves;
1510ac04eaf7SToby Isaac     PetscInt  *numbering_orig, *numbering_new;
1511367003a6SStefano Zampini 
1512b9f40539SMichael Lange     /* Get the original point numbering */
15139566063dSJacob Faibussowitsch     PetscCall(DMPlexCreatePointNumbering(dm, &isOriginal));
15149566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreateIS(isOriginal, &ltogOriginal));
15159566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetSize(ltogOriginal, &size));
15169566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingGetBlockIndices(ltogOriginal, (const PetscInt**)&numbering_orig));
1517b9f40539SMichael Lange     /* Convert to positive global numbers */
1518b9f40539SMichael Lange     for (n=0; n<size; n++) {if (numbering_orig[n] < 0) numbering_orig[n] = -(numbering_orig[n]+1);}
1519b9f40539SMichael Lange     /* Derive the new local-to-global mapping from the old one */
15209566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sf, NULL, &nleaves, NULL, NULL));
15219566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &numbering_new));
15229566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sf, MPIU_INT, numbering_orig, numbering_new,MPI_REPLACE));
15239566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sf, MPIU_INT, numbering_orig, numbering_new,MPI_REPLACE));
15249566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(comm, 1, nleaves, numbering_new, PETSC_OWN_POINTER, &ltogMigration));
15259566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingRestoreIndices(ltogOriginal, (const PetscInt**)&numbering_orig));
15269566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isOriginal));
152715078cd4SMichael Lange   } else {
1528bfb0467fSMichael Lange     /* One-to-all distribution pattern: We can derive LToG from SF */
15299566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreateSF(sf, 0, &ltogMigration));
153015078cd4SMichael Lange   }
15319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-partition_view", &flg));
1532b9f40539SMichael Lange   if (flg) {
15339566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(comm, "Point renumbering for DM migration:\n"));
15349566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingView(ltogMigration, NULL));
1535b9f40539SMichael Lange   }
153615078cd4SMichael Lange   /* Migrate DM data to target DM */
15379566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeCones(dm, sf, ltogOriginal, ltogMigration, targetDM));
15389566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeLabels(dm, sf, targetDM));
15399566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeCoordinates(dm, sf, targetDM));
15409566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeSetupTree(dm, sf, ltogOriginal, ltogMigration, targetDM));
15419566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&ltogOriginal));
15429566063dSJacob Faibussowitsch   PetscCall(ISLocalToGlobalMappingDestroy(&ltogMigration));
15439566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_Migrate, dm, 0, 0, 0));
154415078cd4SMichael Lange   PetscFunctionReturn(0);
154515078cd4SMichael Lange }
154615078cd4SMichael Lange 
15473b8f15a2SMatthew G. Knepley /*@C
154870034214SMatthew G. Knepley   DMPlexDistribute - Distributes the mesh and any associated sections.
154970034214SMatthew G. Knepley 
1550d083f849SBarry Smith   Collective on dm
155170034214SMatthew G. Knepley 
1552064ec1f3Sprj-   Input Parameters:
155370034214SMatthew G. Knepley + dm  - The original DMPlex object
155470034214SMatthew G. Knepley - overlap - The overlap of partitions, 0 is the default
155570034214SMatthew G. Knepley 
1556064ec1f3Sprj-   Output Parameters:
1557f4ae5380SJed Brown + sf - The PetscSF used for point distribution, or NULL if not needed
1558f4ae5380SJed Brown - dmParallel - The distributed DMPlex object
155970034214SMatthew G. Knepley 
1560f4ae5380SJed Brown   Note: If the mesh was not distributed, the output dmParallel will be NULL.
156170034214SMatthew G. Knepley 
1562b0441da4SMatthew G. Knepley   The user can control the definition of adjacency for the mesh using DMSetAdjacency(). They should choose the combination appropriate for the function
156370034214SMatthew G. Knepley   representation on the mesh.
156470034214SMatthew G. Knepley 
156570034214SMatthew G. Knepley   Level: intermediate
156670034214SMatthew G. Knepley 
1567db781477SPatrick Sanan .seealso: `DMPlexCreate()`, `DMSetAdjacency()`, `DMPlexGetOverlap()`
156870034214SMatthew G. Knepley @*/
156980cf41d5SMatthew G. Knepley PetscErrorCode DMPlexDistribute(DM dm, PetscInt overlap, PetscSF *sf, DM *dmParallel)
157070034214SMatthew G. Knepley {
157170034214SMatthew G. Knepley   MPI_Comm               comm;
157215078cd4SMichael Lange   PetscPartitioner       partitioner;
1573f8987ae8SMichael Lange   IS                     cellPart;
1574f8987ae8SMichael Lange   PetscSection           cellPartSection;
1575cf86098cSMatthew G. Knepley   DM                     dmCoord;
1576f8987ae8SMichael Lange   DMLabel                lblPartition, lblMigration;
1577874ddda9SLisandro Dalcin   PetscSF                sfMigration, sfStratified, sfPoint;
157898ba2d7fSLawrence Mitchell   PetscBool              flg, balance;
1579874ddda9SLisandro Dalcin   PetscMPIInt            rank, size;
158070034214SMatthew G. Knepley 
158170034214SMatthew G. Knepley   PetscFunctionBegin;
158270034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1583d5c515a1SMatthew G. Knepley   PetscValidLogicalCollectiveInt(dm, overlap, 2);
1584d5c515a1SMatthew G. Knepley   if (sf) PetscValidPointer(sf,3);
1585d5c515a1SMatthew G. Knepley   PetscValidPointer(dmParallel,4);
158670034214SMatthew G. Knepley 
15870c86c063SLisandro Dalcin   if (sf) *sf = NULL;
15880c86c063SLisandro Dalcin   *dmParallel = NULL;
15899566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
15909566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
15919566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
15920c86c063SLisandro Dalcin   if (size == 1) PetscFunctionReturn(0);
159370034214SMatthew G. Knepley 
15949566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0));
159515078cd4SMichael Lange   /* Create cell partition */
15969566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_Partition,dm,0,0,0));
15979566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &cellPartSection));
15989566063dSJacob Faibussowitsch   PetscCall(DMPlexGetPartitioner(dm, &partitioner));
15999566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerDMPlexPartition(partitioner, dm, NULL, cellPartSection, &cellPart));
16009566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_PartSelf,dm,0,0,0));
1601f8987ae8SMichael Lange   {
1602f8987ae8SMichael Lange     /* Convert partition to DMLabel */
1603825f8a23SLisandro Dalcin     IS             is;
1604825f8a23SLisandro Dalcin     PetscHSetI     ht;
1605f8987ae8SMichael Lange     const PetscInt *points;
16068e330a33SStefano Zampini     PetscInt       *iranks;
16078e330a33SStefano Zampini     PetscInt       pStart, pEnd, proc, npoints, poff = 0, nranks;
1608825f8a23SLisandro Dalcin 
16099566063dSJacob Faibussowitsch     PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Point Partition", &lblPartition));
1610825f8a23SLisandro Dalcin     /* Preallocate strata */
16119566063dSJacob Faibussowitsch     PetscCall(PetscHSetICreate(&ht));
16129566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(cellPartSection, &pStart, &pEnd));
1613825f8a23SLisandro Dalcin     for (proc = pStart; proc < pEnd; proc++) {
16149566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(cellPartSection, proc, &npoints));
16159566063dSJacob Faibussowitsch       if (npoints) PetscCall(PetscHSetIAdd(ht, proc));
1616825f8a23SLisandro Dalcin     }
16179566063dSJacob Faibussowitsch     PetscCall(PetscHSetIGetSize(ht, &nranks));
16189566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nranks, &iranks));
16199566063dSJacob Faibussowitsch     PetscCall(PetscHSetIGetElems(ht, &poff, iranks));
16209566063dSJacob Faibussowitsch     PetscCall(PetscHSetIDestroy(&ht));
16219566063dSJacob Faibussowitsch     PetscCall(DMLabelAddStrata(lblPartition, nranks, iranks));
16229566063dSJacob Faibussowitsch     PetscCall(PetscFree(iranks));
1623825f8a23SLisandro Dalcin     /* Inline DMPlexPartitionLabelClosure() */
16249566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(cellPart, &points));
16259566063dSJacob Faibussowitsch     PetscCall(PetscSectionGetChart(cellPartSection, &pStart, &pEnd));
1626f8987ae8SMichael Lange     for (proc = pStart; proc < pEnd; proc++) {
16279566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetDof(cellPartSection, proc, &npoints));
1628825f8a23SLisandro Dalcin       if (!npoints) continue;
16299566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetOffset(cellPartSection, proc, &poff));
16309566063dSJacob Faibussowitsch       PetscCall(DMPlexClosurePoints_Private(dm, npoints, points+poff, &is));
16319566063dSJacob Faibussowitsch       PetscCall(DMLabelSetStratumIS(lblPartition, proc, is));
16329566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&is));
1633f8987ae8SMichael Lange     }
16349566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(cellPart, &points));
1635f8987ae8SMichael Lange   }
16369566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_PartSelf,dm,0,0,0));
16376e203dd9SStefano Zampini 
16389566063dSJacob Faibussowitsch   PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Point migration", &lblMigration));
16399566063dSJacob Faibussowitsch   PetscCall(DMPlexPartitionLabelInvert(dm, lblPartition, NULL, lblMigration));
16409566063dSJacob Faibussowitsch   PetscCall(DMPlexPartitionLabelCreateSF(dm, lblMigration, &sfMigration));
16419566063dSJacob Faibussowitsch   PetscCall(DMPlexStratifyMigrationSF(dm, sfMigration, &sfStratified));
16429566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfMigration));
164343f7d02bSMichael Lange   sfMigration = sfStratified;
16449566063dSJacob Faibussowitsch   PetscCall(PetscSFSetUp(sfMigration));
16459566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_Partition,dm,0,0,0));
16469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(((PetscObject) dm)->options,((PetscObject) dm)->prefix, "-partition_view", &flg));
164770034214SMatthew G. Knepley   if (flg) {
16489566063dSJacob Faibussowitsch     PetscCall(DMLabelView(lblPartition, PETSC_VIEWER_STDOUT_(comm)));
16499566063dSJacob Faibussowitsch     PetscCall(PetscSFView(sfMigration, PETSC_VIEWER_STDOUT_(comm)));
165070034214SMatthew G. Knepley   }
1651f8987ae8SMichael Lange 
165215078cd4SMichael Lange   /* Create non-overlapping parallel DM and migrate internal data */
16539566063dSJacob Faibussowitsch   PetscCall(DMPlexCreate(comm, dmParallel));
16549566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh"));
16559566063dSJacob Faibussowitsch   PetscCall(DMPlexMigrate(dm, sfMigration, *dmParallel));
165670034214SMatthew G. Knepley 
1657a157612eSMichael Lange   /* Build the point SF without overlap */
16589566063dSJacob Faibussowitsch   PetscCall(DMPlexGetPartitionBalance(dm, &balance));
16599566063dSJacob Faibussowitsch   PetscCall(DMPlexSetPartitionBalance(*dmParallel, balance));
16609566063dSJacob Faibussowitsch   PetscCall(DMPlexCreatePointSF(*dmParallel, sfMigration, PETSC_TRUE, &sfPoint));
16619566063dSJacob Faibussowitsch   PetscCall(DMSetPointSF(*dmParallel, sfPoint));
16629566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(*dmParallel, &dmCoord));
16639566063dSJacob Faibussowitsch   if (dmCoord) PetscCall(DMSetPointSF(dmCoord, sfPoint));
16649566063dSJacob Faibussowitsch   if (flg) PetscCall(PetscSFView(sfPoint, NULL));
166570034214SMatthew G. Knepley 
1666a157612eSMichael Lange   if (overlap > 0) {
166715078cd4SMichael Lange     DM                 dmOverlap;
166857f290daSMatthew G. Knepley     PetscInt           nroots, nleaves, noldleaves, l;
166957f290daSMatthew G. Knepley     const PetscInt    *oldLeaves;
167057f290daSMatthew G. Knepley     PetscSFNode       *newRemote, *permRemote;
167115078cd4SMichael Lange     const PetscSFNode *oldRemote;
167215078cd4SMichael Lange     PetscSF            sfOverlap, sfOverlapPoint;
1673524e35f8SStefano Zampini 
1674a157612eSMichael Lange     /* Add the partition overlap to the distributed DM */
16759566063dSJacob Faibussowitsch     PetscCall(DMPlexDistributeOverlap(*dmParallel, overlap, &sfOverlap, &dmOverlap));
16769566063dSJacob Faibussowitsch     PetscCall(DMDestroy(dmParallel));
1677a157612eSMichael Lange     *dmParallel = dmOverlap;
1678c389ea9fSToby Isaac     if (flg) {
16799566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(comm, "Overlap Migration SF:\n"));
16809566063dSJacob Faibussowitsch       PetscCall(PetscSFView(sfOverlap, NULL));
1681c389ea9fSToby Isaac     }
168243331d4aSMichael Lange 
1683f8987ae8SMichael Lange     /* Re-map the migration SF to establish the full migration pattern */
16849566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sfMigration, &nroots, &noldleaves, &oldLeaves, &oldRemote));
16859566063dSJacob Faibussowitsch     PetscCall(PetscSFGetGraph(sfOverlap, NULL, &nleaves, NULL, NULL));
16869566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nleaves, &newRemote));
168757f290daSMatthew G. Knepley     /* oldRemote: original root point mapping to original leaf point
168857f290daSMatthew G. Knepley        newRemote: original leaf point mapping to overlapped leaf point */
168957f290daSMatthew G. Knepley     if (oldLeaves) {
169073e69a6aSMatthew G. Knepley       /* After stratification, the migration remotes may not be in root (canonical) order, so we reorder using the leaf numbering */
16919566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(noldleaves, &permRemote));
169257f290daSMatthew G. Knepley       for (l = 0; l < noldleaves; ++l) permRemote[oldLeaves[l]] = oldRemote[l];
169357f290daSMatthew G. Knepley       oldRemote = permRemote;
169457f290daSMatthew G. Knepley     }
16959566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastBegin(sfOverlap, MPIU_2INT, oldRemote, newRemote,MPI_REPLACE));
16969566063dSJacob Faibussowitsch     PetscCall(PetscSFBcastEnd(sfOverlap, MPIU_2INT, oldRemote, newRemote,MPI_REPLACE));
16979566063dSJacob Faibussowitsch     if (oldLeaves) PetscCall(PetscFree(oldRemote));
16989566063dSJacob Faibussowitsch     PetscCall(PetscSFCreate(comm, &sfOverlapPoint));
16999566063dSJacob Faibussowitsch     PetscCall(PetscSFSetGraph(sfOverlapPoint, nroots, nleaves, NULL, PETSC_OWN_POINTER, newRemote, PETSC_OWN_POINTER));
17009566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfOverlap));
17019566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&sfMigration));
170215078cd4SMichael Lange     sfMigration = sfOverlapPoint;
1703c389ea9fSToby Isaac   }
1704bf5244c7SMatthew G. Knepley   /* Cleanup Partition */
17059566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&lblPartition));
17069566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&lblMigration));
17079566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&cellPartSection));
17089566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&cellPart));
170945480ffeSMatthew G. Knepley   /* Copy discretization */
17109566063dSJacob Faibussowitsch   PetscCall(DMCopyDisc(dm, *dmParallel));
171166fe0bfeSMatthew G. Knepley   /* Create sfNatural */
171266fe0bfeSMatthew G. Knepley   if (dm->useNatural) {
171366fe0bfeSMatthew G. Knepley     PetscSection section;
171466fe0bfeSMatthew G. Knepley 
17159566063dSJacob Faibussowitsch     PetscCall(DMGetLocalSection(dm, &section));
17169566063dSJacob Faibussowitsch     PetscCall(DMPlexCreateGlobalToNaturalSF(*dmParallel, section, sfMigration, &(*dmParallel)->sfNatural));
17179566063dSJacob Faibussowitsch     PetscCall(DMSetUseNatural(*dmParallel, PETSC_TRUE));
171866fe0bfeSMatthew G. Knepley   }
1719*5de52c6dSVaclav Hapla   PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_FALSE, *dmParallel));
1720721cbd76SMatthew G. Knepley   /* Cleanup */
1721f8987ae8SMichael Lange   if (sf) {*sf = sfMigration;}
17229566063dSJacob Faibussowitsch   else    PetscCall(PetscSFDestroy(&sfMigration));
17239566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfPoint));
17249566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0));
172570034214SMatthew G. Knepley   PetscFunctionReturn(0);
172670034214SMatthew G. Knepley }
1727a157612eSMichael Lange 
1728a157612eSMichael Lange /*@C
172924cc2ca5SMatthew G. Knepley   DMPlexDistributeOverlap - Add partition overlap to a distributed non-overlapping DM.
1730a157612eSMichael Lange 
1731d083f849SBarry Smith   Collective on dm
1732a157612eSMichael Lange 
1733064ec1f3Sprj-   Input Parameters:
1734064ec1f3Sprj- + dm  - The non-overlapping distributed DMPlex object
173557fe9a49SVaclav Hapla - overlap - The overlap of partitions (the same on all ranks)
1736a157612eSMichael Lange 
1737064ec1f3Sprj-   Output Parameters:
1738a157612eSMichael Lange + sf - The PetscSF used for point distribution
1739a157612eSMichael Lange - dmOverlap - The overlapping distributed DMPlex object, or NULL
1740a157612eSMichael Lange 
1741dccdeccaSVaclav Hapla   Notes:
1742dccdeccaSVaclav Hapla   If the mesh was not distributed, the return value is NULL.
1743a157612eSMichael Lange 
1744b0441da4SMatthew G. Knepley   The user can control the definition of adjacency for the mesh using DMSetAdjacency(). They should choose the combination appropriate for the function
1745a157612eSMichael Lange   representation on the mesh.
1746a157612eSMichael Lange 
1747dccdeccaSVaclav Hapla   Level: advanced
1748a157612eSMichael Lange 
1749db781477SPatrick Sanan .seealso: `DMPlexCreate()`, `DMSetAdjacency()`, `DMPlexDistribute()`, `DMPlexCreateOverlapLabel()`, `DMPlexGetOverlap()`
1750a157612eSMichael Lange @*/
1751b9f40539SMichael Lange PetscErrorCode DMPlexDistributeOverlap(DM dm, PetscInt overlap, PetscSF *sf, DM *dmOverlap)
1752a157612eSMichael Lange {
1753a157612eSMichael Lange   MPI_Comm               comm;
17543567eaeeSMatthew G. Knepley   PetscMPIInt            size, rank;
1755a157612eSMichael Lange   PetscSection           rootSection, leafSection;
1756a157612eSMichael Lange   IS                     rootrank, leafrank;
1757cf86098cSMatthew G. Knepley   DM                     dmCoord;
1758a9f1d5b2SMichael Lange   DMLabel                lblOverlap;
1759f5bf2dbfSMichael Lange   PetscSF                sfOverlap, sfStratified, sfPoint;
1760a157612eSMichael Lange 
1761a157612eSMichael Lange   PetscFunctionBegin;
1762a157612eSMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
176357fe9a49SVaclav Hapla   PetscValidLogicalCollectiveInt(dm, overlap, 2);
1764a157612eSMichael Lange   if (sf) PetscValidPointer(sf, 3);
1765a157612eSMichael Lange   PetscValidPointer(dmOverlap, 4);
1766a157612eSMichael Lange 
17670c86c063SLisandro Dalcin   if (sf) *sf = NULL;
17680c86c063SLisandro Dalcin   *dmOverlap  = NULL;
17699566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
17709566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm, &size));
17719566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm, &rank));
17720c86c063SLisandro Dalcin   if (size == 1) PetscFunctionReturn(0);
1773a157612eSMichael Lange 
17749566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_DistributeOverlap, dm, 0, 0, 0));
1775a157612eSMichael Lange   /* Compute point overlap with neighbouring processes on the distributed DM */
17769566063dSJacob Faibussowitsch   PetscCall(PetscLogEventBegin(DMPLEX_Partition,dm,0,0,0));
17779566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &rootSection));
17789566063dSJacob Faibussowitsch   PetscCall(PetscSectionCreate(comm, &leafSection));
17799566063dSJacob Faibussowitsch   PetscCall(DMPlexDistributeOwnership(dm, rootSection, &rootrank, leafSection, &leafrank));
17809566063dSJacob Faibussowitsch   PetscCall(DMPlexCreateOverlapLabel(dm, overlap, rootSection, rootrank, leafSection, leafrank, &lblOverlap));
1781a9f1d5b2SMichael Lange   /* Convert overlap label to stratified migration SF */
17829566063dSJacob Faibussowitsch   PetscCall(DMPlexPartitionLabelCreateSF(dm, lblOverlap, &sfOverlap));
17839566063dSJacob Faibussowitsch   PetscCall(DMPlexStratifyMigrationSF(dm, sfOverlap, &sfStratified));
17849566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfOverlap));
1785a9f1d5b2SMichael Lange   sfOverlap = sfStratified;
17869566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) sfOverlap, "Overlap SF"));
17879566063dSJacob Faibussowitsch   PetscCall(PetscSFSetFromOptions(sfOverlap));
1788a9f1d5b2SMichael Lange 
17899566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&rootSection));
17909566063dSJacob Faibussowitsch   PetscCall(PetscSectionDestroy(&leafSection));
17919566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&rootrank));
17929566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&leafrank));
17939566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_Partition,dm,0,0,0));
1794a157612eSMichael Lange 
1795a157612eSMichael Lange   /* Build the overlapping DM */
17969566063dSJacob Faibussowitsch   PetscCall(DMPlexCreate(comm, dmOverlap));
17979566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *dmOverlap, "Parallel Mesh"));
17989566063dSJacob Faibussowitsch   PetscCall(DMPlexMigrate(dm, sfOverlap, *dmOverlap));
1799cb54e036SVaclav Hapla   /* Store the overlap in the new DM */
1800cb54e036SVaclav Hapla   ((DM_Plex*)(*dmOverlap)->data)->overlap = overlap + ((DM_Plex*)dm->data)->overlap;
1801f5bf2dbfSMichael Lange   /* Build the new point SF */
18029566063dSJacob Faibussowitsch   PetscCall(DMPlexCreatePointSF(*dmOverlap, sfOverlap, PETSC_FALSE, &sfPoint));
18039566063dSJacob Faibussowitsch   PetscCall(DMSetPointSF(*dmOverlap, sfPoint));
18049566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(*dmOverlap, &dmCoord));
18059566063dSJacob Faibussowitsch   if (dmCoord) PetscCall(DMSetPointSF(dmCoord, sfPoint));
18069566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfPoint));
1807a157612eSMichael Lange   /* Cleanup overlap partition */
18089566063dSJacob Faibussowitsch   PetscCall(DMLabelDestroy(&lblOverlap));
1809a9f1d5b2SMichael Lange   if (sf) *sf = sfOverlap;
18109566063dSJacob Faibussowitsch   else    PetscCall(PetscSFDestroy(&sfOverlap));
18119566063dSJacob Faibussowitsch   PetscCall(PetscLogEventEnd(DMPLEX_DistributeOverlap, dm, 0, 0, 0));
1812a157612eSMichael Lange   PetscFunctionReturn(0);
1813a157612eSMichael Lange }
18146462276dSToby Isaac 
1815cb54e036SVaclav Hapla PetscErrorCode DMPlexGetOverlap_Plex(DM dm, PetscInt *overlap)
1816cb54e036SVaclav Hapla {
1817cb54e036SVaclav Hapla   DM_Plex        *mesh  = (DM_Plex*) dm->data;
1818cb54e036SVaclav Hapla 
1819cb54e036SVaclav Hapla   PetscFunctionBegin;
1820cb54e036SVaclav Hapla   *overlap = mesh->overlap;
1821cb54e036SVaclav Hapla   PetscFunctionReturn(0);
1822cb54e036SVaclav Hapla }
1823cb54e036SVaclav Hapla 
1824cb54e036SVaclav Hapla /*@
1825cb54e036SVaclav Hapla   DMPlexGetOverlap - Get the DMPlex partition overlap.
1826cb54e036SVaclav Hapla 
1827cb54e036SVaclav Hapla   Not collective
1828cb54e036SVaclav Hapla 
1829cb54e036SVaclav Hapla   Input Parameter:
1830cb54e036SVaclav Hapla . dm - The DM
1831cb54e036SVaclav Hapla 
1832064ec1f3Sprj-   Output Parameter:
1833cb54e036SVaclav Hapla . overlap - The overlap of this DM
1834cb54e036SVaclav Hapla 
1835cb54e036SVaclav Hapla   Level: intermediate
1836cb54e036SVaclav Hapla 
1837db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexDistributeOverlap()`, `DMPlexCreateOverlapLabel()`
1838cb54e036SVaclav Hapla @*/
1839cb54e036SVaclav Hapla PetscErrorCode DMPlexGetOverlap(DM dm, PetscInt *overlap)
1840cb54e036SVaclav Hapla {
1841cb54e036SVaclav Hapla   PetscFunctionBegin;
1842cb54e036SVaclav Hapla   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1843cac4c232SBarry Smith   PetscUseMethod(dm,"DMPlexGetOverlap_C",(DM,PetscInt*),(dm,overlap));
1844cb54e036SVaclav Hapla   PetscFunctionReturn(0);
1845cb54e036SVaclav Hapla }
1846cb54e036SVaclav Hapla 
1847e600fa54SMatthew G. Knepley PetscErrorCode DMPlexDistributeSetDefault_Plex(DM dm, PetscBool dist)
1848e600fa54SMatthew G. Knepley {
1849e600fa54SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
1850e600fa54SMatthew G. Knepley 
1851e600fa54SMatthew G. Knepley   PetscFunctionBegin;
1852e600fa54SMatthew G. Knepley   mesh->distDefault = dist;
1853e600fa54SMatthew G. Knepley   PetscFunctionReturn(0);
1854e600fa54SMatthew G. Knepley }
1855e600fa54SMatthew G. Knepley 
1856e600fa54SMatthew G. Knepley /*@
1857e600fa54SMatthew G. Knepley   DMPlexDistributeSetDefault - Set flag indicating whether the DM should be distributed by default
1858e600fa54SMatthew G. Knepley 
1859e600fa54SMatthew G. Knepley   Logically collective
1860e600fa54SMatthew G. Knepley 
1861e600fa54SMatthew G. Knepley   Input Parameters:
1862e600fa54SMatthew G. Knepley + dm   - The DM
1863e600fa54SMatthew G. Knepley - dist - Flag for distribution
1864e600fa54SMatthew G. Knepley 
1865e600fa54SMatthew G. Knepley   Level: intermediate
1866e600fa54SMatthew G. Knepley 
1867db781477SPatrick Sanan .seealso: `DMDistributeGetDefault()`, `DMPlexDistribute()`
1868e600fa54SMatthew G. Knepley @*/
1869e600fa54SMatthew G. Knepley PetscErrorCode DMPlexDistributeSetDefault(DM dm, PetscBool dist)
1870e600fa54SMatthew G. Knepley {
1871e600fa54SMatthew G. Knepley   PetscFunctionBegin;
1872e600fa54SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1873e600fa54SMatthew G. Knepley   PetscValidLogicalCollectiveBool(dm, dist, 2);
1874cac4c232SBarry Smith   PetscTryMethod(dm,"DMPlexDistributeSetDefault_C",(DM,PetscBool),(dm,dist));
1875e600fa54SMatthew G. Knepley   PetscFunctionReturn(0);
1876e600fa54SMatthew G. Knepley }
1877e600fa54SMatthew G. Knepley 
1878e600fa54SMatthew G. Knepley PetscErrorCode DMPlexDistributeGetDefault_Plex(DM dm, PetscBool *dist)
1879e600fa54SMatthew G. Knepley {
1880e600fa54SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
1881e600fa54SMatthew G. Knepley 
1882e600fa54SMatthew G. Knepley   PetscFunctionBegin;
1883e600fa54SMatthew G. Knepley   *dist = mesh->distDefault;
1884e600fa54SMatthew G. Knepley   PetscFunctionReturn(0);
1885e600fa54SMatthew G. Knepley }
1886e600fa54SMatthew G. Knepley 
1887e600fa54SMatthew G. Knepley /*@
1888e600fa54SMatthew G. Knepley   DMPlexDistributeGetDefault - Get flag indicating whether the DM should be distributed by default
1889e600fa54SMatthew G. Knepley 
1890e600fa54SMatthew G. Knepley   Not collective
1891e600fa54SMatthew G. Knepley 
1892e600fa54SMatthew G. Knepley   Input Parameter:
1893e600fa54SMatthew G. Knepley . dm   - The DM
1894e600fa54SMatthew G. Knepley 
1895e600fa54SMatthew G. Knepley   Output Parameter:
1896e600fa54SMatthew G. Knepley . dist - Flag for distribution
1897e600fa54SMatthew G. Knepley 
1898e600fa54SMatthew G. Knepley   Level: intermediate
1899e600fa54SMatthew G. Knepley 
1900db781477SPatrick Sanan .seealso: `DMDistributeSetDefault()`, `DMPlexDistribute()`
1901e600fa54SMatthew G. Knepley @*/
1902e600fa54SMatthew G. Knepley PetscErrorCode DMPlexDistributeGetDefault(DM dm, PetscBool *dist)
1903e600fa54SMatthew G. Knepley {
1904e600fa54SMatthew G. Knepley   PetscFunctionBegin;
1905e600fa54SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1906e600fa54SMatthew G. Knepley   PetscValidBoolPointer(dist, 2);
1907cac4c232SBarry Smith   PetscUseMethod(dm,"DMPlexDistributeGetDefault_C",(DM,PetscBool*),(dm,dist));
1908e600fa54SMatthew G. Knepley   PetscFunctionReturn(0);
1909e600fa54SMatthew G. Knepley }
1910e600fa54SMatthew G. Knepley 
19116462276dSToby Isaac /*@C
19126462276dSToby Isaac   DMPlexGetGatherDM - Get a copy of the DMPlex that gathers all points on the
19136462276dSToby Isaac   root process of the original's communicator.
19146462276dSToby Isaac 
1915d083f849SBarry Smith   Collective on dm
191683655b49SVáclav Hapla 
1917064ec1f3Sprj-   Input Parameter:
19186462276dSToby Isaac . dm - the original DMPlex object
19196462276dSToby Isaac 
19206462276dSToby Isaac   Output Parameters:
1921a13df41bSStefano Zampini + sf - the PetscSF used for point distribution (optional)
1922a13df41bSStefano Zampini - gatherMesh - the gathered DM object, or NULL
19236462276dSToby Isaac 
19246462276dSToby Isaac   Level: intermediate
19256462276dSToby Isaac 
1926db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexGetRedundantDM()`
19276462276dSToby Isaac @*/
1928a13df41bSStefano Zampini PetscErrorCode DMPlexGetGatherDM(DM dm, PetscSF *sf, DM *gatherMesh)
19296462276dSToby Isaac {
19306462276dSToby Isaac   MPI_Comm       comm;
19316462276dSToby Isaac   PetscMPIInt    size;
19326462276dSToby Isaac   PetscPartitioner oldPart, gatherPart;
19336462276dSToby Isaac 
19346462276dSToby Isaac   PetscFunctionBegin;
19356462276dSToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1936064a246eSJacob Faibussowitsch   PetscValidPointer(gatherMesh,3);
19370c86c063SLisandro Dalcin   *gatherMesh = NULL;
1938a13df41bSStefano Zampini   if (sf) *sf = NULL;
19396462276dSToby Isaac   comm = PetscObjectComm((PetscObject)dm);
19409566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
19416462276dSToby Isaac   if (size == 1) PetscFunctionReturn(0);
19429566063dSJacob Faibussowitsch   PetscCall(DMPlexGetPartitioner(dm,&oldPart));
19439566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)oldPart));
19449566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerCreate(comm,&gatherPart));
19459566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerSetType(gatherPart,PETSCPARTITIONERGATHER));
19469566063dSJacob Faibussowitsch   PetscCall(DMPlexSetPartitioner(dm,gatherPart));
19479566063dSJacob Faibussowitsch   PetscCall(DMPlexDistribute(dm,0,sf,gatherMesh));
1948a13df41bSStefano Zampini 
19499566063dSJacob Faibussowitsch   PetscCall(DMPlexSetPartitioner(dm,oldPart));
19509566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerDestroy(&gatherPart));
19519566063dSJacob Faibussowitsch   PetscCall(PetscPartitionerDestroy(&oldPart));
19526462276dSToby Isaac   PetscFunctionReturn(0);
19536462276dSToby Isaac }
19546462276dSToby Isaac 
19556462276dSToby Isaac /*@C
19566462276dSToby Isaac   DMPlexGetRedundantDM - Get a copy of the DMPlex that is completely copied on each process.
19576462276dSToby Isaac 
1958d083f849SBarry Smith   Collective on dm
195983655b49SVáclav Hapla 
1960064ec1f3Sprj-   Input Parameter:
19616462276dSToby Isaac . dm - the original DMPlex object
19626462276dSToby Isaac 
19636462276dSToby Isaac   Output Parameters:
1964a13df41bSStefano Zampini + sf - the PetscSF used for point distribution (optional)
1965a13df41bSStefano Zampini - redundantMesh - the redundant DM object, or NULL
19666462276dSToby Isaac 
19676462276dSToby Isaac   Level: intermediate
19686462276dSToby Isaac 
1969db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexGetGatherDM()`
19706462276dSToby Isaac @*/
1971a13df41bSStefano Zampini PetscErrorCode DMPlexGetRedundantDM(DM dm, PetscSF *sf, DM *redundantMesh)
19726462276dSToby Isaac {
19736462276dSToby Isaac   MPI_Comm       comm;
19746462276dSToby Isaac   PetscMPIInt    size, rank;
19756462276dSToby Isaac   PetscInt       pStart, pEnd, p;
19766462276dSToby Isaac   PetscInt       numPoints = -1;
1977a13df41bSStefano Zampini   PetscSF        migrationSF, sfPoint, gatherSF;
19786462276dSToby Isaac   DM             gatherDM, dmCoord;
19796462276dSToby Isaac   PetscSFNode    *points;
19806462276dSToby Isaac 
19816462276dSToby Isaac   PetscFunctionBegin;
19826462276dSToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1983064a246eSJacob Faibussowitsch   PetscValidPointer(redundantMesh,3);
19840c86c063SLisandro Dalcin   *redundantMesh = NULL;
19856462276dSToby Isaac   comm = PetscObjectComm((PetscObject)dm);
19869566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
198768dbc166SMatthew G. Knepley   if (size == 1) {
19889566063dSJacob Faibussowitsch     PetscCall(PetscObjectReference((PetscObject) dm));
198968dbc166SMatthew G. Knepley     *redundantMesh = dm;
1990a13df41bSStefano Zampini     if (sf) *sf = NULL;
199168dbc166SMatthew G. Knepley     PetscFunctionReturn(0);
199268dbc166SMatthew G. Knepley   }
19939566063dSJacob Faibussowitsch   PetscCall(DMPlexGetGatherDM(dm,&gatherSF,&gatherDM));
19946462276dSToby Isaac   if (!gatherDM) PetscFunctionReturn(0);
19959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(comm,&rank));
19969566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(gatherDM,&pStart,&pEnd));
19976462276dSToby Isaac   numPoints = pEnd - pStart;
19989566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Bcast(&numPoints,1,MPIU_INT,0,comm));
19999566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(numPoints,&points));
20009566063dSJacob Faibussowitsch   PetscCall(PetscSFCreate(comm,&migrationSF));
20016462276dSToby Isaac   for (p = 0; p < numPoints; p++) {
20026462276dSToby Isaac     points[p].index = p;
20036462276dSToby Isaac     points[p].rank  = 0;
20046462276dSToby Isaac   }
20059566063dSJacob Faibussowitsch   PetscCall(PetscSFSetGraph(migrationSF,pEnd-pStart,numPoints,NULL,PETSC_OWN_POINTER,points,PETSC_OWN_POINTER));
20069566063dSJacob Faibussowitsch   PetscCall(DMPlexCreate(comm, redundantMesh));
20079566063dSJacob Faibussowitsch   PetscCall(PetscObjectSetName((PetscObject) *redundantMesh, "Redundant Mesh"));
20089566063dSJacob Faibussowitsch   PetscCall(DMPlexMigrate(gatherDM, migrationSF, *redundantMesh));
20099566063dSJacob Faibussowitsch   PetscCall(DMPlexCreatePointSF(*redundantMesh, migrationSF, PETSC_FALSE, &sfPoint));
20109566063dSJacob Faibussowitsch   PetscCall(DMSetPointSF(*redundantMesh, sfPoint));
20119566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(*redundantMesh, &dmCoord));
20129566063dSJacob Faibussowitsch   if (dmCoord) PetscCall(DMSetPointSF(dmCoord, sfPoint));
20139566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&sfPoint));
2014a13df41bSStefano Zampini   if (sf) {
2015a13df41bSStefano Zampini     PetscSF tsf;
2016a13df41bSStefano Zampini 
20179566063dSJacob Faibussowitsch     PetscCall(PetscSFCompose(gatherSF,migrationSF,&tsf));
20189566063dSJacob Faibussowitsch     PetscCall(DMPlexStratifyMigrationSF(dm, tsf, sf));
20199566063dSJacob Faibussowitsch     PetscCall(PetscSFDestroy(&tsf));
2020a13df41bSStefano Zampini   }
20219566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&migrationSF));
20229566063dSJacob Faibussowitsch   PetscCall(PetscSFDestroy(&gatherSF));
20239566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&gatherDM));
20249566063dSJacob Faibussowitsch   PetscCall(DMCopyDisc(dm, *redundantMesh));
2025*5de52c6dSVaclav Hapla   PetscCall(DMPlexCopy_Internal(dm, PETSC_TRUE, PETSC_FALSE, *redundantMesh));
20266462276dSToby Isaac   PetscFunctionReturn(0);
20276462276dSToby Isaac }
20285fa78c88SVaclav Hapla 
20295fa78c88SVaclav Hapla /*@
20305fa78c88SVaclav Hapla   DMPlexIsDistributed - Find out whether this DM is distributed, i.e. more than one rank owns some points.
20315fa78c88SVaclav Hapla 
20325fa78c88SVaclav Hapla   Collective
20335fa78c88SVaclav Hapla 
20345fa78c88SVaclav Hapla   Input Parameter:
20355fa78c88SVaclav Hapla . dm      - The DM object
20365fa78c88SVaclav Hapla 
20375fa78c88SVaclav Hapla   Output Parameter:
20385fa78c88SVaclav Hapla . distributed - Flag whether the DM is distributed
20395fa78c88SVaclav Hapla 
20405fa78c88SVaclav Hapla   Level: intermediate
20415fa78c88SVaclav Hapla 
20425fa78c88SVaclav Hapla   Notes:
20435fa78c88SVaclav Hapla   This currently finds out whether at least two ranks have any DAG points.
20445fa78c88SVaclav Hapla   This involves MPI_Allreduce() with one integer.
20455fa78c88SVaclav Hapla   The result is currently not stashed so every call to this routine involves this global communication.
20465fa78c88SVaclav Hapla 
2047db781477SPatrick Sanan .seealso: `DMPlexDistribute()`, `DMPlexGetOverlap()`, `DMPlexIsInterpolated()`
20485fa78c88SVaclav Hapla @*/
20495fa78c88SVaclav Hapla PetscErrorCode DMPlexIsDistributed(DM dm, PetscBool *distributed)
20505fa78c88SVaclav Hapla {
20515fa78c88SVaclav Hapla   PetscInt          pStart, pEnd, count;
20525fa78c88SVaclav Hapla   MPI_Comm          comm;
205335d70e31SStefano Zampini   PetscMPIInt       size;
20545fa78c88SVaclav Hapla 
20555fa78c88SVaclav Hapla   PetscFunctionBegin;
20565fa78c88SVaclav Hapla   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2057dadcf809SJacob Faibussowitsch   PetscValidBoolPointer(distributed,2);
20589566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
20599566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(comm,&size));
206035d70e31SStefano Zampini   if (size == 1) { *distributed = PETSC_FALSE; PetscFunctionReturn(0); }
20619566063dSJacob Faibussowitsch   PetscCall(DMPlexGetChart(dm, &pStart, &pEnd));
206235d70e31SStefano Zampini   count = (pEnd - pStart) > 0 ? 1 : 0;
20639566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(MPI_IN_PLACE, &count, 1, MPIU_INT, MPI_SUM, comm));
20645fa78c88SVaclav Hapla   *distributed = count > 1 ? PETSC_TRUE : PETSC_FALSE;
20655fa78c88SVaclav Hapla   PetscFunctionReturn(0);
20665fa78c88SVaclav Hapla }
2067