xref: /petsc/src/dm/impls/plex/plexdistribute.c (revision 854ce69b31eecaec02632e265405c6229e3f3488)
170034214SMatthew G. Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
270034214SMatthew G. Knepley #include <petscsf.h>
370034214SMatthew G. Knepley 
470034214SMatthew G. Knepley #undef __FUNCT__
570034214SMatthew G. Knepley #define __FUNCT__ "DMPlexSetAdjacencyUseCone"
670034214SMatthew G. Knepley /*@
770034214SMatthew G. Knepley   DMPlexSetAdjacencyUseCone - Define adjacency in the mesh using either the cone or the support first
870034214SMatthew G. Knepley 
970034214SMatthew G. Knepley   Input Parameters:
1070034214SMatthew G. Knepley + dm      - The DM object
1170034214SMatthew G. Knepley - useCone - Flag to use the cone first
1270034214SMatthew G. Knepley 
1370034214SMatthew G. Knepley   Level: intermediate
1470034214SMatthew G. Knepley 
1570034214SMatthew G. Knepley   Notes:
1670034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
174b6b44bdSMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
1870034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
1970034214SMatthew G. Knepley 
2070034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
2170034214SMatthew G. Knepley @*/
2270034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseCone(DM dm, PetscBool useCone)
2370034214SMatthew G. Knepley {
2470034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
2570034214SMatthew G. Knepley 
2670034214SMatthew G. Knepley   PetscFunctionBegin;
2770034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2870034214SMatthew G. Knepley   mesh->useCone = useCone;
2970034214SMatthew G. Knepley   PetscFunctionReturn(0);
3070034214SMatthew G. Knepley }
3170034214SMatthew G. Knepley 
3270034214SMatthew G. Knepley #undef __FUNCT__
3370034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacencyUseCone"
3470034214SMatthew G. Knepley /*@
3570034214SMatthew G. Knepley   DMPlexGetAdjacencyUseCone - Query whether adjacency in the mesh uses the cone or the support first
3670034214SMatthew G. Knepley 
3770034214SMatthew G. Knepley   Input Parameter:
3870034214SMatthew G. Knepley . dm      - The DM object
3970034214SMatthew G. Knepley 
4070034214SMatthew G. Knepley   Output Parameter:
4170034214SMatthew G. Knepley . useCone - Flag to use the cone first
4270034214SMatthew G. Knepley 
4370034214SMatthew G. Knepley   Level: intermediate
4470034214SMatthew G. Knepley 
4570034214SMatthew G. Knepley   Notes:
4670034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
474b6b44bdSMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
4870034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
4970034214SMatthew G. Knepley 
5070034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
5170034214SMatthew G. Knepley @*/
5270034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseCone(DM dm, PetscBool *useCone)
5370034214SMatthew G. Knepley {
5470034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
5570034214SMatthew G. Knepley 
5670034214SMatthew G. Knepley   PetscFunctionBegin;
5770034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5870034214SMatthew G. Knepley   PetscValidIntPointer(useCone, 2);
5970034214SMatthew G. Knepley   *useCone = mesh->useCone;
6070034214SMatthew G. Knepley   PetscFunctionReturn(0);
6170034214SMatthew G. Knepley }
6270034214SMatthew G. Knepley 
6370034214SMatthew G. Knepley #undef __FUNCT__
6470034214SMatthew G. Knepley #define __FUNCT__ "DMPlexSetAdjacencyUseClosure"
6570034214SMatthew G. Knepley /*@
6670034214SMatthew G. Knepley   DMPlexSetAdjacencyUseClosure - Define adjacency in the mesh using the transitive closure
6770034214SMatthew G. Knepley 
6870034214SMatthew G. Knepley   Input Parameters:
6970034214SMatthew G. Knepley + dm      - The DM object
7070034214SMatthew G. Knepley - useClosure - Flag to use the closure
7170034214SMatthew G. Knepley 
7270034214SMatthew G. Knepley   Level: intermediate
7370034214SMatthew G. Knepley 
7470034214SMatthew G. Knepley   Notes:
7570034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
764b6b44bdSMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
7770034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
7870034214SMatthew G. Knepley 
7970034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
8070034214SMatthew G. Knepley @*/
8170034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseClosure(DM dm, PetscBool useClosure)
8270034214SMatthew G. Knepley {
8370034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
8470034214SMatthew G. Knepley 
8570034214SMatthew G. Knepley   PetscFunctionBegin;
8670034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8770034214SMatthew G. Knepley   mesh->useClosure = useClosure;
8870034214SMatthew G. Knepley   PetscFunctionReturn(0);
8970034214SMatthew G. Knepley }
9070034214SMatthew G. Knepley 
9170034214SMatthew G. Knepley #undef __FUNCT__
9270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacencyUseClosure"
9370034214SMatthew G. Knepley /*@
9470034214SMatthew G. Knepley   DMPlexGetAdjacencyUseClosure - Query whether adjacency in the mesh uses the transitive closure
9570034214SMatthew G. Knepley 
9670034214SMatthew G. Knepley   Input Parameter:
9770034214SMatthew G. Knepley . dm      - The DM object
9870034214SMatthew G. Knepley 
9970034214SMatthew G. Knepley   Output Parameter:
10070034214SMatthew G. Knepley . useClosure - Flag to use the closure
10170034214SMatthew G. Knepley 
10270034214SMatthew G. Knepley   Level: intermediate
10370034214SMatthew G. Knepley 
10470034214SMatthew G. Knepley   Notes:
10570034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
1064b6b44bdSMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
10770034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
10870034214SMatthew G. Knepley 
10970034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
11070034214SMatthew G. Knepley @*/
11170034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseClosure(DM dm, PetscBool *useClosure)
11270034214SMatthew G. Knepley {
11370034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
11470034214SMatthew G. Knepley 
11570034214SMatthew G. Knepley   PetscFunctionBegin;
11670034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11770034214SMatthew G. Knepley   PetscValidIntPointer(useClosure, 2);
11870034214SMatthew G. Knepley   *useClosure = mesh->useClosure;
11970034214SMatthew G. Knepley   PetscFunctionReturn(0);
12070034214SMatthew G. Knepley }
12170034214SMatthew G. Knepley 
12270034214SMatthew G. Knepley #undef __FUNCT__
123a17985deSToby Isaac #define __FUNCT__ "DMPlexSetAdjacencyUseAnchors"
1248b0b4c70SToby Isaac /*@
125a17985deSToby Isaac   DMPlexSetAdjacencyUseAnchors - Define adjacency in the mesh using the point-to-point constraints.
1268b0b4c70SToby Isaac 
1278b0b4c70SToby Isaac   Input Parameters:
1288b0b4c70SToby Isaac + dm      - The DM object
1295b317d89SToby 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.
1308b0b4c70SToby Isaac 
1318b0b4c70SToby Isaac   Level: intermediate
1328b0b4c70SToby Isaac 
133a17985deSToby Isaac .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
1348b0b4c70SToby Isaac @*/
1355b317d89SToby Isaac PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM dm, PetscBool useAnchors)
1368b0b4c70SToby Isaac {
1378b0b4c70SToby Isaac   DM_Plex *mesh = (DM_Plex *) dm->data;
1388b0b4c70SToby Isaac 
1398b0b4c70SToby Isaac   PetscFunctionBegin;
1408b0b4c70SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1415b317d89SToby Isaac   mesh->useAnchors = useAnchors;
1428b0b4c70SToby Isaac   PetscFunctionReturn(0);
1438b0b4c70SToby Isaac }
1448b0b4c70SToby Isaac 
1458b0b4c70SToby Isaac #undef __FUNCT__
146a17985deSToby Isaac #define __FUNCT__ "DMPlexGetAdjacencyUseAnchors"
1478b0b4c70SToby Isaac /*@
148a17985deSToby Isaac   DMPlexGetAdjacencyUseAnchors - Query whether adjacency in the mesh uses the point-to-point constraints.
1498b0b4c70SToby Isaac 
1508b0b4c70SToby Isaac   Input Parameter:
1518b0b4c70SToby Isaac . dm      - The DM object
1528b0b4c70SToby Isaac 
1538b0b4c70SToby Isaac   Output Parameter:
1545b317d89SToby 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.
1558b0b4c70SToby Isaac 
1568b0b4c70SToby Isaac   Level: intermediate
1578b0b4c70SToby Isaac 
158a17985deSToby Isaac .seealso: DMPlexSetAdjacencyUseAnchors(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
1598b0b4c70SToby Isaac @*/
1605b317d89SToby Isaac PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM dm, PetscBool *useAnchors)
1618b0b4c70SToby Isaac {
1628b0b4c70SToby Isaac   DM_Plex *mesh = (DM_Plex *) dm->data;
1638b0b4c70SToby Isaac 
1648b0b4c70SToby Isaac   PetscFunctionBegin;
1658b0b4c70SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1665b317d89SToby Isaac   PetscValidIntPointer(useAnchors, 2);
1675b317d89SToby Isaac   *useAnchors = mesh->useAnchors;
1688b0b4c70SToby Isaac   PetscFunctionReturn(0);
1698b0b4c70SToby Isaac }
1708b0b4c70SToby Isaac 
1718b0b4c70SToby Isaac #undef __FUNCT__
17270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Cone_Internal"
17370034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Cone_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
17470034214SMatthew G. Knepley {
17570034214SMatthew G. Knepley   const PetscInt *cone = NULL;
17670034214SMatthew G. Knepley   PetscInt        numAdj = 0, maxAdjSize = *adjSize, coneSize, c;
17770034214SMatthew G. Knepley   PetscErrorCode  ierr;
17870034214SMatthew G. Knepley 
17970034214SMatthew G. Knepley   PetscFunctionBeginHot;
18070034214SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
18170034214SMatthew G. Knepley   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
1824b6b44bdSMatthew G. Knepley   for (c = 0; c <= coneSize; ++c) {
1834b6b44bdSMatthew G. Knepley     const PetscInt  point   = !c ? p : cone[c-1];
18470034214SMatthew G. Knepley     const PetscInt *support = NULL;
18570034214SMatthew G. Knepley     PetscInt        supportSize, s, q;
18670034214SMatthew G. Knepley 
1874b6b44bdSMatthew G. Knepley     ierr = DMPlexGetSupportSize(dm, point, &supportSize);CHKERRQ(ierr);
1884b6b44bdSMatthew G. Knepley     ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr);
18970034214SMatthew G. Knepley     for (s = 0; s < supportSize; ++s) {
19070034214SMatthew G. Knepley       for (q = 0; q < numAdj || (adj[numAdj++] = support[s],0); ++q) {
19170034214SMatthew G. Knepley         if (support[s] == adj[q]) break;
19270034214SMatthew G. Knepley       }
19370034214SMatthew G. Knepley       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
19470034214SMatthew G. Knepley     }
19570034214SMatthew G. Knepley   }
19670034214SMatthew G. Knepley   *adjSize = numAdj;
19770034214SMatthew G. Knepley   PetscFunctionReturn(0);
19870034214SMatthew G. Knepley }
19970034214SMatthew G. Knepley 
20070034214SMatthew G. Knepley #undef __FUNCT__
20170034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Support_Internal"
20270034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Support_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
20370034214SMatthew G. Knepley {
20470034214SMatthew G. Knepley   const PetscInt *support = NULL;
20570034214SMatthew G. Knepley   PetscInt        numAdj   = 0, maxAdjSize = *adjSize, supportSize, s;
20670034214SMatthew G. Knepley   PetscErrorCode  ierr;
20770034214SMatthew G. Knepley 
20870034214SMatthew G. Knepley   PetscFunctionBeginHot;
20970034214SMatthew G. Knepley   ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr);
21070034214SMatthew G. Knepley   ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr);
2114b6b44bdSMatthew G. Knepley   for (s = 0; s <= supportSize; ++s) {
2124b6b44bdSMatthew G. Knepley     const PetscInt  point = !s ? p : support[s-1];
21370034214SMatthew G. Knepley     const PetscInt *cone  = NULL;
21470034214SMatthew G. Knepley     PetscInt        coneSize, c, q;
21570034214SMatthew G. Knepley 
2164b6b44bdSMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
2174b6b44bdSMatthew G. Knepley     ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
21870034214SMatthew G. Knepley     for (c = 0; c < coneSize; ++c) {
21970034214SMatthew G. Knepley       for (q = 0; q < numAdj || (adj[numAdj++] = cone[c],0); ++q) {
22070034214SMatthew G. Knepley         if (cone[c] == adj[q]) break;
22170034214SMatthew G. Knepley       }
22270034214SMatthew G. Knepley       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
22370034214SMatthew G. Knepley     }
22470034214SMatthew G. Knepley   }
22570034214SMatthew G. Knepley   *adjSize = numAdj;
22670034214SMatthew G. Knepley   PetscFunctionReturn(0);
22770034214SMatthew G. Knepley }
22870034214SMatthew G. Knepley 
22970034214SMatthew G. Knepley #undef __FUNCT__
23070034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Transitive_Internal"
23170034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Transitive_Internal(DM dm, PetscInt p, PetscBool useClosure, PetscInt *adjSize, PetscInt adj[])
23270034214SMatthew G. Knepley {
23370034214SMatthew G. Knepley   PetscInt      *star = NULL;
23470034214SMatthew G. Knepley   PetscInt       numAdj = 0, maxAdjSize = *adjSize, starSize, s;
23570034214SMatthew G. Knepley   PetscErrorCode ierr;
23670034214SMatthew G. Knepley 
23770034214SMatthew G. Knepley   PetscFunctionBeginHot;
23870034214SMatthew G. Knepley   ierr = DMPlexGetTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr);
23970034214SMatthew G. Knepley   for (s = 0; s < starSize*2; s += 2) {
24070034214SMatthew G. Knepley     const PetscInt *closure = NULL;
24170034214SMatthew G. Knepley     PetscInt        closureSize, c, q;
24270034214SMatthew G. Knepley 
24370034214SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
24470034214SMatthew G. Knepley     for (c = 0; c < closureSize*2; c += 2) {
24570034214SMatthew G. Knepley       for (q = 0; q < numAdj || (adj[numAdj++] = closure[c],0); ++q) {
24670034214SMatthew G. Knepley         if (closure[c] == adj[q]) break;
24770034214SMatthew G. Knepley       }
24870034214SMatthew G. Knepley       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
24970034214SMatthew G. Knepley     }
25070034214SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
25170034214SMatthew G. Knepley   }
25270034214SMatthew G. Knepley   ierr = DMPlexRestoreTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr);
25370034214SMatthew G. Knepley   *adjSize = numAdj;
25470034214SMatthew G. Knepley   PetscFunctionReturn(0);
25570034214SMatthew G. Knepley }
25670034214SMatthew G. Knepley 
25770034214SMatthew G. Knepley #undef __FUNCT__
25870034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Internal"
2595b317d89SToby Isaac PetscErrorCode DMPlexGetAdjacency_Internal(DM dm, PetscInt p, PetscBool useCone, PetscBool useTransitiveClosure, PetscBool useAnchors, PetscInt *adjSize, PetscInt *adj[])
26070034214SMatthew G. Knepley {
26179bad331SMatthew G. Knepley   static PetscInt asiz = 0;
2628b0b4c70SToby Isaac   PetscInt maxAnchors = 1;
2638b0b4c70SToby Isaac   PetscInt aStart = -1, aEnd = -1;
2648b0b4c70SToby Isaac   PetscInt maxAdjSize;
2658b0b4c70SToby Isaac   PetscSection aSec = NULL;
2668b0b4c70SToby Isaac   IS aIS = NULL;
2678b0b4c70SToby Isaac   const PetscInt *anchors;
26870034214SMatthew G. Knepley   PetscErrorCode  ierr;
26970034214SMatthew G. Knepley 
27070034214SMatthew G. Knepley   PetscFunctionBeginHot;
2715b317d89SToby Isaac   if (useAnchors) {
272a17985deSToby Isaac     ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr);
2738b0b4c70SToby Isaac     if (aSec) {
2748b0b4c70SToby Isaac       ierr = PetscSectionGetMaxDof(aSec,&maxAnchors);CHKERRQ(ierr);
27524c766afSToby Isaac       maxAnchors = PetscMax(1,maxAnchors);
2768b0b4c70SToby Isaac       ierr = PetscSectionGetChart(aSec,&aStart,&aEnd);CHKERRQ(ierr);
2778b0b4c70SToby Isaac       ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr);
2788b0b4c70SToby Isaac     }
2798b0b4c70SToby Isaac   }
28079bad331SMatthew G. Knepley   if (!*adj) {
28124c766afSToby Isaac     PetscInt depth, coneSeries, supportSeries, maxC, maxS, pStart, pEnd;
28279bad331SMatthew G. Knepley 
28324c766afSToby Isaac     ierr  = DMPlexGetChart(dm, &pStart,&pEnd);CHKERRQ(ierr);
28479bad331SMatthew G. Knepley     ierr  = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
28524c766afSToby Isaac     ierr  = DMPlexGetMaxSizes(dm, &maxC, &maxS);CHKERRQ(ierr);
28624c766afSToby Isaac     coneSeries    = (maxC > 1) ? ((PetscPowInt(maxC,depth+1)-1)/(maxC-1)) : depth+1;
28724c766afSToby Isaac     supportSeries = (maxS > 1) ? ((PetscPowInt(maxS,depth+1)-1)/(maxS-1)) : depth+1;
28824c766afSToby Isaac     asiz  = PetscMax(PetscPowInt(maxS,depth)*coneSeries,PetscPowInt(maxC,depth)*supportSeries);
2898b0b4c70SToby Isaac     asiz *= maxAnchors;
29024c766afSToby Isaac     asiz  = PetscMin(asiz,pEnd-pStart);
29179bad331SMatthew G. Knepley     ierr  = PetscMalloc1(asiz,adj);CHKERRQ(ierr);
29279bad331SMatthew G. Knepley   }
29379bad331SMatthew G. Knepley   if (*adjSize < 0) *adjSize = asiz;
2948b0b4c70SToby Isaac   maxAdjSize = *adjSize;
29570034214SMatthew G. Knepley   if (useTransitiveClosure) {
29679bad331SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Transitive_Internal(dm, p, useCone, adjSize, *adj);CHKERRQ(ierr);
29770034214SMatthew G. Knepley   } else if (useCone) {
29879bad331SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Cone_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
29970034214SMatthew G. Knepley   } else {
30079bad331SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Support_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
30170034214SMatthew G. Knepley   }
3025b317d89SToby Isaac   if (useAnchors && aSec) {
3038b0b4c70SToby Isaac     PetscInt origSize = *adjSize;
3048b0b4c70SToby Isaac     PetscInt numAdj = origSize;
3058b0b4c70SToby Isaac     PetscInt i = 0, j;
3068b0b4c70SToby Isaac     PetscInt *orig = *adj;
3078b0b4c70SToby Isaac 
3088b0b4c70SToby Isaac     while (i < origSize) {
3098b0b4c70SToby Isaac       PetscInt p = orig[i];
3108b0b4c70SToby Isaac       PetscInt aDof = 0;
3118b0b4c70SToby Isaac 
3128b0b4c70SToby Isaac       if (p >= aStart && p < aEnd) {
3138b0b4c70SToby Isaac         ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr);
3148b0b4c70SToby Isaac       }
3158b0b4c70SToby Isaac       if (aDof) {
3168b0b4c70SToby Isaac         PetscInt aOff;
3178b0b4c70SToby Isaac         PetscInt s, q;
3188b0b4c70SToby Isaac 
3198b0b4c70SToby Isaac         for (j = i + 1; j < numAdj; j++) {
3208b0b4c70SToby Isaac           orig[j - 1] = orig[j];
3218b0b4c70SToby Isaac         }
3228b0b4c70SToby Isaac         origSize--;
3238b0b4c70SToby Isaac         numAdj--;
3248b0b4c70SToby Isaac         ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr);
3258b0b4c70SToby Isaac         for (s = 0; s < aDof; ++s) {
3268b0b4c70SToby Isaac           for (q = 0; q < numAdj || (orig[numAdj++] = anchors[aOff+s],0); ++q) {
3278b0b4c70SToby Isaac             if (anchors[aOff+s] == orig[q]) break;
3288b0b4c70SToby Isaac           }
3298b0b4c70SToby Isaac           if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
3308b0b4c70SToby Isaac         }
3318b0b4c70SToby Isaac       }
3328b0b4c70SToby Isaac       else {
3338b0b4c70SToby Isaac         i++;
3348b0b4c70SToby Isaac       }
3358b0b4c70SToby Isaac     }
3368b0b4c70SToby Isaac     *adjSize = numAdj;
3378b0b4c70SToby Isaac     ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr);
3388b0b4c70SToby Isaac   }
33970034214SMatthew G. Knepley   PetscFunctionReturn(0);
34070034214SMatthew G. Knepley }
34170034214SMatthew G. Knepley 
34270034214SMatthew G. Knepley #undef __FUNCT__
34370034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency"
34470034214SMatthew G. Knepley /*@
34570034214SMatthew G. Knepley   DMPlexGetAdjacency - Return all points adjacent to the given point
34670034214SMatthew G. Knepley 
34770034214SMatthew G. Knepley   Input Parameters:
34870034214SMatthew G. Knepley + dm - The DM object
34970034214SMatthew G. Knepley . p  - The point
35070034214SMatthew G. Knepley . adjSize - The maximum size of adj if it is non-NULL, or PETSC_DETERMINE
35170034214SMatthew G. Knepley - adj - Either NULL so that the array is allocated, or an existing array with size adjSize
35270034214SMatthew G. Knepley 
35370034214SMatthew G. Knepley   Output Parameters:
35470034214SMatthew G. Knepley + adjSize - The number of adjacent points
35570034214SMatthew G. Knepley - adj - The adjacent points
35670034214SMatthew G. Knepley 
35770034214SMatthew G. Knepley   Level: advanced
35870034214SMatthew G. Knepley 
35970034214SMatthew G. Knepley   Notes: The user must PetscFree the adj array if it was not passed in.
36070034214SMatthew G. Knepley 
36170034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexDistribute(), DMCreateMatrix(), DMPlexPreallocateOperator()
36270034214SMatthew G. Knepley @*/
36370034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacency(DM dm, PetscInt p, PetscInt *adjSize, PetscInt *adj[])
36470034214SMatthew G. Knepley {
36570034214SMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
36670034214SMatthew G. Knepley   PetscErrorCode ierr;
36770034214SMatthew G. Knepley 
36870034214SMatthew G. Knepley   PetscFunctionBeginHot;
36970034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
37070034214SMatthew G. Knepley   PetscValidPointer(adjSize,3);
37170034214SMatthew G. Knepley   PetscValidPointer(adj,4);
3725b317d89SToby Isaac   ierr = DMPlexGetAdjacency_Internal(dm, p, mesh->useCone, mesh->useClosure, mesh->useAnchors, adjSize, adj);CHKERRQ(ierr);
37370034214SMatthew G. Knepley   PetscFunctionReturn(0);
37470034214SMatthew G. Knepley }
375b0a623aaSMatthew G. Knepley #undef __FUNCT__
376b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexCreateTwoSidedProcessSF"
377b0a623aaSMatthew G. Knepley /*@
378b0a623aaSMatthew G. Knepley   DMPlexCreateTwoSidedProcessSF - Create an SF which just has process connectivity
379b0a623aaSMatthew G. Knepley 
380b0a623aaSMatthew G. Knepley   Collective on DM
381b0a623aaSMatthew G. Knepley 
382b0a623aaSMatthew G. Knepley   Input Parameters:
383b0a623aaSMatthew G. Knepley + dm      - The DM
384b0a623aaSMatthew G. Knepley - sfPoint - The PetscSF which encodes point connectivity
385b0a623aaSMatthew G. Knepley 
386b0a623aaSMatthew G. Knepley   Output Parameters:
387b0a623aaSMatthew G. Knepley + processRanks - A list of process neighbors, or NULL
388b0a623aaSMatthew G. Knepley - sfProcess    - An SF encoding the two-sided process connectivity, or NULL
389b0a623aaSMatthew G. Knepley 
390b0a623aaSMatthew G. Knepley   Level: developer
391b0a623aaSMatthew G. Knepley 
392b0a623aaSMatthew G. Knepley .seealso: PetscSFCreate(), DMPlexCreateProcessSF()
393b0a623aaSMatthew G. Knepley @*/
394b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM dm, PetscSF sfPoint, PetscSection rootRankSection, IS rootRanks, PetscSection leafRankSection, IS leafRanks, IS *processRanks, PetscSF *sfProcess)
395b0a623aaSMatthew G. Knepley {
396b0a623aaSMatthew G. Knepley   const PetscSFNode *remotePoints;
397b0a623aaSMatthew G. Knepley   PetscInt          *localPointsNew;
398b0a623aaSMatthew G. Knepley   PetscSFNode       *remotePointsNew;
399b0a623aaSMatthew G. Knepley   const PetscInt    *nranks;
400b0a623aaSMatthew G. Knepley   PetscInt          *ranksNew;
401b0a623aaSMatthew G. Knepley   PetscBT            neighbors;
402b0a623aaSMatthew G. Knepley   PetscInt           pStart, pEnd, p, numLeaves, l, numNeighbors, n;
403b0a623aaSMatthew G. Knepley   PetscMPIInt        numProcs, proc, rank;
404b0a623aaSMatthew G. Knepley   PetscErrorCode     ierr;
405b0a623aaSMatthew G. Knepley 
406b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
407b0a623aaSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
408b0a623aaSMatthew G. Knepley   PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2);
409b0a623aaSMatthew G. Knepley   if (processRanks) {PetscValidPointer(processRanks, 3);}
410b0a623aaSMatthew G. Knepley   if (sfProcess)    {PetscValidPointer(sfProcess, 4);}
411b0a623aaSMatthew G. Knepley   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr);
412b0a623aaSMatthew G. Knepley   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
413b0a623aaSMatthew G. Knepley   ierr = PetscSFGetGraph(sfPoint, NULL, &numLeaves, NULL, &remotePoints);CHKERRQ(ierr);
414b0a623aaSMatthew G. Knepley   ierr = PetscBTCreate(numProcs, &neighbors);CHKERRQ(ierr);
415b0a623aaSMatthew G. Knepley   ierr = PetscBTMemzero(numProcs, neighbors);CHKERRQ(ierr);
416b0a623aaSMatthew G. Knepley   /* Compute root-to-leaf process connectivity */
417b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetChart(rootRankSection, &pStart, &pEnd);CHKERRQ(ierr);
418b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(rootRanks, &nranks);CHKERRQ(ierr);
419b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
420b0a623aaSMatthew G. Knepley     PetscInt ndof, noff, n;
421b0a623aaSMatthew G. Knepley 
422b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetDof(rootRankSection, p, &ndof);CHKERRQ(ierr);
423b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(rootRankSection, p, &noff);CHKERRQ(ierr);
424b0a623aaSMatthew G. Knepley     for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);}
425b0a623aaSMatthew G. Knepley   }
426b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(rootRanks, &nranks);CHKERRQ(ierr);
427b0a623aaSMatthew G. Knepley   /* Compute leaf-to-neighbor process connectivity */
428b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetChart(leafRankSection, &pStart, &pEnd);CHKERRQ(ierr);
429b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(leafRanks, &nranks);CHKERRQ(ierr);
430b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
431b0a623aaSMatthew G. Knepley     PetscInt ndof, noff, n;
432b0a623aaSMatthew G. Knepley 
433b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetDof(leafRankSection, p, &ndof);CHKERRQ(ierr);
434b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(leafRankSection, p, &noff);CHKERRQ(ierr);
435b0a623aaSMatthew G. Knepley     for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);}
436b0a623aaSMatthew G. Knepley   }
437b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(leafRanks, &nranks);CHKERRQ(ierr);
438b0a623aaSMatthew G. Knepley   /* Compute leaf-to-root process connectivity */
439b0a623aaSMatthew G. Knepley   for (l = 0; l < numLeaves; ++l) {PetscBTSet(neighbors, remotePoints[l].rank);}
440b0a623aaSMatthew G. Knepley   /* Calculate edges */
441b0a623aaSMatthew G. Knepley   PetscBTClear(neighbors, rank);
442b0a623aaSMatthew G. Knepley   for(proc = 0, numNeighbors = 0; proc < numProcs; ++proc) {if (PetscBTLookup(neighbors, proc)) ++numNeighbors;}
443b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(numNeighbors, &ranksNew);CHKERRQ(ierr);
444b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(numNeighbors, &localPointsNew);CHKERRQ(ierr);
445b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(numNeighbors, &remotePointsNew);CHKERRQ(ierr);
446b0a623aaSMatthew G. Knepley   for(proc = 0, n = 0; proc < numProcs; ++proc) {
447b0a623aaSMatthew G. Knepley     if (PetscBTLookup(neighbors, proc)) {
448b0a623aaSMatthew G. Knepley       ranksNew[n]              = proc;
449b0a623aaSMatthew G. Knepley       localPointsNew[n]        = proc;
450b0a623aaSMatthew G. Knepley       remotePointsNew[n].index = rank;
451b0a623aaSMatthew G. Knepley       remotePointsNew[n].rank  = proc;
452b0a623aaSMatthew G. Knepley       ++n;
453b0a623aaSMatthew G. Knepley     }
454b0a623aaSMatthew G. Knepley   }
455b0a623aaSMatthew G. Knepley   ierr = PetscBTDestroy(&neighbors);CHKERRQ(ierr);
456b0a623aaSMatthew G. Knepley   if (processRanks) {ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), numNeighbors, ranksNew, PETSC_OWN_POINTER, processRanks);CHKERRQ(ierr);}
457b0a623aaSMatthew G. Knepley   else              {ierr = PetscFree(ranksNew);CHKERRQ(ierr);}
458b0a623aaSMatthew G. Knepley   if (sfProcess) {
459b0a623aaSMatthew G. Knepley     ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess);CHKERRQ(ierr);
460b0a623aaSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) *sfProcess, "Two-Sided Process SF");CHKERRQ(ierr);
461b0a623aaSMatthew G. Knepley     ierr = PetscSFSetFromOptions(*sfProcess);CHKERRQ(ierr);
462b0a623aaSMatthew G. Knepley     ierr = PetscSFSetGraph(*sfProcess, numProcs, numNeighbors, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr);
463b0a623aaSMatthew G. Knepley   }
464b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
465b0a623aaSMatthew G. Knepley }
466b0a623aaSMatthew G. Knepley 
467b0a623aaSMatthew G. Knepley #undef __FUNCT__
468b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeOwnership"
469b0a623aaSMatthew G. Knepley /*@
470b0a623aaSMatthew G. Knepley   DMPlexDistributeOwnership - Compute owner information for shared points. This basically gets two-sided for an SF.
471b0a623aaSMatthew G. Knepley 
472b0a623aaSMatthew G. Knepley   Collective on DM
473b0a623aaSMatthew G. Knepley 
474b0a623aaSMatthew G. Knepley   Input Parameter:
475b0a623aaSMatthew G. Knepley . dm - The DM
476b0a623aaSMatthew G. Knepley 
477b0a623aaSMatthew G. Knepley   Output Parameters:
478b0a623aaSMatthew G. Knepley + rootSection - The number of leaves for a given root point
479b0a623aaSMatthew G. Knepley . rootrank    - The rank of each edge into the root point
480b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point
481b0a623aaSMatthew G. Knepley - leafrank    - The rank of each process sharing a leaf point
482b0a623aaSMatthew G. Knepley 
483b0a623aaSMatthew G. Knepley   Level: developer
484b0a623aaSMatthew G. Knepley 
485b0a623aaSMatthew G. Knepley .seealso: DMPlexCreateOverlap()
486b0a623aaSMatthew G. Knepley @*/
487b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexDistributeOwnership(DM dm, PetscSection rootSection, IS *rootrank, PetscSection leafSection, IS *leafrank)
488b0a623aaSMatthew G. Knepley {
489b0a623aaSMatthew G. Knepley   MPI_Comm        comm;
490b0a623aaSMatthew G. Knepley   PetscSF         sfPoint;
491b0a623aaSMatthew G. Knepley   const PetscInt *rootdegree;
492b0a623aaSMatthew G. Knepley   PetscInt       *myrank, *remoterank;
493b0a623aaSMatthew G. Knepley   PetscInt        pStart, pEnd, p, nedges;
494b0a623aaSMatthew G. Knepley   PetscMPIInt     rank;
495b0a623aaSMatthew G. Knepley   PetscErrorCode  ierr;
496b0a623aaSMatthew G. Knepley 
497b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
498b0a623aaSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
499b0a623aaSMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
500b0a623aaSMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
501b0a623aaSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
502b0a623aaSMatthew G. Knepley   /* Compute number of leaves for each root */
503b0a623aaSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) rootSection, "Root Section");CHKERRQ(ierr);
504b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetChart(rootSection, pStart, pEnd);CHKERRQ(ierr);
505b0a623aaSMatthew G. Knepley   ierr = PetscSFComputeDegreeBegin(sfPoint, &rootdegree);CHKERRQ(ierr);
506b0a623aaSMatthew G. Knepley   ierr = PetscSFComputeDegreeEnd(sfPoint, &rootdegree);CHKERRQ(ierr);
507b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {ierr = PetscSectionSetDof(rootSection, p, rootdegree[p-pStart]);CHKERRQ(ierr);}
508b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
509b0a623aaSMatthew G. Knepley   /* Gather rank of each leaf to root */
510b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(rootSection, &nedges);CHKERRQ(ierr);
511b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(pEnd-pStart, &myrank);CHKERRQ(ierr);
512b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(nedges,  &remoterank);CHKERRQ(ierr);
513b0a623aaSMatthew G. Knepley   for (p = 0; p < pEnd-pStart; ++p) myrank[p] = rank;
514b0a623aaSMatthew G. Knepley   ierr = PetscSFGatherBegin(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr);
515b0a623aaSMatthew G. Knepley   ierr = PetscSFGatherEnd(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr);
516b0a623aaSMatthew G. Knepley   ierr = PetscFree(myrank);CHKERRQ(ierr);
517b0a623aaSMatthew G. Knepley   ierr = ISCreateGeneral(comm, nedges, remoterank, PETSC_OWN_POINTER, rootrank);CHKERRQ(ierr);
518b0a623aaSMatthew G. Knepley   /* Distribute remote ranks to leaves */
519b0a623aaSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) leafSection, "Leaf Section");CHKERRQ(ierr);
520b0a623aaSMatthew G. Knepley   ierr = DMPlexDistributeFieldIS(dm, sfPoint, rootSection, *rootrank, leafSection, leafrank);CHKERRQ(ierr);
521b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
522b0a623aaSMatthew G. Knepley }
523b0a623aaSMatthew G. Knepley 
524b0a623aaSMatthew G. Knepley #undef __FUNCT__
525b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexCreateOverlap"
526278397a0SMatthew G. Knepley /*@C
527b0a623aaSMatthew G. Knepley   DMPlexCreateOverlap - Compute owner information for shared points. This basically gets two-sided for an SF.
528b0a623aaSMatthew G. Knepley 
529b0a623aaSMatthew G. Knepley   Collective on DM
530b0a623aaSMatthew G. Knepley 
531b0a623aaSMatthew G. Knepley   Input Parameters:
532b0a623aaSMatthew G. Knepley + dm          - The DM
533b0a623aaSMatthew G. Knepley . rootSection - The number of leaves for a given root point
534b0a623aaSMatthew G. Knepley . rootrank    - The rank of each edge into the root point
535b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point
536b0a623aaSMatthew G. Knepley - leafrank    - The rank of each process sharing a leaf point
537b0a623aaSMatthew G. Knepley 
538b0a623aaSMatthew G. Knepley   Output Parameters:
539b0a623aaSMatthew G. Knepley + ovRootSection - The number of new overlap points for each neighboring process
540b0a623aaSMatthew G. Knepley . ovRootPoints  - The new overlap points for each neighboring process
541b0a623aaSMatthew G. Knepley . ovLeafSection - The number of new overlap points from each neighboring process
542b0a623aaSMatthew G. Knepley - ovLeafPoints  - The new overlap points from each neighboring process
543b0a623aaSMatthew G. Knepley 
544b0a623aaSMatthew G. Knepley   Level: developer
545b0a623aaSMatthew G. Knepley 
546b0a623aaSMatthew G. Knepley .seealso: DMPlexDistributeOwnership()
547b0a623aaSMatthew G. Knepley @*/
548e540f424SMichael Lange PetscErrorCode DMPlexCreateOverlap(DM dm, PetscSection rootSection, IS rootrank, PetscSection leafSection, IS leafrank, PetscSF *overlapSF)
549b0a623aaSMatthew G. Knepley {
550e540f424SMichael Lange   MPI_Comm           comm;
551b0a623aaSMatthew G. Knepley   DMLabel            ovAdjByRank; /* A DMLabel containing all points adjacent to shared points, separated by rank (value in label) */
552b0a623aaSMatthew G. Knepley   PetscSF            sfPoint, sfProc;
553b0a623aaSMatthew G. Knepley   IS                 valueIS;
554e540f424SMichael Lange   DMLabel            ovLeafLabel;
555b0a623aaSMatthew G. Knepley   const PetscSFNode *remote;
556b0a623aaSMatthew G. Knepley   const PetscInt    *local;
557b0a623aaSMatthew G. Knepley   const PetscInt    *nrank, *rrank, *neighbors;
558e540f424SMichael Lange   PetscSFNode       *ovRootPoints, *ovLeafPoints, *remotePoints;
559e540f424SMichael Lange   PetscSection       ovRootSection, ovLeafSection;
560b0a623aaSMatthew G. Knepley   PetscInt          *adj = NULL;
561b0a623aaSMatthew G. Knepley   PetscInt           pStart, pEnd, p, sStart, sEnd, nleaves, l, numNeighbors, n, ovSize;
562e540f424SMichael Lange   PetscInt           idx, numRemote;
563b0a623aaSMatthew G. Knepley   PetscMPIInt        rank, numProcs;
56426a7d390SMatthew G. Knepley   PetscBool          useCone, useClosure, flg;
565b0a623aaSMatthew G. Knepley   PetscErrorCode     ierr;
566b0a623aaSMatthew G. Knepley 
567b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
568e540f424SMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
569e540f424SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
570e540f424SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
571b0a623aaSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
572b0a623aaSMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
573b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetChart(leafSection, &sStart, &sEnd);CHKERRQ(ierr);
574b0a623aaSMatthew G. Knepley   ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr);
575b0a623aaSMatthew G. Knepley   ierr = DMLabelCreate("Overlap adjacency", &ovAdjByRank);CHKERRQ(ierr);
576b0a623aaSMatthew G. Knepley   /* Handle leaves: shared with the root point */
577b0a623aaSMatthew G. Knepley   for (l = 0; l < nleaves; ++l) {
578b0a623aaSMatthew G. Knepley     PetscInt adjSize = PETSC_DETERMINE, a;
579b0a623aaSMatthew G. Knepley 
580b0a623aaSMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, local[l], &adjSize, &adj);CHKERRQ(ierr);
581b0a623aaSMatthew G. Knepley     for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remote[l].rank);CHKERRQ(ierr);}
582b0a623aaSMatthew G. Knepley   }
583b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(rootrank, &rrank);CHKERRQ(ierr);
584b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(leafrank, &nrank);CHKERRQ(ierr);
585b0a623aaSMatthew G. Knepley   /* Handle roots */
586b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
587b0a623aaSMatthew G. Knepley     PetscInt adjSize = PETSC_DETERMINE, neighbors = 0, noff, n, a;
588b0a623aaSMatthew G. Knepley 
589b0a623aaSMatthew G. Knepley     if ((p >= sStart) && (p < sEnd)) {
590b0a623aaSMatthew G. Knepley       /* Some leaves share a root with other leaves on different processes */
591b0a623aaSMatthew G. Knepley       ierr = PetscSectionGetDof(leafSection, p, &neighbors);CHKERRQ(ierr);
592b0a623aaSMatthew G. Knepley       if (neighbors) {
593b0a623aaSMatthew G. Knepley         ierr = PetscSectionGetOffset(leafSection, p, &noff);CHKERRQ(ierr);
594b0a623aaSMatthew G. Knepley         ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
595b0a623aaSMatthew G. Knepley         for (n = 0; n < neighbors; ++n) {
596b0a623aaSMatthew G. Knepley           const PetscInt remoteRank = nrank[noff+n];
597b0a623aaSMatthew G. Knepley 
598b0a623aaSMatthew G. Knepley           if (remoteRank == rank) continue;
599b0a623aaSMatthew G. Knepley           for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);}
600b0a623aaSMatthew G. Knepley         }
601b0a623aaSMatthew G. Knepley       }
602b0a623aaSMatthew G. Knepley     }
603b0a623aaSMatthew G. Knepley     /* Roots are shared with leaves */
604b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetDof(rootSection, p, &neighbors);CHKERRQ(ierr);
605b0a623aaSMatthew G. Knepley     if (!neighbors) continue;
606b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(rootSection, p, &noff);CHKERRQ(ierr);
607b0a623aaSMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
608b0a623aaSMatthew G. Knepley     for (n = 0; n < neighbors; ++n) {
609b0a623aaSMatthew G. Knepley       const PetscInt remoteRank = rrank[noff+n];
610b0a623aaSMatthew G. Knepley 
611b0a623aaSMatthew G. Knepley       if (remoteRank == rank) continue;
612b0a623aaSMatthew G. Knepley       for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);}
613b0a623aaSMatthew G. Knepley     }
614b0a623aaSMatthew G. Knepley   }
615b0a623aaSMatthew G. Knepley   ierr = PetscFree(adj);CHKERRQ(ierr);
616b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(rootrank, &rrank);CHKERRQ(ierr);
617b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(leafrank, &nrank);CHKERRQ(ierr);
61826a7d390SMatthew G. Knepley   /* We require the closure in the overlap */
61926a7d390SMatthew G. Knepley   ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
62026a7d390SMatthew G. Knepley   ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr);
62126a7d390SMatthew G. Knepley   if (useCone || !useClosure) {
62226a7d390SMatthew G. Knepley     IS              rankIS,   pointIS;
62326a7d390SMatthew G. Knepley     const PetscInt *ranks,   *points;
62426a7d390SMatthew G. Knepley     PetscInt        numRanks, numPoints, r, p;
62526a7d390SMatthew G. Knepley 
62626a7d390SMatthew G. Knepley     ierr = DMLabelGetValueIS(ovAdjByRank, &rankIS);CHKERRQ(ierr);
62726a7d390SMatthew G. Knepley     ierr = ISGetLocalSize(rankIS, &numRanks);CHKERRQ(ierr);
62826a7d390SMatthew G. Knepley     ierr = ISGetIndices(rankIS, &ranks);CHKERRQ(ierr);
62926a7d390SMatthew G. Knepley     for (r = 0; r < numRanks; ++r) {
63026a7d390SMatthew G. Knepley       const PetscInt rank = ranks[r];
63126a7d390SMatthew G. Knepley 
63226a7d390SMatthew G. Knepley       ierr = DMLabelGetStratumIS(ovAdjByRank, rank, &pointIS);CHKERRQ(ierr);
63326a7d390SMatthew G. Knepley       ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
63426a7d390SMatthew G. Knepley       ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
63526a7d390SMatthew G. Knepley       for (p = 0; p < numPoints; ++p) {
63626a7d390SMatthew G. Knepley         PetscInt *closure = NULL, closureSize, c;
63726a7d390SMatthew G. Knepley 
63826a7d390SMatthew G. Knepley         ierr = DMPlexGetTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
63926a7d390SMatthew G. Knepley         for (c = 0; c < closureSize*2; c += 2) {ierr = DMLabelSetValue(ovAdjByRank, closure[c], rank);CHKERRQ(ierr);}
64026a7d390SMatthew G. Knepley         ierr = DMPlexRestoreTransitiveClosure(dm, points[p], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
64126a7d390SMatthew G. Knepley       }
64226a7d390SMatthew G. Knepley       ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
64326a7d390SMatthew G. Knepley       ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
64426a7d390SMatthew G. Knepley     }
64526a7d390SMatthew G. Knepley     ierr = ISRestoreIndices(rankIS, &ranks);CHKERRQ(ierr);
64626a7d390SMatthew G. Knepley     ierr = ISDestroy(&rankIS);CHKERRQ(ierr);
64726a7d390SMatthew G. Knepley   }
648e540f424SMichael Lange   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-overlap_view", &flg);CHKERRQ(ierr);
649e540f424SMichael Lange   if (flg) {
650b0a623aaSMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr);
651b0a623aaSMatthew G. Knepley     ierr = DMLabelView(ovAdjByRank, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
652b0a623aaSMatthew G. Knepley   }
653b0a623aaSMatthew G. Knepley   /* Convert to (point, rank) and use actual owners */
654e540f424SMichael Lange   ierr = PetscSectionCreate(comm, &ovRootSection);CHKERRQ(ierr);
655e540f424SMichael Lange   ierr = PetscObjectSetName((PetscObject) ovRootSection, "Overlap Root Section");CHKERRQ(ierr);
656b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetChart(ovRootSection, 0, numProcs);CHKERRQ(ierr);
657b0a623aaSMatthew G. Knepley   ierr = DMLabelGetValueIS(ovAdjByRank, &valueIS);CHKERRQ(ierr);
658b0a623aaSMatthew G. Knepley   ierr = ISGetLocalSize(valueIS, &numNeighbors);CHKERRQ(ierr);
659b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(valueIS, &neighbors);CHKERRQ(ierr);
660b0a623aaSMatthew G. Knepley   for (n = 0; n < numNeighbors; ++n) {
661b0a623aaSMatthew G. Knepley     PetscInt numPoints;
662b0a623aaSMatthew G. Knepley 
663b0a623aaSMatthew G. Knepley     ierr = DMLabelGetStratumSize(ovAdjByRank, neighbors[n], &numPoints);CHKERRQ(ierr);
664b0a623aaSMatthew G. Knepley     ierr = PetscSectionAddDof(ovRootSection, neighbors[n], numPoints);CHKERRQ(ierr);
665b0a623aaSMatthew G. Knepley   }
666b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetUp(ovRootSection);CHKERRQ(ierr);
667b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(ovRootSection, &ovSize);CHKERRQ(ierr);
668e540f424SMichael Lange   ierr = PetscMalloc1(ovSize, &ovRootPoints);CHKERRQ(ierr);
669b0a623aaSMatthew G. Knepley   ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr);
670b0a623aaSMatthew G. Knepley   for (n = 0; n < numNeighbors; ++n) {
671b0a623aaSMatthew G. Knepley     IS              pointIS;
672b0a623aaSMatthew G. Knepley     const PetscInt *points;
673b0a623aaSMatthew G. Knepley     PetscInt        off, numPoints, p;
674b0a623aaSMatthew G. Knepley 
675b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(ovRootSection, neighbors[n], &off);CHKERRQ(ierr);
676b0a623aaSMatthew G. Knepley     ierr = DMLabelGetStratumIS(ovAdjByRank, neighbors[n], &pointIS);CHKERRQ(ierr);
677b0a623aaSMatthew G. Knepley     ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
678b0a623aaSMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
679b0a623aaSMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
680b0a623aaSMatthew G. Knepley       ierr = PetscFindInt(points[p], nleaves, local, &l);CHKERRQ(ierr);
681e540f424SMichael Lange       if (l >= 0) {ovRootPoints[off+p] = remote[l];}
682e540f424SMichael Lange       else        {ovRootPoints[off+p].index = points[p]; ovRootPoints[off+p].rank = rank;}
683b0a623aaSMatthew G. Knepley     }
684b0a623aaSMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
685b0a623aaSMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
686b0a623aaSMatthew G. Knepley   }
687b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(valueIS, &neighbors);CHKERRQ(ierr);
688b0a623aaSMatthew G. Knepley   ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
689b0a623aaSMatthew G. Knepley   ierr = DMLabelDestroy(&ovAdjByRank);CHKERRQ(ierr);
690b0a623aaSMatthew G. Knepley   /* Make process SF */
691b0a623aaSMatthew G. Knepley   ierr = DMPlexCreateTwoSidedProcessSF(dm, sfPoint, rootSection, rootrank, leafSection, leafrank, NULL, &sfProc);CHKERRQ(ierr);
692e540f424SMichael Lange   if (flg) {
693b0a623aaSMatthew G. Knepley     ierr = PetscSFView(sfProc, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
694b0a623aaSMatthew G. Knepley   }
695b0a623aaSMatthew G. Knepley   /* Communicate overlap */
696e540f424SMichael Lange   ierr = PetscSectionCreate(comm, &ovLeafSection);CHKERRQ(ierr);
697e540f424SMichael Lange   ierr = PetscObjectSetName((PetscObject) ovLeafSection, "Overlap Leaf Section");CHKERRQ(ierr);
698e540f424SMichael Lange   ierr = DMPlexDistributeData(dm, sfProc, ovRootSection, MPIU_2INT, ovRootPoints, ovLeafSection, (void**) &ovLeafPoints);CHKERRQ(ierr);
699e540f424SMichael Lange   /* Filter remote contributions (ovLeafPoints) into the overlapSF */
700e540f424SMichael Lange   ierr = DMLabelCreate("ovLeafs", &ovLeafLabel);CHKERRQ(ierr);
701e540f424SMichael Lange   ierr = PetscSectionGetStorageSize(ovLeafSection, &ovSize);CHKERRQ(ierr);
702e540f424SMichael Lange   for (p = 0; p < ovSize; p++) {
703e540f424SMichael Lange     /* Don't import points from yourself */
704e540f424SMichael Lange     if (ovLeafPoints[p].rank == rank) continue;
705e540f424SMichael Lange     ierr = DMLabelSetValue(ovLeafLabel, ovLeafPoints[p].index, ovLeafPoints[p].rank);CHKERRQ(ierr);
706e540f424SMichael Lange   }
707e540f424SMichael Lange   /* Don't import points already in the pointSF */
708e540f424SMichael Lange   for (l = 0; l < nleaves; ++l) {
709e540f424SMichael Lange     ierr = DMLabelClearValue(ovLeafLabel, remote[l].index, remote[l].rank);CHKERRQ(ierr);
710e540f424SMichael Lange   }
711e540f424SMichael Lange   for (numRemote = 0, n = 0; n < numProcs; ++n) {
712e540f424SMichael Lange     PetscInt numPoints;
713e540f424SMichael Lange     ierr = DMLabelGetStratumSize(ovLeafLabel, n, &numPoints);CHKERRQ(ierr);
714e540f424SMichael Lange     numRemote += numPoints;
715e540f424SMichael Lange   }
716e540f424SMichael Lange   ierr = PetscMalloc1(numRemote, &remotePoints);CHKERRQ(ierr);
717e540f424SMichael Lange   for (idx = 0, n = 0; n < numProcs; ++n) {
718e540f424SMichael Lange     IS remoteRootIS;
719e540f424SMichael Lange     PetscInt numPoints;
720e540f424SMichael Lange     const PetscInt *remoteRoots;
721e540f424SMichael Lange     ierr = DMLabelGetStratumSize(ovLeafLabel, n, &numPoints);CHKERRQ(ierr);
722e540f424SMichael Lange     if (numPoints <= 0) continue;
723e540f424SMichael Lange     ierr = DMLabelGetStratumIS(ovLeafLabel, n, &remoteRootIS);CHKERRQ(ierr);
724e540f424SMichael Lange     ierr = ISGetIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
725e540f424SMichael Lange     for (p = 0; p < numPoints; p++) {
726e540f424SMichael Lange       remotePoints[idx].index = remoteRoots[p];
727e540f424SMichael Lange       remotePoints[idx].rank = n;
728e540f424SMichael Lange       idx++;
729e540f424SMichael Lange     }
730e540f424SMichael Lange     ierr = ISRestoreIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
73115fff7beSMatthew G. Knepley     ierr = ISDestroy(&remoteRootIS);CHKERRQ(ierr);
732e540f424SMichael Lange   }
73315fff7beSMatthew G. Knepley   ierr = DMLabelDestroy(&ovLeafLabel);CHKERRQ(ierr);
734e540f424SMichael Lange   ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), overlapSF);CHKERRQ(ierr);
735e540f424SMichael Lange   ierr = PetscObjectSetName((PetscObject) *overlapSF, "Overlap SF");CHKERRQ(ierr);
736e540f424SMichael Lange   ierr = PetscSFSetFromOptions(*overlapSF);CHKERRQ(ierr);
737e540f424SMichael Lange   ierr = PetscSFSetGraph(*overlapSF, pEnd-pStart, numRemote, NULL, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
738e540f424SMichael Lange   if (flg) {
739e540f424SMichael Lange     ierr = PetscPrintf(comm, "Overlap SF\n");CHKERRQ(ierr);
740e540f424SMichael Lange     ierr = PetscSFView(*overlapSF, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
741e540f424SMichael Lange   }
742e540f424SMichael Lange   /* Clean up */
743b0a623aaSMatthew G. Knepley   ierr = PetscSFDestroy(&sfProc);CHKERRQ(ierr);
744e540f424SMichael Lange   ierr = PetscFree(ovRootPoints);CHKERRQ(ierr);
745e540f424SMichael Lange   ierr = PetscSectionDestroy(&ovRootSection);CHKERRQ(ierr);
746e540f424SMichael Lange   ierr = PetscFree(ovLeafPoints);CHKERRQ(ierr);
747e540f424SMichael Lange   ierr = PetscSectionDestroy(&ovLeafSection);CHKERRQ(ierr);
748b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
749b0a623aaSMatthew G. Knepley }
75070034214SMatthew G. Knepley 
75170034214SMatthew G. Knepley #undef __FUNCT__
75246f9b1c3SMichael Lange #define __FUNCT__ "DMPlexCreateOverlapMigrationSF"
75346f9b1c3SMichael Lange PetscErrorCode DMPlexCreateOverlapMigrationSF(DM dm, PetscSF overlapSF, PetscSF *migrationSF)
75446f9b1c3SMichael Lange {
75546f9b1c3SMichael Lange   MPI_Comm           comm;
75646f9b1c3SMichael Lange   PetscMPIInt        rank, numProcs;
75746f9b1c3SMichael Lange   PetscInt           d, dim, p, pStart, pEnd, nroots, nleaves, newLeaves, point, numSharedPoints;
75846f9b1c3SMichael Lange   PetscInt          *pointDepths, *remoteDepths, *ilocal;
75946f9b1c3SMichael Lange   PetscInt          *depthRecv, *depthShift, *depthIdx;
76046f9b1c3SMichael Lange   PetscSFNode       *iremote;
76146f9b1c3SMichael Lange   PetscSF            pointSF;
76246f9b1c3SMichael Lange   const PetscInt    *sharedLocal;
76346f9b1c3SMichael Lange   const PetscSFNode *overlapRemote, *sharedRemote;
76446f9b1c3SMichael Lange   PetscErrorCode     ierr;
76546f9b1c3SMichael Lange 
76646f9b1c3SMichael Lange   PetscFunctionBegin;
76746f9b1c3SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
76846f9b1c3SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
76946f9b1c3SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
77046f9b1c3SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
77146f9b1c3SMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
77246f9b1c3SMichael Lange 
77346f9b1c3SMichael Lange   /* Before building the migration SF we need to know the new stratum offsets */
77446f9b1c3SMichael Lange   ierr = PetscSFGetGraph(overlapSF, &nroots, &nleaves, NULL, &overlapRemote);CHKERRQ(ierr);
77546f9b1c3SMichael Lange   ierr = PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);CHKERRQ(ierr);
77646f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
77746f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
77846f9b1c3SMichael Lange     for (p=pStart; p<pEnd; p++) pointDepths[p] = d;
77946f9b1c3SMichael Lange   }
78046f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) remoteDepths[p] = -1;
78146f9b1c3SMichael Lange   ierr = PetscSFBcastBegin(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
78246f9b1c3SMichael Lange   ierr = PetscSFBcastEnd(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
78346f9b1c3SMichael Lange 
78446f9b1c3SMichael Lange   /* Count recevied points in each stratum and compute the internal strata shift */
78546f9b1c3SMichael Lange   ierr = PetscMalloc3(dim+1, &depthRecv, dim+1, &depthShift, dim+1, &depthIdx);CHKERRQ(ierr);
78646f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) depthRecv[d]=0;
78746f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) depthRecv[remoteDepths[p]]++;
78846f9b1c3SMichael Lange   depthShift[dim] = 0;
78946f9b1c3SMichael Lange   for (d=0; d<dim; d++) depthShift[d] = depthRecv[dim];
79046f9b1c3SMichael Lange   for (d=1; d<dim; d++) depthShift[d] += depthRecv[0];
79146f9b1c3SMichael Lange   for (d=dim-2; d>0; d--) depthShift[d] += depthRecv[d+1];
79246f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
79346f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
79446f9b1c3SMichael Lange     depthIdx[d] = pStart + depthShift[d];
79546f9b1c3SMichael Lange   }
79646f9b1c3SMichael Lange 
79746f9b1c3SMichael Lange   /* Form the overlap SF build an SF that describes the full overlap migration SF */
79846f9b1c3SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
79946f9b1c3SMichael Lange   newLeaves = pEnd - pStart + nleaves;
80009b7985cSMatthew G. Knepley   ierr = PetscMalloc1(newLeaves, &ilocal);CHKERRQ(ierr);
80109b7985cSMatthew G. Knepley   ierr = PetscMalloc1(newLeaves, &iremote);CHKERRQ(ierr);
80246f9b1c3SMichael Lange   /* First map local points to themselves */
80346f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
80446f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
80546f9b1c3SMichael Lange     for (p=pStart; p<pEnd; p++) {
80646f9b1c3SMichael Lange       point = p + depthShift[d];
80746f9b1c3SMichael Lange       ilocal[point] = point;
80846f9b1c3SMichael Lange       iremote[point].index = p;
80946f9b1c3SMichael Lange       iremote[point].rank = rank;
81046f9b1c3SMichael Lange       depthIdx[d]++;
81146f9b1c3SMichael Lange     }
81246f9b1c3SMichael Lange   }
81346f9b1c3SMichael Lange 
81446f9b1c3SMichael Lange   /* Add in the remote roots for currently shared points */
81546f9b1c3SMichael Lange   ierr = DMGetPointSF(dm, &pointSF);CHKERRQ(ierr);
81646f9b1c3SMichael Lange   ierr = PetscSFGetGraph(pointSF, NULL, &numSharedPoints, &sharedLocal, &sharedRemote);CHKERRQ(ierr);
81746f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
81846f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
81946f9b1c3SMichael Lange     for (p=0; p<numSharedPoints; p++) {
82046f9b1c3SMichael Lange       if (pStart <= sharedLocal[p] && sharedLocal[p] < pEnd) {
82146f9b1c3SMichael Lange         point = sharedLocal[p] + depthShift[d];
82246f9b1c3SMichael Lange         iremote[point].index = sharedRemote[p].index;
82346f9b1c3SMichael Lange         iremote[point].rank = sharedRemote[p].rank;
82446f9b1c3SMichael Lange       }
82546f9b1c3SMichael Lange     }
82646f9b1c3SMichael Lange   }
82746f9b1c3SMichael Lange 
82846f9b1c3SMichael Lange   /* Now add the incoming overlap points */
82946f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) {
83046f9b1c3SMichael Lange     point = depthIdx[remoteDepths[p]];
83146f9b1c3SMichael Lange     ilocal[point] = point;
83246f9b1c3SMichael Lange     iremote[point].index = overlapRemote[p].index;
83346f9b1c3SMichael Lange     iremote[point].rank = overlapRemote[p].rank;
83446f9b1c3SMichael Lange     depthIdx[remoteDepths[p]]++;
83546f9b1c3SMichael Lange   }
83615fff7beSMatthew G. Knepley   ierr = PetscFree2(pointDepths,remoteDepths);CHKERRQ(ierr);
83746f9b1c3SMichael Lange 
83846f9b1c3SMichael Lange   ierr = PetscSFCreate(comm, migrationSF);CHKERRQ(ierr);
83946f9b1c3SMichael Lange   ierr = PetscObjectSetName((PetscObject) *migrationSF, "Overlap Migration SF");CHKERRQ(ierr);
84046f9b1c3SMichael Lange   ierr = PetscSFSetFromOptions(*migrationSF);CHKERRQ(ierr);
84146f9b1c3SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
84246f9b1c3SMichael Lange   ierr = PetscSFSetGraph(*migrationSF, pEnd-pStart, newLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
84346f9b1c3SMichael Lange 
84446f9b1c3SMichael Lange   ierr = PetscFree3(depthRecv, depthShift, depthIdx);CHKERRQ(ierr);
84570034214SMatthew G. Knepley   PetscFunctionReturn(0);
84670034214SMatthew G. Knepley }
84770034214SMatthew G. Knepley 
84870034214SMatthew G. Knepley #undef __FUNCT__
84970034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeField"
85070034214SMatthew G. Knepley /*@
85170034214SMatthew G. Knepley   DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
85270034214SMatthew G. Knepley 
85370034214SMatthew G. Knepley   Collective on DM
85470034214SMatthew G. Knepley 
85570034214SMatthew G. Knepley   Input Parameters:
85670034214SMatthew G. Knepley + dm - The DMPlex object
85770034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
85870034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
85970034214SMatthew G. Knepley - originalVec - The existing data
86070034214SMatthew G. Knepley 
86170034214SMatthew G. Knepley   Output Parameters:
86270034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout
86370034214SMatthew G. Knepley - newVec - The new data
86470034214SMatthew G. Knepley 
86570034214SMatthew G. Knepley   Level: developer
86670034214SMatthew G. Knepley 
86770034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeFieldIS(), DMPlexDistributeData()
86870034214SMatthew G. Knepley @*/
86970034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec)
87070034214SMatthew G. Knepley {
87170034214SMatthew G. Knepley   PetscSF        fieldSF;
87270034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
87370034214SMatthew G. Knepley   PetscScalar   *originalValues, *newValues;
87470034214SMatthew G. Knepley   PetscErrorCode ierr;
87570034214SMatthew G. Knepley 
87670034214SMatthew G. Knepley   PetscFunctionBegin;
87770034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
87870034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
87970034214SMatthew G. Knepley 
88070034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
88170034214SMatthew G. Knepley   ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr);
88270034214SMatthew G. Knepley   ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr);
88370034214SMatthew G. Knepley 
88470034214SMatthew G. Knepley   ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr);
88570034214SMatthew G. Knepley   ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr);
88670034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
88770034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
88870034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
88970034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
89070034214SMatthew G. Knepley   ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr);
89170034214SMatthew G. Knepley   ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr);
89270034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
89370034214SMatthew G. Knepley   PetscFunctionReturn(0);
89470034214SMatthew G. Knepley }
89570034214SMatthew G. Knepley 
89670034214SMatthew G. Knepley #undef __FUNCT__
89770034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeFieldIS"
89870034214SMatthew G. Knepley /*@
89970034214SMatthew G. Knepley   DMPlexDistributeFieldIS - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
90070034214SMatthew G. Knepley 
90170034214SMatthew G. Knepley   Collective on DM
90270034214SMatthew G. Knepley 
90370034214SMatthew G. Knepley   Input Parameters:
90470034214SMatthew G. Knepley + dm - The DMPlex object
90570034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
90670034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
90770034214SMatthew G. Knepley - originalIS - The existing data
90870034214SMatthew G. Knepley 
90970034214SMatthew G. Knepley   Output Parameters:
91070034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout
91170034214SMatthew G. Knepley - newIS - The new data
91270034214SMatthew G. Knepley 
91370034214SMatthew G. Knepley   Level: developer
91470034214SMatthew G. Knepley 
91570034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexDistributeData()
91670034214SMatthew G. Knepley @*/
91770034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeFieldIS(DM dm, PetscSF pointSF, PetscSection originalSection, IS originalIS, PetscSection newSection, IS *newIS)
91870034214SMatthew G. Knepley {
91970034214SMatthew G. Knepley   PetscSF         fieldSF;
92070034214SMatthew G. Knepley   PetscInt       *newValues, *remoteOffsets, fieldSize;
92170034214SMatthew G. Knepley   const PetscInt *originalValues;
92270034214SMatthew G. Knepley   PetscErrorCode  ierr;
92370034214SMatthew G. Knepley 
92470034214SMatthew G. Knepley   PetscFunctionBegin;
92570034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
92670034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
92770034214SMatthew G. Knepley 
92870034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
929*854ce69bSBarry Smith   ierr = PetscMalloc1(fieldSize, &newValues);CHKERRQ(ierr);
93070034214SMatthew G. Knepley 
93170034214SMatthew G. Knepley   ierr = ISGetIndices(originalIS, &originalValues);CHKERRQ(ierr);
93270034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
933aaf8c182SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr);
934aaf8c182SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr);
93570034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
93670034214SMatthew G. Knepley   ierr = ISRestoreIndices(originalIS, &originalValues);CHKERRQ(ierr);
93770034214SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) pointSF), fieldSize, newValues, PETSC_OWN_POINTER, newIS);CHKERRQ(ierr);
93870034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
93970034214SMatthew G. Knepley   PetscFunctionReturn(0);
94070034214SMatthew G. Knepley }
94170034214SMatthew G. Knepley 
94270034214SMatthew G. Knepley #undef __FUNCT__
94370034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeData"
94470034214SMatthew G. Knepley /*@
94570034214SMatthew G. Knepley   DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
94670034214SMatthew G. Knepley 
94770034214SMatthew G. Knepley   Collective on DM
94870034214SMatthew G. Knepley 
94970034214SMatthew G. Knepley   Input Parameters:
95070034214SMatthew G. Knepley + dm - The DMPlex object
95170034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
95270034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
95370034214SMatthew G. Knepley . datatype - The type of data
95470034214SMatthew G. Knepley - originalData - The existing data
95570034214SMatthew G. Knepley 
95670034214SMatthew G. Knepley   Output Parameters:
95770034214SMatthew G. Knepley + newSection - The PetscSection describing the new data layout
95870034214SMatthew G. Knepley - newData - The new data
95970034214SMatthew G. Knepley 
96070034214SMatthew G. Knepley   Level: developer
96170034214SMatthew G. Knepley 
96270034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField()
96370034214SMatthew G. Knepley @*/
96470034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData)
96570034214SMatthew G. Knepley {
96670034214SMatthew G. Knepley   PetscSF        fieldSF;
96770034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
96870034214SMatthew G. Knepley   PetscMPIInt    dataSize;
96970034214SMatthew G. Knepley   PetscErrorCode ierr;
97070034214SMatthew G. Knepley 
97170034214SMatthew G. Knepley   PetscFunctionBegin;
97270034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
97370034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
97470034214SMatthew G. Knepley 
97570034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
97670034214SMatthew G. Knepley   ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr);
97770034214SMatthew G. Knepley   ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr);
97870034214SMatthew G. Knepley 
97970034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
98070034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
98170034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
98270034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
98370034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
98470034214SMatthew G. Knepley   PetscFunctionReturn(0);
98570034214SMatthew G. Knepley }
98670034214SMatthew G. Knepley 
98770034214SMatthew G. Knepley #undef __FUNCT__
988cc71bff1SMichael Lange #define __FUNCT__ "DMPlexDistributeCones"
989cc71bff1SMichael Lange PetscErrorCode DMPlexDistributeCones(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
990cc71bff1SMichael Lange {
991cc71bff1SMichael Lange   DM_Plex               *pmesh   = (DM_Plex*) (dmParallel)->data;
992cc71bff1SMichael Lange   MPI_Comm               comm;
993cc71bff1SMichael Lange   PetscSF                coneSF;
994cc71bff1SMichael Lange   PetscSection           originalConeSection, newConeSection;
995cc71bff1SMichael Lange   PetscInt              *remoteOffsets, *cones, *newCones, newConesSize;
996cc71bff1SMichael Lange   PetscBool              flg;
997cc71bff1SMichael Lange   PetscErrorCode         ierr;
998cc71bff1SMichael Lange 
999cc71bff1SMichael Lange   PetscFunctionBegin;
1000cc71bff1SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1001cc71bff1SMichael Lange   PetscValidPointer(dmParallel,4);
1002cc71bff1SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
1003cc71bff1SMichael Lange 
1004cc71bff1SMichael Lange   /* Distribute cone section */
1005cc71bff1SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
1006cc71bff1SMichael Lange   ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr);
1007cc71bff1SMichael Lange   ierr = DMPlexGetConeSection(dmParallel, &newConeSection);CHKERRQ(ierr);
1008cc71bff1SMichael Lange   ierr = PetscSFDistributeSection(migrationSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr);
1009cc71bff1SMichael Lange   ierr = DMSetUp(dmParallel);CHKERRQ(ierr);
1010cc71bff1SMichael Lange   {
1011cc71bff1SMichael Lange     PetscInt pStart, pEnd, p;
1012cc71bff1SMichael Lange 
1013cc71bff1SMichael Lange     ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr);
1014cc71bff1SMichael Lange     for (p = pStart; p < pEnd; ++p) {
1015cc71bff1SMichael Lange       PetscInt coneSize;
1016cc71bff1SMichael Lange       ierr               = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr);
1017cc71bff1SMichael Lange       pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize);
1018cc71bff1SMichael Lange     }
1019cc71bff1SMichael Lange   }
1020cc71bff1SMichael Lange   /* Communicate and renumber cones */
1021cc71bff1SMichael Lange   ierr = PetscSFCreateSectionSF(migrationSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr);
1022cc71bff1SMichael Lange   ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
1023cc71bff1SMichael Lange   ierr = DMPlexGetCones(dmParallel, &newCones);CHKERRQ(ierr);
1024cc71bff1SMichael Lange   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
1025cc71bff1SMichael Lange   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
1026cc71bff1SMichael Lange   ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr);
1027cc71bff1SMichael Lange   ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr);
10283533c52bSMatthew G. Knepley #if PETSC_USE_DEBUG
10293533c52bSMatthew G. Knepley   {
10303533c52bSMatthew G. Knepley     PetscInt  p;
10313533c52bSMatthew G. Knepley     PetscBool valid = PETSC_TRUE;
10323533c52bSMatthew G. Knepley     for (p = 0; p < newConesSize; ++p) {
10333533c52bSMatthew G. Knepley       if (newCones[p] < 0) {valid = PETSC_FALSE; ierr = PetscPrintf(PETSC_COMM_SELF, "Point %d not in overlap SF\n", p);CHKERRQ(ierr);}
10343533c52bSMatthew G. Knepley     }
10353533c52bSMatthew G. Knepley     if (!valid) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map");
10363533c52bSMatthew G. Knepley   }
10373533c52bSMatthew G. Knepley #endif
1038cc71bff1SMichael Lange   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr);
1039cc71bff1SMichael Lange   if (flg) {
1040cc71bff1SMichael Lange     ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr);
1041cc71bff1SMichael Lange     ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1042cc71bff1SMichael Lange     ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr);
1043cc71bff1SMichael Lange     ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1044cc71bff1SMichael Lange     ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr);
1045cc71bff1SMichael Lange   }
1046cc71bff1SMichael Lange   ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr);
1047cc71bff1SMichael Lange   ierr = DMPlexGetConeOrientations(dmParallel, &newCones);CHKERRQ(ierr);
1048cc71bff1SMichael Lange   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
1049cc71bff1SMichael Lange   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
1050cc71bff1SMichael Lange   ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr);
1051cc71bff1SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
1052cc71bff1SMichael Lange   /* Create supports and stratify sieve */
1053cc71bff1SMichael Lange   {
1054cc71bff1SMichael Lange     PetscInt pStart, pEnd;
1055cc71bff1SMichael Lange 
1056cc71bff1SMichael Lange     ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr);
1057cc71bff1SMichael Lange     ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr);
1058cc71bff1SMichael Lange   }
1059cc71bff1SMichael Lange   ierr = DMPlexSymmetrize(dmParallel);CHKERRQ(ierr);
1060cc71bff1SMichael Lange   ierr = DMPlexStratify(dmParallel);CHKERRQ(ierr);
1061cc71bff1SMichael Lange   PetscFunctionReturn(0);
1062cc71bff1SMichael Lange }
1063cc71bff1SMichael Lange 
1064cc71bff1SMichael Lange #undef __FUNCT__
10650df0e737SMichael Lange #define __FUNCT__ "DMPlexDistributeCoordinates"
10660df0e737SMichael Lange PetscErrorCode DMPlexDistributeCoordinates(DM dm, PetscSF migrationSF, DM dmParallel)
10670df0e737SMichael Lange {
10680df0e737SMichael Lange   MPI_Comm         comm;
10690df0e737SMichael Lange   PetscSection     originalCoordSection, newCoordSection;
10700df0e737SMichael Lange   Vec              originalCoordinates, newCoordinates;
10710df0e737SMichael Lange   PetscInt         bs;
10720df0e737SMichael Lange   const char      *name;
10730df0e737SMichael Lange   const PetscReal *maxCell, *L;
10740df0e737SMichael Lange   PetscErrorCode   ierr;
10750df0e737SMichael Lange 
10760df0e737SMichael Lange   PetscFunctionBegin;
10770df0e737SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10780df0e737SMichael Lange   PetscValidPointer(dmParallel, 3);
10790df0e737SMichael Lange 
10800df0e737SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
10810df0e737SMichael Lange   ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr);
10820df0e737SMichael Lange   ierr = DMGetCoordinateSection(dmParallel, &newCoordSection);CHKERRQ(ierr);
10830df0e737SMichael Lange   ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr);
10840df0e737SMichael Lange   if (originalCoordinates) {
10850df0e737SMichael Lange     ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr);
10860df0e737SMichael Lange     ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr);
10870df0e737SMichael Lange     ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr);
10880df0e737SMichael Lange 
10890df0e737SMichael Lange     ierr = DMPlexDistributeField(dm, migrationSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr);
10900df0e737SMichael Lange     ierr = DMSetCoordinatesLocal(dmParallel, newCoordinates);CHKERRQ(ierr);
10910df0e737SMichael Lange     ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr);
10920df0e737SMichael Lange     ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr);
10930df0e737SMichael Lange     ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr);
10940df0e737SMichael Lange   }
10950df0e737SMichael Lange   ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr);
10960df0e737SMichael Lange   if (L) {ierr = DMSetPeriodicity(dmParallel, maxCell, L);CHKERRQ(ierr);}
10970df0e737SMichael Lange   PetscFunctionReturn(0);
10980df0e737SMichael Lange }
10990df0e737SMichael Lange 
11000df0e737SMichael Lange #undef __FUNCT__
11010df0e737SMichael Lange #define __FUNCT__ "DMPlexDistributeLabels"
1102d995df53SMatthew G. Knepley /* Here we are assuming that process 0 always has everything */
11032eb1fa14SMichael Lange PetscErrorCode DMPlexDistributeLabels(DM dm, PetscSF migrationSF, DM dmParallel)
11040df0e737SMichael Lange {
11050df0e737SMichael Lange   MPI_Comm       comm;
11060df0e737SMichael Lange   PetscMPIInt    rank;
1107d995df53SMatthew G. Knepley   PetscInt       numLabels, numLocalLabels, l;
1108d995df53SMatthew G. Knepley   PetscBool      hasLabels = PETSC_FALSE;
11090df0e737SMichael Lange   PetscErrorCode ierr;
11100df0e737SMichael Lange 
11110df0e737SMichael Lange   PetscFunctionBegin;
11120df0e737SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11132eb1fa14SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
11140df0e737SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
11150df0e737SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
11160df0e737SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
11170df0e737SMichael Lange 
1118d995df53SMatthew G. Knepley   /* Everyone must have either the same number of labels, or none */
1119d995df53SMatthew G. Knepley   ierr = DMPlexGetNumLabels(dm, &numLocalLabels);CHKERRQ(ierr);
1120d995df53SMatthew G. Knepley   numLabels = numLocalLabels;
11210df0e737SMichael Lange   ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1122627847f0SMatthew G. Knepley   if (numLabels == numLocalLabels) hasLabels = PETSC_TRUE;
1123bdd2d751SMichael Lange   for (l = numLabels-1; l >= 0; --l) {
1124bdd2d751SMichael Lange     DMLabel     label = NULL, labelNew = NULL;
11250df0e737SMichael Lange     PetscBool   isdepth;
11260df0e737SMichael Lange 
1127d995df53SMatthew G. Knepley     if (hasLabels) {
1128bdd2d751SMichael Lange       ierr = DMPlexGetLabelByNum(dm, l, &label);CHKERRQ(ierr);
11290df0e737SMichael Lange       /* Skip "depth" because it is recreated */
1130bdd2d751SMichael Lange       ierr = PetscStrcmp(label->name, "depth", &isdepth);CHKERRQ(ierr);
1131bdd2d751SMichael Lange     }
11320df0e737SMichael Lange     ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
1133bdd2d751SMichael Lange     if (isdepth) continue;
1134bdd2d751SMichael Lange     ierr = DMLabelDistribute(label, migrationSF, &labelNew);CHKERRQ(ierr);
1135bdd2d751SMichael Lange     ierr = DMPlexAddLabel(dmParallel, labelNew);CHKERRQ(ierr);
11360df0e737SMichael Lange   }
11370df0e737SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
11380df0e737SMichael Lange   PetscFunctionReturn(0);
11390df0e737SMichael Lange }
11400df0e737SMichael Lange 
11419734c634SMichael Lange #undef __FUNCT__
11429734c634SMichael Lange #define __FUNCT__ "DMPlexDistributeSetupHybrid"
11439734c634SMichael Lange PetscErrorCode DMPlexDistributeSetupHybrid(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
11449734c634SMichael Lange {
11459734c634SMichael Lange   DM_Plex        *mesh  = (DM_Plex*) dm->data;
11469734c634SMichael Lange   DM_Plex        *pmesh = (DM_Plex*) (dmParallel)->data;
11479734c634SMichael Lange   MPI_Comm        comm;
11489734c634SMichael Lange   const PetscInt *gpoints;
11499734c634SMichael Lange   PetscInt        dim, depth, n, d;
11509734c634SMichael Lange   PetscErrorCode  ierr;
11519734c634SMichael Lange 
11529734c634SMichael Lange   PetscFunctionBegin;
11539734c634SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11549734c634SMichael Lange   PetscValidPointer(dmParallel, 4);
11559734c634SMichael Lange 
11569734c634SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
11579734c634SMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
11589734c634SMichael Lange 
11599734c634SMichael Lange   /* Setup hybrid structure */
11609734c634SMichael Lange   for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];}
11619734c634SMichael Lange   ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr);
11629734c634SMichael Lange   ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr);
11639734c634SMichael Lange   ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr);
11649734c634SMichael Lange   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
11659734c634SMichael Lange   for (d = 0; d <= dim; ++d) {
11669734c634SMichael Lange     PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p;
11679734c634SMichael Lange 
11689734c634SMichael Lange     if (pmax < 0) continue;
11699734c634SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr);
11709734c634SMichael Lange     ierr = DMPlexGetDepthStratum(dmParallel, d, NULL, &pEnd);CHKERRQ(ierr);
11719734c634SMichael Lange     ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr);
11729734c634SMichael Lange     for (p = 0; p < n; ++p) {
11739734c634SMichael Lange       const PetscInt point = gpoints[p];
11749734c634SMichael Lange 
11759734c634SMichael Lange       if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax;
11769734c634SMichael Lange     }
11779734c634SMichael Lange     if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax;
11789734c634SMichael Lange     else            pmesh->hybridPointMax[d] = -1;
11799734c634SMichael Lange   }
11809734c634SMichael Lange   ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr);
11819734c634SMichael Lange   PetscFunctionReturn(0);
11829734c634SMichael Lange }
11830df0e737SMichael Lange 
1184a6f36705SMichael Lange #undef __FUNCT__
1185a6f36705SMichael Lange #define __FUNCT__ "DMPlexDistributeSetupTree"
1186a6f36705SMichael Lange PetscErrorCode DMPlexDistributeSetupTree(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
1187a6f36705SMichael Lange {
1188a6f36705SMichael Lange   MPI_Comm        comm;
1189a6f36705SMichael Lange   DM              refTree;
1190a6f36705SMichael Lange   PetscSection    origParentSection, newParentSection;
1191a6f36705SMichael Lange   PetscInt        *origParents, *origChildIDs;
1192a6f36705SMichael Lange   PetscBool       flg;
1193a6f36705SMichael Lange   PetscErrorCode  ierr;
1194a6f36705SMichael Lange 
1195a6f36705SMichael Lange   PetscFunctionBegin;
1196a6f36705SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1197a6f36705SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 4);
1198a6f36705SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
1199a6f36705SMichael Lange 
1200a6f36705SMichael Lange   /* Set up tree */
1201a6f36705SMichael Lange   ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr);
1202a6f36705SMichael Lange   ierr = DMPlexSetReferenceTree(dmParallel,refTree);CHKERRQ(ierr);
1203a6f36705SMichael Lange   ierr = DMPlexGetTree(dm,&origParentSection,&origParents,&origChildIDs,NULL,NULL);CHKERRQ(ierr);
1204a6f36705SMichael Lange   if (origParentSection) {
1205a6f36705SMichael Lange     PetscInt        pStart, pEnd;
1206a6f36705SMichael Lange     PetscInt        *newParents, *newChildIDs;
1207a6f36705SMichael Lange     PetscInt        *remoteOffsetsParents, newParentSize;
1208a6f36705SMichael Lange     PetscSF         parentSF;
1209a6f36705SMichael Lange 
1210a6f36705SMichael Lange     ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
1211a6f36705SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dmParallel),&newParentSection);CHKERRQ(ierr);
1212a6f36705SMichael Lange     ierr = PetscSectionSetChart(newParentSection,pStart,pEnd);CHKERRQ(ierr);
1213a6f36705SMichael Lange     ierr = PetscSFDistributeSection(migrationSF, origParentSection, &remoteOffsetsParents, newParentSection);CHKERRQ(ierr);
1214a6f36705SMichael Lange     ierr = PetscSFCreateSectionSF(migrationSF, origParentSection, remoteOffsetsParents, newParentSection, &parentSF);CHKERRQ(ierr);
1215a6f36705SMichael Lange     ierr = PetscSectionGetStorageSize(newParentSection,&newParentSize);CHKERRQ(ierr);
1216a6f36705SMichael Lange     ierr = PetscMalloc2(newParentSize,&newParents,newParentSize,&newChildIDs);CHKERRQ(ierr);
1217a6f36705SMichael Lange     ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
1218a6f36705SMichael Lange     ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
1219a6f36705SMichael Lange     ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1220a6f36705SMichael Lange     ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1221a6f36705SMichael Lange     ierr = ISGlobalToLocalMappingApplyBlock(renumbering,IS_GTOLM_MASK, newParentSize, newParents, NULL, newParents);CHKERRQ(ierr);
1222a6f36705SMichael Lange     ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-parents_view", &flg);CHKERRQ(ierr);
1223a6f36705SMichael Lange     if (flg) {
1224a6f36705SMichael Lange       ierr = PetscPrintf(comm, "Serial Parent Section: \n");CHKERRQ(ierr);
1225a6f36705SMichael Lange       ierr = PetscSectionView(origParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1226a6f36705SMichael Lange       ierr = PetscPrintf(comm, "Parallel Parent Section: \n");CHKERRQ(ierr);
1227a6f36705SMichael Lange       ierr = PetscSectionView(newParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1228a6f36705SMichael Lange       ierr = PetscSFView(parentSF, NULL);CHKERRQ(ierr);
1229a6f36705SMichael Lange     }
1230a6f36705SMichael Lange     ierr = DMPlexSetTree(dmParallel,newParentSection,newParents,newChildIDs);CHKERRQ(ierr);
1231a6f36705SMichael Lange     ierr = PetscSectionDestroy(&newParentSection);CHKERRQ(ierr);
1232a6f36705SMichael Lange     ierr = PetscFree2(newParents,newChildIDs);CHKERRQ(ierr);
1233a6f36705SMichael Lange     ierr = PetscSFDestroy(&parentSF);CHKERRQ(ierr);
1234a6f36705SMichael Lange   }
1235a6f36705SMichael Lange   PetscFunctionReturn(0);
1236a6f36705SMichael Lange }
12370df0e737SMichael Lange 
12380df0e737SMichael Lange #undef __FUNCT__
12394eca1733SMichael Lange #define __FUNCT__ "DMPlexDistributeSF"
12404eca1733SMichael Lange PetscErrorCode DMPlexDistributeSF(DM dm, PetscSF migrationSF, PetscSection partSection, IS part, PetscSection origPartSection, IS origPart, DM dmParallel)
12414eca1733SMichael Lange {
12424eca1733SMichael Lange   DM_Plex               *mesh  = (DM_Plex*) dm->data;
12434eca1733SMichael Lange   DM_Plex               *pmesh = (DM_Plex*) (dmParallel)->data;
12444eca1733SMichael Lange   PetscMPIInt            rank, numProcs;
12454eca1733SMichael Lange   MPI_Comm               comm;
12464eca1733SMichael Lange   PetscErrorCode         ierr;
12474eca1733SMichael Lange 
12484eca1733SMichael Lange   PetscFunctionBegin;
12494eca1733SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12504eca1733SMichael Lange   PetscValidPointer(dmParallel,7);
12514eca1733SMichael Lange 
12524eca1733SMichael Lange   /* Create point SF for parallel mesh */
12534eca1733SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
12544eca1733SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
12554eca1733SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
12564eca1733SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
12574eca1733SMichael Lange   {
12584eca1733SMichael Lange     const PetscInt *leaves;
12594eca1733SMichael Lange     PetscSFNode    *remotePoints, *rowners, *lowners;
12604eca1733SMichael Lange     PetscInt        numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints;
12614eca1733SMichael Lange     PetscInt        pStart, pEnd;
12624eca1733SMichael Lange 
12634eca1733SMichael Lange     ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
12644eca1733SMichael Lange     ierr = PetscSFGetGraph(migrationSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr);
12654eca1733SMichael Lange     ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr);
12664eca1733SMichael Lange     for (p=0; p<numRoots; p++) {
12674eca1733SMichael Lange       rowners[p].rank  = -1;
12684eca1733SMichael Lange       rowners[p].index = -1;
12694eca1733SMichael Lange     }
12704eca1733SMichael Lange     if (origPart) {
12714eca1733SMichael Lange       /* Make sure points in the original partition are not assigned to other procs */
12724eca1733SMichael Lange       const PetscInt *origPoints;
12734eca1733SMichael Lange 
12744eca1733SMichael Lange       ierr = ISGetIndices(origPart, &origPoints);CHKERRQ(ierr);
12754eca1733SMichael Lange       for (p = 0; p < numProcs; ++p) {
12764eca1733SMichael Lange         PetscInt dof, off, d;
12774eca1733SMichael Lange 
12784eca1733SMichael Lange         ierr = PetscSectionGetDof(origPartSection, p, &dof);CHKERRQ(ierr);
12794eca1733SMichael Lange         ierr = PetscSectionGetOffset(origPartSection, p, &off);CHKERRQ(ierr);
12804eca1733SMichael Lange         for (d = off; d < off+dof; ++d) {
12814eca1733SMichael Lange           rowners[origPoints[d]].rank = p;
12824eca1733SMichael Lange         }
12834eca1733SMichael Lange       }
12844eca1733SMichael Lange       ierr = ISRestoreIndices(origPart, &origPoints);CHKERRQ(ierr);
12854eca1733SMichael Lange     }
12864eca1733SMichael Lange     ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12874eca1733SMichael Lange     ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12884eca1733SMichael Lange     for (p = 0; p < numLeaves; ++p) {
12894eca1733SMichael Lange       if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */
12904eca1733SMichael Lange         lowners[p].rank  = rank;
12914eca1733SMichael Lange         lowners[p].index = leaves ? leaves[p] : p;
12924eca1733SMichael Lange       } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */
12934eca1733SMichael Lange         lowners[p].rank  = -2;
12944eca1733SMichael Lange         lowners[p].index = -2;
12954eca1733SMichael Lange       }
12964eca1733SMichael Lange     }
12974eca1733SMichael Lange     for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */
12984eca1733SMichael Lange       rowners[p].rank  = -3;
12994eca1733SMichael Lange       rowners[p].index = -3;
13004eca1733SMichael Lange     }
13014eca1733SMichael Lange     ierr = PetscSFReduceBegin(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
13024eca1733SMichael Lange     ierr = PetscSFReduceEnd(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
13034eca1733SMichael Lange     ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
13044eca1733SMichael Lange     ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
13054eca1733SMichael Lange     for (p = 0; p < numLeaves; ++p) {
13064eca1733SMichael Lange       if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed");
13074eca1733SMichael Lange       if (lowners[p].rank != rank) ++numGhostPoints;
13084eca1733SMichael Lange     }
13094eca1733SMichael Lange     ierr = PetscMalloc1(numGhostPoints, &ghostPoints);CHKERRQ(ierr);
13104eca1733SMichael Lange     ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr);
13114eca1733SMichael Lange     for (p = 0, gp = 0; p < numLeaves; ++p) {
13124eca1733SMichael Lange       if (lowners[p].rank != rank) {
13134eca1733SMichael Lange         ghostPoints[gp]        = leaves ? leaves[p] : p;
13144eca1733SMichael Lange         remotePoints[gp].rank  = lowners[p].rank;
13154eca1733SMichael Lange         remotePoints[gp].index = lowners[p].index;
13164eca1733SMichael Lange         ++gp;
13174eca1733SMichael Lange       }
13184eca1733SMichael Lange     }
13194eca1733SMichael Lange     ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr);
13204eca1733SMichael Lange     ierr = PetscSFSetGraph((dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
13214eca1733SMichael Lange     ierr = PetscSFSetFromOptions((dmParallel)->sf);CHKERRQ(ierr);
13224eca1733SMichael Lange   }
13234eca1733SMichael Lange   pmesh->useCone    = mesh->useCone;
13244eca1733SMichael Lange   pmesh->useClosure = mesh->useClosure;
13254eca1733SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
13264eca1733SMichael Lange   PetscFunctionReturn(0);
13274eca1733SMichael Lange }
13284eca1733SMichael Lange 
13294eca1733SMichael Lange #undef __FUNCT__
133070034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistribute"
133156117bffSMatthew G. Knepley /*@
133270034214SMatthew G. Knepley   DMPlexDistribute - Distributes the mesh and any associated sections.
133370034214SMatthew G. Knepley 
133470034214SMatthew G. Knepley   Not Collective
133570034214SMatthew G. Knepley 
133670034214SMatthew G. Knepley   Input Parameter:
133770034214SMatthew G. Knepley + dm  - The original DMPlex object
133870034214SMatthew G. Knepley - overlap - The overlap of partitions, 0 is the default
133970034214SMatthew G. Knepley 
134070034214SMatthew G. Knepley   Output Parameter:
134170034214SMatthew G. Knepley + sf - The PetscSF used for point distribution
134270034214SMatthew G. Knepley - parallelMesh - The distributed DMPlex object, or NULL
134370034214SMatthew G. Knepley 
134470034214SMatthew G. Knepley   Note: If the mesh was not distributed, the return value is NULL.
134570034214SMatthew G. Knepley 
134670034214SMatthew G. Knepley   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
134770034214SMatthew G. Knepley   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
134870034214SMatthew G. Knepley   representation on the mesh.
134970034214SMatthew G. Knepley 
135070034214SMatthew G. Knepley   Level: intermediate
135170034214SMatthew G. Knepley 
135270034214SMatthew G. Knepley .keywords: mesh, elements
135370034214SMatthew G. Knepley .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
135470034214SMatthew G. Knepley @*/
135580cf41d5SMatthew G. Knepley PetscErrorCode DMPlexDistribute(DM dm, PetscInt overlap, PetscSF *sf, DM *dmParallel)
135670034214SMatthew G. Knepley {
135770034214SMatthew G. Knepley   DM_Plex               *mesh   = (DM_Plex*) dm->data, *pmesh;
135870034214SMatthew G. Knepley   MPI_Comm               comm;
135943331d4aSMichael Lange   PetscInt               dim, numRemoteRanks, nroots, nleaves;
1360a157612eSMichael Lange   DM                     dmOverlap;
1361a157612eSMichael Lange   IS                     cellPart,        part;
1362a157612eSMichael Lange   PetscSection           cellPartSection, partSection;
136343331d4aSMichael Lange   PetscSFNode           *remoteRanks, *newRemote;
136443331d4aSMichael Lange   const PetscSFNode     *oldRemote;
136543331d4aSMichael Lange   PetscSF                partSF, pointSF, overlapPointSF, overlapSF;
136670034214SMatthew G. Knepley   ISLocalToGlobalMapping renumbering;
136770034214SMatthew G. Knepley   PetscBool              flg;
136870034214SMatthew G. Knepley   PetscMPIInt            rank, numProcs, p;
136970034214SMatthew G. Knepley   PetscErrorCode         ierr;
137070034214SMatthew G. Knepley 
137170034214SMatthew G. Knepley   PetscFunctionBegin;
137270034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
137370034214SMatthew G. Knepley   if (sf) PetscValidPointer(sf,4);
137470034214SMatthew G. Knepley   PetscValidPointer(dmParallel,5);
137570034214SMatthew G. Knepley 
137670034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
137770034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
137870034214SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
137970034214SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
138070034214SMatthew G. Knepley 
138170034214SMatthew G. Knepley   *dmParallel = NULL;
138270034214SMatthew G. Knepley   if (numProcs == 1) PetscFunctionReturn(0);
138370034214SMatthew G. Knepley 
138470034214SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
138570034214SMatthew G. Knepley   /* Create cell partition - We need to rewrite to use IS, use the MatPartition stuff */
138677623264SMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
138777623264SMatthew G. Knepley   if (overlap > 1) SETERRQ(comm, PETSC_ERR_SUP, "Overlap > 1 not yet implemented");
138877623264SMatthew G. Knepley   ierr = PetscSectionCreate(comm, &cellPartSection);CHKERRQ(ierr);
1389a157612eSMichael Lange   ierr = PetscPartitionerPartition(mesh->partitioner, dm, PETSC_FALSE, cellPartSection, &cellPart, NULL, NULL);CHKERRQ(ierr);
139070034214SMatthew G. Knepley   /* Create SF assuming a serial partition for all processes: Could check for IS length here */
139170034214SMatthew G. Knepley   if (!rank) numRemoteRanks = numProcs;
139270034214SMatthew G. Knepley   else       numRemoteRanks = 0;
139370034214SMatthew G. Knepley   ierr = PetscMalloc1(numRemoteRanks, &remoteRanks);CHKERRQ(ierr);
139470034214SMatthew G. Knepley   for (p = 0; p < numRemoteRanks; ++p) {
139570034214SMatthew G. Knepley     remoteRanks[p].rank  = p;
139670034214SMatthew G. Knepley     remoteRanks[p].index = 0;
139770034214SMatthew G. Knepley   }
139870034214SMatthew G. Knepley   ierr = PetscSFCreate(comm, &partSF);CHKERRQ(ierr);
139970034214SMatthew G. Knepley   ierr = PetscSFSetGraph(partSF, 1, numRemoteRanks, NULL, PETSC_OWN_POINTER, remoteRanks, PETSC_OWN_POINTER);CHKERRQ(ierr);
140070034214SMatthew G. Knepley   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
140170034214SMatthew G. Knepley   if (flg) {
1402a157612eSMichael Lange     ierr = PetscPrintf(comm, "Original Cell Partition:\n");CHKERRQ(ierr);
140370034214SMatthew G. Knepley     ierr = PetscSectionView(cellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
140470034214SMatthew G. Knepley     ierr = ISView(cellPart, NULL);CHKERRQ(ierr);
140570034214SMatthew G. Knepley     ierr = PetscSFView(partSF, NULL);CHKERRQ(ierr);
140670034214SMatthew G. Knepley   }
140770034214SMatthew G. Knepley   /* Close the partition over the mesh */
140870034214SMatthew G. Knepley   ierr = DMPlexCreatePartitionClosure(dm, cellPartSection, cellPart, &partSection, &part);CHKERRQ(ierr);
140970034214SMatthew G. Knepley   /* Create new mesh */
141070034214SMatthew G. Knepley   ierr  = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr);
141170034214SMatthew G. Knepley   ierr  = DMSetDimension(*dmParallel, dim);CHKERRQ(ierr);
141270034214SMatthew G. Knepley   ierr  = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr);
141370034214SMatthew G. Knepley   pmesh = (DM_Plex*) (*dmParallel)->data;
14144eca1733SMichael Lange   pmesh->useAnchors = mesh->useAnchors;
14154eca1733SMichael Lange 
141670034214SMatthew G. Knepley   /* Distribute sieve points and the global point numbering (replaces creating remote bases) */
141770034214SMatthew G. Knepley   ierr = PetscSFConvertPartition(partSF, partSection, part, &renumbering, &pointSF);CHKERRQ(ierr);
141870034214SMatthew G. Knepley   if (flg) {
141970034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Point Partition:\n");CHKERRQ(ierr);
142070034214SMatthew G. Knepley     ierr = PetscSectionView(partSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
142170034214SMatthew G. Knepley     ierr = ISView(part, NULL);CHKERRQ(ierr);
142270034214SMatthew G. Knepley     ierr = PetscSFView(pointSF, NULL);CHKERRQ(ierr);
142370034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Point Renumbering after partition:\n");CHKERRQ(ierr);
142470034214SMatthew G. Knepley     ierr = ISLocalToGlobalMappingView(renumbering, NULL);CHKERRQ(ierr);
142570034214SMatthew G. Knepley   }
142677623264SMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
142770034214SMatthew G. Knepley 
1428a157612eSMichael Lange   /* Migrate data to a non-overlapping parallel DM */
1429cc71bff1SMichael Lange   ierr = DMPlexDistributeCones(dm, pointSF, renumbering, *dmParallel);CHKERRQ(ierr);
14300df0e737SMichael Lange   ierr = DMPlexDistributeCoordinates(dm, pointSF, *dmParallel);CHKERRQ(ierr);
14312eb1fa14SMichael Lange   ierr = DMPlexDistributeLabels(dm, pointSF, *dmParallel);CHKERRQ(ierr);
14329734c634SMichael Lange   ierr = DMPlexDistributeSetupHybrid(dm, pointSF, renumbering, *dmParallel);CHKERRQ(ierr);
1433a6f36705SMichael Lange   ierr = DMPlexDistributeSetupTree(dm, pointSF, renumbering, *dmParallel);CHKERRQ(ierr);
143470034214SMatthew G. Knepley 
1435a157612eSMichael Lange   /* Build the point SF without overlap */
1436a157612eSMichael Lange   ierr = DMPlexDistributeSF(dm, pointSF, partSection, part, NULL, NULL, *dmParallel);CHKERRQ(ierr);
143770034214SMatthew G. Knepley   if (flg) {
143815fff7beSMatthew G. Knepley     ierr = PetscPrintf(comm, "Point SF:\n");CHKERRQ(ierr);
143915fff7beSMatthew G. Knepley     ierr = PetscSFView((*dmParallel)->sf, NULL);CHKERRQ(ierr);
144070034214SMatthew G. Knepley   }
144170034214SMatthew G. Knepley 
1442a157612eSMichael Lange   if (overlap > 0) {
14433d822a50SMichael Lange     ierr = PetscLogEventBegin(DMPLEX_DistributeOverlap,dm,0,0,0);CHKERRQ(ierr);
1444a157612eSMichael Lange     /* Add the partition overlap to the distributed DM */
1445a157612eSMichael Lange     ierr = DMPlexDistributeOverlap(*dmParallel, overlap, renumbering, &overlapSF, &dmOverlap);CHKERRQ(ierr);
1446a157612eSMichael Lange     ierr = DMDestroy(dmParallel);CHKERRQ(ierr);
1447a157612eSMichael Lange     *dmParallel = dmOverlap;
1448c389ea9fSToby Isaac     if (flg) {
14493d822a50SMichael Lange       ierr = PetscPrintf(comm, "Overlap Migration SF:\n");CHKERRQ(ierr);
14503d822a50SMichael Lange       ierr = PetscSFView(overlapSF, NULL);CHKERRQ(ierr);
1451c389ea9fSToby Isaac     }
145243331d4aSMichael Lange 
145343331d4aSMichael Lange     /* Re-map the pointSF to establish the full migration pattern */
145443331d4aSMichael Lange     ierr = PetscSFGetGraph(pointSF, &nroots, NULL, NULL, &oldRemote);CHKERRQ(ierr);
145543331d4aSMichael Lange     ierr = PetscSFGetGraph(overlapSF, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr);
145643331d4aSMichael Lange     ierr = PetscMalloc1(nleaves, &newRemote);CHKERRQ(ierr);
145743331d4aSMichael Lange     ierr = PetscSFBcastBegin(overlapSF, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr);
145843331d4aSMichael Lange     ierr = PetscSFBcastEnd(overlapSF, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr);
145943331d4aSMichael Lange     ierr = PetscSFCreate(comm, &overlapPointSF);CHKERRQ(ierr);
146043331d4aSMichael Lange     ierr = PetscSFSetGraph(overlapPointSF, nroots, nleaves, NULL, PETSC_OWN_POINTER, newRemote, PETSC_OWN_POINTER);CHKERRQ(ierr);
146115fff7beSMatthew G. Knepley     ierr = PetscSFDestroy(&overlapSF);CHKERRQ(ierr);
146243331d4aSMichael Lange     ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr);
146343331d4aSMichael Lange     pointSF = overlapPointSF;
14643d822a50SMichael Lange     ierr = PetscLogEventEnd(DMPLEX_DistributeOverlap,dm,0,0,0);CHKERRQ(ierr);
1465c389ea9fSToby Isaac   }
1466bf5244c7SMatthew G. Knepley   /* Cleanup Partition */
1467bf5244c7SMatthew G. Knepley   ierr = ISLocalToGlobalMappingDestroy(&renumbering);CHKERRQ(ierr);
1468bf5244c7SMatthew G. Knepley   ierr = PetscSFDestroy(&partSF);CHKERRQ(ierr);
1469bf5244c7SMatthew G. Knepley   ierr = PetscSectionDestroy(&partSection);CHKERRQ(ierr);
1470bf5244c7SMatthew G. Knepley   ierr = ISDestroy(&part);CHKERRQ(ierr);
14714eca1733SMichael Lange   ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr);
14724eca1733SMichael Lange   ierr = ISDestroy(&cellPart);CHKERRQ(ierr);
147386719b60SMatthew G. Knepley   /* Copy BC */
147486719b60SMatthew G. Knepley   ierr = DMPlexCopyBoundary(dm, *dmParallel);CHKERRQ(ierr);
147570034214SMatthew G. Knepley   /* Cleanup */
147670034214SMatthew G. Knepley   if (sf) {*sf = pointSF;}
147770034214SMatthew G. Knepley   else    {ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr);}
147870034214SMatthew G. Knepley   ierr = DMSetFromOptions(*dmParallel);CHKERRQ(ierr);
147970034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
148070034214SMatthew G. Knepley   PetscFunctionReturn(0);
148170034214SMatthew G. Knepley }
1482a157612eSMichael Lange 
1483a157612eSMichael Lange #undef __FUNCT__
1484a157612eSMichael Lange #define __FUNCT__ "DMPlexDistributeOverlap"
1485a157612eSMichael Lange /*@C
1486a157612eSMichael Lange   DMPlexDistribute - Add partition overlap to a distributed non-overlapping DM.
1487a157612eSMichael Lange 
1488a157612eSMichael Lange   Not Collective
1489a157612eSMichael Lange 
1490a157612eSMichael Lange   Input Parameter:
1491a157612eSMichael Lange + dm  - The non-overlapping distrbuted DMPlex object
1492a157612eSMichael Lange - overlap - The overlap of partitions, 0 is the default
1493a157612eSMichael Lange 
1494a157612eSMichael Lange   Output Parameter:
1495a157612eSMichael Lange + sf - The PetscSF used for point distribution
1496a157612eSMichael Lange - dmOverlap - The overlapping distributed DMPlex object, or NULL
1497a157612eSMichael Lange 
1498a157612eSMichael Lange   Note: If the mesh was not distributed, the return value is NULL.
1499a157612eSMichael Lange 
1500a157612eSMichael Lange   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
1501a157612eSMichael Lange   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
1502a157612eSMichael Lange   representation on the mesh.
1503a157612eSMichael Lange 
1504a157612eSMichael Lange   Level: intermediate
1505a157612eSMichael Lange 
1506a157612eSMichael Lange .keywords: mesh, elements
1507a157612eSMichael Lange .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
1508a157612eSMichael Lange @*/
1509a157612eSMichael Lange PetscErrorCode DMPlexDistributeOverlap(DM dm, PetscInt overlap, ISLocalToGlobalMapping renumbering, PetscSF *sf, DM *dmOverlap)
1510a157612eSMichael Lange {
1511a157612eSMichael Lange   MPI_Comm               comm;
1512a157612eSMichael Lange   PetscMPIInt            rank;
1513a157612eSMichael Lange   PetscSection           rootSection, leafSection;
1514a157612eSMichael Lange   IS                     rootrank, leafrank;
1515a157612eSMichael Lange   PetscSection           coneSection;
1516a157612eSMichael Lange   PetscSF                overlapSF, migrationSF, pointSF, newPointSF;
1517a157612eSMichael Lange   PetscSFNode           *ghostRemote;
1518a157612eSMichael Lange   const PetscSFNode     *overlapRemote;
1519a157612eSMichael Lange   ISLocalToGlobalMapping overlapRenumbering;
1520a157612eSMichael Lange   const PetscInt        *renumberingArray, *overlapLocal;
1521a157612eSMichael Lange   PetscInt               dim, p, pStart, pEnd, conesSize, idx;
1522a157612eSMichael Lange   PetscInt               numGhostPoints, numOverlapPoints, numSharedPoints, overlapLeaves;
1523a157612eSMichael Lange   PetscInt              *cones, *ghostLocal, *overlapRenumberingArray, *pointIDs, *recvPointIDs;
1524a157612eSMichael Lange   PetscErrorCode         ierr;
1525a157612eSMichael Lange 
1526a157612eSMichael Lange   PetscFunctionBegin;
1527a157612eSMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1528a157612eSMichael Lange   if (sf) PetscValidPointer(sf, 3);
1529a157612eSMichael Lange   PetscValidPointer(dmOverlap, 4);
1530a157612eSMichael Lange 
1531a157612eSMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1532a157612eSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1533a157612eSMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1534a157612eSMichael Lange 
1535a157612eSMichael Lange   /* Compute point overlap with neighbouring processes on the distributed DM */
1536a157612eSMichael Lange   ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1537a157612eSMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
1538a157612eSMichael Lange   ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr);
1539a157612eSMichael Lange   ierr = DMPlexDistributeOwnership(dm, rootSection, &rootrank, leafSection, &leafrank);CHKERRQ(ierr);
1540a157612eSMichael Lange   ierr = DMPlexCreateOverlap(dm, rootSection, rootrank, leafSection, leafrank, &overlapSF);CHKERRQ(ierr);
154115fff7beSMatthew G. Knepley   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
154215fff7beSMatthew G. Knepley   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
154315fff7beSMatthew G. Knepley   ierr = ISDestroy(&rootrank);CHKERRQ(ierr);
154415fff7beSMatthew G. Knepley   ierr = ISDestroy(&leafrank);CHKERRQ(ierr);
1545a157612eSMichael Lange   ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1546a157612eSMichael Lange 
1547a157612eSMichael Lange   /* Build dense migration SF that maps the non-overlapping partition to the overlapping one */
1548a157612eSMichael Lange   ierr = DMPlexCreateOverlapMigrationSF(dm, overlapSF, &migrationSF);CHKERRQ(ierr);
1549a157612eSMichael Lange 
1550a157612eSMichael Lange   /* Convert cones to global numbering before migrating them */
1551a157612eSMichael Lange   ierr = DMPlexGetConeSection(dm, &coneSection);CHKERRQ(ierr);
1552a157612eSMichael Lange   ierr = PetscSectionGetStorageSize(coneSection, &conesSize);CHKERRQ(ierr);
1553a157612eSMichael Lange   ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
1554a157612eSMichael Lange   ierr = ISLocalToGlobalMappingApplyBlock(renumbering, conesSize, cones, cones);CHKERRQ(ierr);
1555a157612eSMichael Lange 
1556a157612eSMichael Lange   /* Derive the new local-to-global mapping from the old one */
1557a157612eSMichael Lange   ierr = PetscSFGetGraph(migrationSF, NULL, &overlapLeaves, &overlapLocal, &overlapRemote);CHKERRQ(ierr);
1558a157612eSMichael Lange   ierr = PetscMalloc1(overlapLeaves, &overlapRenumberingArray);CHKERRQ(ierr);
1559a157612eSMichael Lange   ierr = ISLocalToGlobalMappingGetBlockIndices(renumbering, &renumberingArray);CHKERRQ(ierr);
1560729b3788SMatthew G. Knepley   ierr = PetscSFBcastBegin(migrationSF, MPIU_INT, (PetscInt *) renumberingArray, overlapRenumberingArray);CHKERRQ(ierr);
1561729b3788SMatthew G. Knepley   ierr = PetscSFBcastEnd(migrationSF, MPIU_INT, (PetscInt *) renumberingArray, overlapRenumberingArray);CHKERRQ(ierr);
1562a157612eSMichael Lange   ierr = ISLocalToGlobalMappingCreate(comm, 1, overlapLeaves, (const PetscInt*) overlapRenumberingArray, PETSC_OWN_POINTER, &overlapRenumbering);CHKERRQ(ierr);
1563a157612eSMichael Lange 
1564a157612eSMichael Lange   /* Build the overlapping DM */
1565a157612eSMichael Lange   ierr = DMPlexCreate(comm, dmOverlap);CHKERRQ(ierr);
1566a157612eSMichael Lange   ierr = DMSetDimension(*dmOverlap, dim);CHKERRQ(ierr);
1567a157612eSMichael Lange   ierr = PetscObjectSetName((PetscObject) *dmOverlap, "Parallel Mesh");CHKERRQ(ierr);
1568a157612eSMichael Lange   ierr = DMPlexDistributeCones(dm, migrationSF, overlapRenumbering, *dmOverlap);CHKERRQ(ierr);
1569a157612eSMichael Lange   ierr = DMPlexDistributeCoordinates(dm, migrationSF, *dmOverlap);CHKERRQ(ierr);
1570a157612eSMichael Lange   ierr = DMPlexDistributeLabels(dm, migrationSF, *dmOverlap);CHKERRQ(ierr);
1571a157612eSMichael Lange   ierr = DMPlexDistributeSetupHybrid(dm, migrationSF, overlapRenumbering, *dmOverlap);CHKERRQ(ierr);
1572a157612eSMichael Lange 
1573a157612eSMichael Lange   /* Build the new point SF by propagating the depthShift generate remote root indices */
1574a157612eSMichael Lange   ierr = DMGetPointSF(dm, &pointSF);CHKERRQ(ierr);
1575a157612eSMichael Lange   ierr = PetscSFGetGraph(pointSF, NULL, &numSharedPoints, NULL, NULL);CHKERRQ(ierr);
1576a157612eSMichael Lange   ierr = PetscSFGetGraph(overlapSF, NULL, &numOverlapPoints, NULL, NULL);CHKERRQ(ierr);
1577a157612eSMichael Lange   numGhostPoints = numSharedPoints + numOverlapPoints;
157809b7985cSMatthew G. Knepley   ierr = PetscMalloc1(numGhostPoints, &ghostLocal);CHKERRQ(ierr);
157909b7985cSMatthew G. Knepley   ierr = PetscMalloc1(numGhostPoints, &ghostRemote);CHKERRQ(ierr);
1580a157612eSMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
1581a157612eSMichael Lange   ierr = PetscMalloc2(pEnd-pStart, &pointIDs, overlapLeaves, &recvPointIDs);CHKERRQ(ierr);
1582a157612eSMichael Lange   for (p=0; p<overlapLeaves; p++) {
1583a157612eSMichael Lange     if (overlapRemote[p].rank == rank) pointIDs[overlapRemote[p].index] = overlapLocal[p];
1584a157612eSMichael Lange   }
1585a157612eSMichael Lange   ierr = PetscSFBcastBegin(migrationSF, MPIU_INT, pointIDs, recvPointIDs);CHKERRQ(ierr);
1586a157612eSMichael Lange   ierr = PetscSFBcastEnd(migrationSF, MPIU_INT, pointIDs, recvPointIDs);CHKERRQ(ierr);
1587a157612eSMichael Lange   for (idx=0, p=0; p<overlapLeaves; p++) {
1588a157612eSMichael Lange     if (overlapRemote[p].rank != rank) {
1589a157612eSMichael Lange       ghostLocal[idx] = overlapLocal[p];
1590a157612eSMichael Lange       ghostRemote[idx].index = recvPointIDs[p];
1591a157612eSMichael Lange       ghostRemote[idx].rank = overlapRemote[p].rank;
1592a157612eSMichael Lange       idx++;
1593a157612eSMichael Lange     }
1594a157612eSMichael Lange   }
1595a157612eSMichael Lange   ierr = DMPlexGetChart(*dmOverlap, &pStart, &pEnd);CHKERRQ(ierr);
1596a157612eSMichael Lange   ierr = PetscSFCreate(comm, &newPointSF);;CHKERRQ(ierr);
1597a157612eSMichael Lange   ierr = PetscSFSetGraph(newPointSF, pEnd - pStart, numGhostPoints, ghostLocal, PETSC_OWN_POINTER, ghostRemote, PETSC_OWN_POINTER);CHKERRQ(ierr);
1598a157612eSMichael Lange   ierr = DMSetPointSF(*dmOverlap, newPointSF);CHKERRQ(ierr);
159915fff7beSMatthew G. Knepley   ierr = PetscSFDestroy(&newPointSF);CHKERRQ(ierr);
1600a157612eSMichael Lange   /* Cleanup overlap partition */
1601a157612eSMichael Lange   ierr = ISLocalToGlobalMappingDestroy(&overlapRenumbering);CHKERRQ(ierr);
1602a157612eSMichael Lange   ierr = PetscSFDestroy(&overlapSF);CHKERRQ(ierr);
1603a157612eSMichael Lange   ierr = PetscFree2(pointIDs, recvPointIDs);CHKERRQ(ierr);
1604a157612eSMichael Lange   if (sf) *sf = migrationSF;
1605a157612eSMichael Lange   else    {ierr = PetscSFDestroy(&migrationSF);CHKERRQ(ierr);}
1606a157612eSMichael Lange   PetscFunctionReturn(0);
1607a157612eSMichael Lange }
1608