xref: /petsc/src/dm/impls/plex/plexdistribute.c (revision 3533c52b29edbc7b9f52f0c179567251e6920429)
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;
564e540f424SMichael Lange   PetscBool          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);
618e540f424SMichael Lange   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-overlap_view", &flg);CHKERRQ(ierr);
619e540f424SMichael Lange   if (flg) {
620b0a623aaSMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr);
621b0a623aaSMatthew G. Knepley     ierr = DMLabelView(ovAdjByRank, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
622b0a623aaSMatthew G. Knepley   }
623b0a623aaSMatthew G. Knepley   /* Convert to (point, rank) and use actual owners */
624e540f424SMichael Lange   ierr = PetscSectionCreate(comm, &ovRootSection);CHKERRQ(ierr);
625e540f424SMichael Lange   ierr = PetscObjectSetName((PetscObject) ovRootSection, "Overlap Root Section");CHKERRQ(ierr);
626b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetChart(ovRootSection, 0, numProcs);CHKERRQ(ierr);
627b0a623aaSMatthew G. Knepley   ierr = DMLabelGetValueIS(ovAdjByRank, &valueIS);CHKERRQ(ierr);
628b0a623aaSMatthew G. Knepley   ierr = ISGetLocalSize(valueIS, &numNeighbors);CHKERRQ(ierr);
629b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(valueIS, &neighbors);CHKERRQ(ierr);
630b0a623aaSMatthew G. Knepley   for (n = 0; n < numNeighbors; ++n) {
631b0a623aaSMatthew G. Knepley     PetscInt numPoints;
632b0a623aaSMatthew G. Knepley 
633b0a623aaSMatthew G. Knepley     ierr = DMLabelGetStratumSize(ovAdjByRank, neighbors[n], &numPoints);CHKERRQ(ierr);
634b0a623aaSMatthew G. Knepley     ierr = PetscSectionAddDof(ovRootSection, neighbors[n], numPoints);CHKERRQ(ierr);
635b0a623aaSMatthew G. Knepley   }
636b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetUp(ovRootSection);CHKERRQ(ierr);
637b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(ovRootSection, &ovSize);CHKERRQ(ierr);
638e540f424SMichael Lange   ierr = PetscMalloc1(ovSize, &ovRootPoints);CHKERRQ(ierr);
639b0a623aaSMatthew G. Knepley   ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr);
640b0a623aaSMatthew G. Knepley   for (n = 0; n < numNeighbors; ++n) {
641b0a623aaSMatthew G. Knepley     IS              pointIS;
642b0a623aaSMatthew G. Knepley     const PetscInt *points;
643b0a623aaSMatthew G. Knepley     PetscInt        off, numPoints, p;
644b0a623aaSMatthew G. Knepley 
645b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(ovRootSection, neighbors[n], &off);CHKERRQ(ierr);
646b0a623aaSMatthew G. Knepley     ierr = DMLabelGetStratumIS(ovAdjByRank, neighbors[n], &pointIS);CHKERRQ(ierr);
647b0a623aaSMatthew G. Knepley     ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
648b0a623aaSMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
649b0a623aaSMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
650b0a623aaSMatthew G. Knepley       ierr = PetscFindInt(points[p], nleaves, local, &l);CHKERRQ(ierr);
651e540f424SMichael Lange       if (l >= 0) {ovRootPoints[off+p] = remote[l];}
652e540f424SMichael Lange       else        {ovRootPoints[off+p].index = points[p]; ovRootPoints[off+p].rank = rank;}
653b0a623aaSMatthew G. Knepley     }
654b0a623aaSMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
655b0a623aaSMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
656b0a623aaSMatthew G. Knepley   }
657b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(valueIS, &neighbors);CHKERRQ(ierr);
658b0a623aaSMatthew G. Knepley   ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
659b0a623aaSMatthew G. Knepley   ierr = DMLabelDestroy(&ovAdjByRank);CHKERRQ(ierr);
660b0a623aaSMatthew G. Knepley   /* Make process SF */
661b0a623aaSMatthew G. Knepley   ierr = DMPlexCreateTwoSidedProcessSF(dm, sfPoint, rootSection, rootrank, leafSection, leafrank, NULL, &sfProc);CHKERRQ(ierr);
662e540f424SMichael Lange   if (flg) {
663b0a623aaSMatthew G. Knepley     ierr = PetscSFView(sfProc, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
664b0a623aaSMatthew G. Knepley   }
665b0a623aaSMatthew G. Knepley   /* Communicate overlap */
666e540f424SMichael Lange   ierr = PetscSectionCreate(comm, &ovLeafSection);CHKERRQ(ierr);
667e540f424SMichael Lange   ierr = PetscObjectSetName((PetscObject) ovLeafSection, "Overlap Leaf Section");CHKERRQ(ierr);
668e540f424SMichael Lange   ierr = DMPlexDistributeData(dm, sfProc, ovRootSection, MPIU_2INT, ovRootPoints, ovLeafSection, (void**) &ovLeafPoints);CHKERRQ(ierr);
669e540f424SMichael Lange   /* Filter remote contributions (ovLeafPoints) into the overlapSF */
670e540f424SMichael Lange   ierr = DMLabelCreate("ovLeafs", &ovLeafLabel);CHKERRQ(ierr);
671e540f424SMichael Lange   ierr = PetscSectionGetStorageSize(ovLeafSection, &ovSize);CHKERRQ(ierr);
672e540f424SMichael Lange   for (p = 0; p < ovSize; p++) {
673e540f424SMichael Lange     /* Don't import points from yourself */
674e540f424SMichael Lange     if (ovLeafPoints[p].rank == rank) continue;
675e540f424SMichael Lange     ierr = DMLabelSetValue(ovLeafLabel, ovLeafPoints[p].index, ovLeafPoints[p].rank);CHKERRQ(ierr);
676e540f424SMichael Lange   }
677e540f424SMichael Lange   /* Don't import points already in the pointSF */
678e540f424SMichael Lange   for (l = 0; l < nleaves; ++l) {
679e540f424SMichael Lange     ierr = DMLabelClearValue(ovLeafLabel, remote[l].index, remote[l].rank);CHKERRQ(ierr);
680e540f424SMichael Lange   }
681e540f424SMichael Lange   for (numRemote = 0, n = 0; n < numProcs; ++n) {
682e540f424SMichael Lange     PetscInt numPoints;
683e540f424SMichael Lange     ierr = DMLabelGetStratumSize(ovLeafLabel, n, &numPoints);CHKERRQ(ierr);
684e540f424SMichael Lange     numRemote += numPoints;
685e540f424SMichael Lange   }
686e540f424SMichael Lange   ierr = PetscMalloc1(numRemote, &remotePoints);CHKERRQ(ierr);
687e540f424SMichael Lange   for (idx = 0, n = 0; n < numProcs; ++n) {
688e540f424SMichael Lange     IS remoteRootIS;
689e540f424SMichael Lange     PetscInt numPoints;
690e540f424SMichael Lange     const PetscInt *remoteRoots;
691e540f424SMichael Lange     ierr = DMLabelGetStratumSize(ovLeafLabel, n, &numPoints);CHKERRQ(ierr);
692e540f424SMichael Lange     if (numPoints <= 0) continue;
693e540f424SMichael Lange     ierr = DMLabelGetStratumIS(ovLeafLabel, n, &remoteRootIS);CHKERRQ(ierr);
694e540f424SMichael Lange     ierr = ISGetIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
695e540f424SMichael Lange     for (p = 0; p < numPoints; p++) {
696e540f424SMichael Lange       remotePoints[idx].index = remoteRoots[p];
697e540f424SMichael Lange       remotePoints[idx].rank = n;
698e540f424SMichael Lange       idx++;
699e540f424SMichael Lange     }
700e540f424SMichael Lange     ierr = ISRestoreIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
70115fff7beSMatthew G. Knepley     ierr = ISDestroy(&remoteRootIS);CHKERRQ(ierr);
702e540f424SMichael Lange   }
70315fff7beSMatthew G. Knepley   ierr = DMLabelDestroy(&ovLeafLabel);CHKERRQ(ierr);
704e540f424SMichael Lange   ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), overlapSF);CHKERRQ(ierr);
705e540f424SMichael Lange   ierr = PetscObjectSetName((PetscObject) *overlapSF, "Overlap SF");CHKERRQ(ierr);
706e540f424SMichael Lange   ierr = PetscSFSetFromOptions(*overlapSF);CHKERRQ(ierr);
707e540f424SMichael Lange   ierr = PetscSFSetGraph(*overlapSF, pEnd-pStart, numRemote, NULL, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
708e540f424SMichael Lange   if (flg) {
709e540f424SMichael Lange     ierr = PetscPrintf(comm, "Overlap SF\n");CHKERRQ(ierr);
710e540f424SMichael Lange     ierr = PetscSFView(*overlapSF, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
711e540f424SMichael Lange   }
712e540f424SMichael Lange   /* Clean up */
713b0a623aaSMatthew G. Knepley   ierr = PetscSFDestroy(&sfProc);CHKERRQ(ierr);
714e540f424SMichael Lange   ierr = PetscFree(ovRootPoints);CHKERRQ(ierr);
715e540f424SMichael Lange   ierr = PetscSectionDestroy(&ovRootSection);CHKERRQ(ierr);
716e540f424SMichael Lange   ierr = PetscFree(ovLeafPoints);CHKERRQ(ierr);
717e540f424SMichael Lange   ierr = PetscSectionDestroy(&ovLeafSection);CHKERRQ(ierr);
718b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
719b0a623aaSMatthew G. Knepley }
72070034214SMatthew G. Knepley 
72170034214SMatthew G. Knepley #undef __FUNCT__
72246f9b1c3SMichael Lange #define __FUNCT__ "DMPlexCreateOverlapMigrationSF"
72346f9b1c3SMichael Lange PetscErrorCode DMPlexCreateOverlapMigrationSF(DM dm, PetscSF overlapSF, PetscSF *migrationSF)
72446f9b1c3SMichael Lange {
72546f9b1c3SMichael Lange   MPI_Comm           comm;
72646f9b1c3SMichael Lange   PetscMPIInt        rank, numProcs;
72746f9b1c3SMichael Lange   PetscInt           d, dim, p, pStart, pEnd, nroots, nleaves, newLeaves, point, numSharedPoints;
72846f9b1c3SMichael Lange   PetscInt          *pointDepths, *remoteDepths, *ilocal;
72946f9b1c3SMichael Lange   PetscInt          *depthRecv, *depthShift, *depthIdx;
73046f9b1c3SMichael Lange   PetscSFNode       *iremote;
73146f9b1c3SMichael Lange   PetscSF            pointSF;
73246f9b1c3SMichael Lange   const PetscInt    *sharedLocal;
73346f9b1c3SMichael Lange   const PetscSFNode *overlapRemote, *sharedRemote;
73446f9b1c3SMichael Lange   PetscErrorCode     ierr;
73546f9b1c3SMichael Lange 
73646f9b1c3SMichael Lange   PetscFunctionBegin;
73746f9b1c3SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
73846f9b1c3SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
73946f9b1c3SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
74046f9b1c3SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
74146f9b1c3SMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
74246f9b1c3SMichael Lange 
74346f9b1c3SMichael Lange   /* Before building the migration SF we need to know the new stratum offsets */
74446f9b1c3SMichael Lange   ierr = PetscSFGetGraph(overlapSF, &nroots, &nleaves, NULL, &overlapRemote);CHKERRQ(ierr);
74546f9b1c3SMichael Lange   ierr = PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);CHKERRQ(ierr);
74646f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
74746f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
74846f9b1c3SMichael Lange     for (p=pStart; p<pEnd; p++) pointDepths[p] = d;
74946f9b1c3SMichael Lange   }
75046f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) remoteDepths[p] = -1;
75146f9b1c3SMichael Lange   ierr = PetscSFBcastBegin(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
75246f9b1c3SMichael Lange   ierr = PetscSFBcastEnd(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
75346f9b1c3SMichael Lange 
75446f9b1c3SMichael Lange   /* Count recevied points in each stratum and compute the internal strata shift */
75546f9b1c3SMichael Lange   ierr = PetscMalloc3(dim+1, &depthRecv, dim+1, &depthShift, dim+1, &depthIdx);CHKERRQ(ierr);
75646f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) depthRecv[d]=0;
75746f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) depthRecv[remoteDepths[p]]++;
75846f9b1c3SMichael Lange   depthShift[dim] = 0;
75946f9b1c3SMichael Lange   for (d=0; d<dim; d++) depthShift[d] = depthRecv[dim];
76046f9b1c3SMichael Lange   for (d=1; d<dim; d++) depthShift[d] += depthRecv[0];
76146f9b1c3SMichael Lange   for (d=dim-2; d>0; d--) depthShift[d] += depthRecv[d+1];
76246f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
76346f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
76446f9b1c3SMichael Lange     depthIdx[d] = pStart + depthShift[d];
76546f9b1c3SMichael Lange   }
76646f9b1c3SMichael Lange 
76746f9b1c3SMichael Lange   /* Form the overlap SF build an SF that describes the full overlap migration SF */
76846f9b1c3SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
76946f9b1c3SMichael Lange   newLeaves = pEnd - pStart + nleaves;
77009b7985cSMatthew G. Knepley   ierr = PetscMalloc1(newLeaves, &ilocal);CHKERRQ(ierr);
77109b7985cSMatthew G. Knepley   ierr = PetscMalloc1(newLeaves, &iremote);CHKERRQ(ierr);
77246f9b1c3SMichael Lange   /* First map local points to themselves */
77346f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
77446f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
77546f9b1c3SMichael Lange     for (p=pStart; p<pEnd; p++) {
77646f9b1c3SMichael Lange       point = p + depthShift[d];
77746f9b1c3SMichael Lange       ilocal[point] = point;
77846f9b1c3SMichael Lange       iremote[point].index = p;
77946f9b1c3SMichael Lange       iremote[point].rank = rank;
78046f9b1c3SMichael Lange       depthIdx[d]++;
78146f9b1c3SMichael Lange     }
78246f9b1c3SMichael Lange   }
78346f9b1c3SMichael Lange 
78446f9b1c3SMichael Lange   /* Add in the remote roots for currently shared points */
78546f9b1c3SMichael Lange   ierr = DMGetPointSF(dm, &pointSF);CHKERRQ(ierr);
78646f9b1c3SMichael Lange   ierr = PetscSFGetGraph(pointSF, NULL, &numSharedPoints, &sharedLocal, &sharedRemote);CHKERRQ(ierr);
78746f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
78846f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
78946f9b1c3SMichael Lange     for (p=0; p<numSharedPoints; p++) {
79046f9b1c3SMichael Lange       if (pStart <= sharedLocal[p] && sharedLocal[p] < pEnd) {
79146f9b1c3SMichael Lange         point = sharedLocal[p] + depthShift[d];
79246f9b1c3SMichael Lange         iremote[point].index = sharedRemote[p].index;
79346f9b1c3SMichael Lange         iremote[point].rank = sharedRemote[p].rank;
79446f9b1c3SMichael Lange       }
79546f9b1c3SMichael Lange     }
79646f9b1c3SMichael Lange   }
79746f9b1c3SMichael Lange 
79846f9b1c3SMichael Lange   /* Now add the incoming overlap points */
79946f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) {
80046f9b1c3SMichael Lange     point = depthIdx[remoteDepths[p]];
80146f9b1c3SMichael Lange     ilocal[point] = point;
80246f9b1c3SMichael Lange     iremote[point].index = overlapRemote[p].index;
80346f9b1c3SMichael Lange     iremote[point].rank = overlapRemote[p].rank;
80446f9b1c3SMichael Lange     depthIdx[remoteDepths[p]]++;
80546f9b1c3SMichael Lange   }
80615fff7beSMatthew G. Knepley   ierr = PetscFree2(pointDepths,remoteDepths);CHKERRQ(ierr);
80746f9b1c3SMichael Lange 
80846f9b1c3SMichael Lange   ierr = PetscSFCreate(comm, migrationSF);CHKERRQ(ierr);
80946f9b1c3SMichael Lange   ierr = PetscObjectSetName((PetscObject) *migrationSF, "Overlap Migration SF");CHKERRQ(ierr);
81046f9b1c3SMichael Lange   ierr = PetscSFSetFromOptions(*migrationSF);CHKERRQ(ierr);
81146f9b1c3SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
81246f9b1c3SMichael Lange   ierr = PetscSFSetGraph(*migrationSF, pEnd-pStart, newLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
81346f9b1c3SMichael Lange 
81446f9b1c3SMichael Lange   ierr = PetscFree3(depthRecv, depthShift, depthIdx);CHKERRQ(ierr);
81546f9b1c3SMichael Lange   PetscFunctionReturn(0);
81646f9b1c3SMichael Lange }
81746f9b1c3SMichael Lange 
81846f9b1c3SMichael Lange #undef __FUNCT__
81970034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeField"
82070034214SMatthew G. Knepley /*@
82170034214SMatthew G. Knepley   DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
82270034214SMatthew G. Knepley 
82370034214SMatthew G. Knepley   Collective on DM
82470034214SMatthew G. Knepley 
82570034214SMatthew G. Knepley   Input Parameters:
82670034214SMatthew G. Knepley + dm - The DMPlex object
82770034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
82870034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
82970034214SMatthew G. Knepley - originalVec - The existing data
83070034214SMatthew G. Knepley 
83170034214SMatthew G. Knepley   Output Parameters:
83270034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout
83370034214SMatthew G. Knepley - newVec - The new data
83470034214SMatthew G. Knepley 
83570034214SMatthew G. Knepley   Level: developer
83670034214SMatthew G. Knepley 
8371e8fc25dSMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeFieldIS(), DMPlexDistributeData()
83870034214SMatthew G. Knepley @*/
83970034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec)
84070034214SMatthew G. Knepley {
84170034214SMatthew G. Knepley   PetscSF        fieldSF;
84270034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
84370034214SMatthew G. Knepley   PetscScalar   *originalValues, *newValues;
84470034214SMatthew G. Knepley   PetscErrorCode ierr;
84570034214SMatthew G. Knepley 
84670034214SMatthew G. Knepley   PetscFunctionBegin;
84770034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
84870034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
84970034214SMatthew G. Knepley 
85070034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
85170034214SMatthew G. Knepley   ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr);
85270034214SMatthew G. Knepley   ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr);
85370034214SMatthew G. Knepley 
85470034214SMatthew G. Knepley   ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr);
85570034214SMatthew G. Knepley   ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr);
85670034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
85770034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
85870034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
85970034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
86070034214SMatthew G. Knepley   ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr);
86170034214SMatthew G. Knepley   ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr);
86270034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
86370034214SMatthew G. Knepley   PetscFunctionReturn(0);
86470034214SMatthew G. Knepley }
86570034214SMatthew G. Knepley 
86670034214SMatthew G. Knepley #undef __FUNCT__
8671e8fc25dSMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeFieldIS"
8681e8fc25dSMatthew G. Knepley /*@
8691e8fc25dSMatthew G. Knepley   DMPlexDistributeFieldIS - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
8701e8fc25dSMatthew G. Knepley 
8711e8fc25dSMatthew G. Knepley   Collective on DM
8721e8fc25dSMatthew G. Knepley 
8731e8fc25dSMatthew G. Knepley   Input Parameters:
8741e8fc25dSMatthew G. Knepley + dm - The DMPlex object
8751e8fc25dSMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
8761e8fc25dSMatthew G. Knepley . originalSection - The PetscSection for existing data layout
8771e8fc25dSMatthew G. Knepley - originalIS - The existing data
8781e8fc25dSMatthew G. Knepley 
8791e8fc25dSMatthew G. Knepley   Output Parameters:
8801e8fc25dSMatthew G. Knepley + newSection - The PetscSF describing the new data layout
8811e8fc25dSMatthew G. Knepley - newIS - The new data
8821e8fc25dSMatthew G. Knepley 
8831e8fc25dSMatthew G. Knepley   Level: developer
8841e8fc25dSMatthew G. Knepley 
8851e8fc25dSMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexDistributeData()
8861e8fc25dSMatthew G. Knepley @*/
8871e8fc25dSMatthew G. Knepley PetscErrorCode DMPlexDistributeFieldIS(DM dm, PetscSF pointSF, PetscSection originalSection, IS originalIS, PetscSection newSection, IS *newIS)
8881e8fc25dSMatthew G. Knepley {
8891e8fc25dSMatthew G. Knepley   PetscSF         fieldSF;
8901e8fc25dSMatthew G. Knepley   PetscInt       *newValues, *remoteOffsets, fieldSize;
8911e8fc25dSMatthew G. Knepley   const PetscInt *originalValues;
8921e8fc25dSMatthew G. Knepley   PetscErrorCode  ierr;
8931e8fc25dSMatthew G. Knepley 
8941e8fc25dSMatthew G. Knepley   PetscFunctionBegin;
8951e8fc25dSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
8961e8fc25dSMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
8971e8fc25dSMatthew G. Knepley 
8981e8fc25dSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
8991e8fc25dSMatthew G. Knepley   ierr = PetscMalloc(fieldSize * sizeof(PetscInt), &newValues);CHKERRQ(ierr);
9001e8fc25dSMatthew G. Knepley 
9011e8fc25dSMatthew G. Knepley   ierr = ISGetIndices(originalIS, &originalValues);CHKERRQ(ierr);
9021e8fc25dSMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
9031e8fc25dSMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, MPIU_INT, originalValues, newValues);CHKERRQ(ierr);
9041e8fc25dSMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, MPIU_INT, originalValues, newValues);CHKERRQ(ierr);
9051e8fc25dSMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
9061e8fc25dSMatthew G. Knepley   ierr = ISRestoreIndices(originalIS, &originalValues);CHKERRQ(ierr);
9071e8fc25dSMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) pointSF), fieldSize, newValues, PETSC_OWN_POINTER, newIS);CHKERRQ(ierr);
9081e8fc25dSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
9091e8fc25dSMatthew G. Knepley   PetscFunctionReturn(0);
9101e8fc25dSMatthew G. Knepley }
9111e8fc25dSMatthew G. Knepley 
9121e8fc25dSMatthew G. Knepley #undef __FUNCT__
91370034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeData"
91470034214SMatthew G. Knepley /*@
91570034214SMatthew G. Knepley   DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
91670034214SMatthew G. Knepley 
91770034214SMatthew G. Knepley   Collective on DM
91870034214SMatthew G. Knepley 
91970034214SMatthew G. Knepley   Input Parameters:
92070034214SMatthew G. Knepley + dm - The DMPlex object
92170034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
92270034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
92370034214SMatthew G. Knepley . datatype - The type of data
92470034214SMatthew G. Knepley - originalData - The existing data
92570034214SMatthew G. Knepley 
92670034214SMatthew G. Knepley   Output Parameters:
9271e8fc25dSMatthew G. Knepley + newSection - The PetscSection describing the new data layout
92870034214SMatthew G. Knepley - newData - The new data
92970034214SMatthew G. Knepley 
93070034214SMatthew G. Knepley   Level: developer
93170034214SMatthew G. Knepley 
93270034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField()
93370034214SMatthew G. Knepley @*/
93470034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData)
93570034214SMatthew G. Knepley {
93670034214SMatthew G. Knepley   PetscSF        fieldSF;
93770034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
93870034214SMatthew G. Knepley   PetscMPIInt    dataSize;
93970034214SMatthew G. Knepley   PetscErrorCode ierr;
94070034214SMatthew G. Knepley 
94170034214SMatthew G. Knepley   PetscFunctionBegin;
94270034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
94370034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
94470034214SMatthew G. Knepley 
94570034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
94670034214SMatthew G. Knepley   ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr);
94770034214SMatthew G. Knepley   ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr);
94870034214SMatthew G. Knepley 
94970034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
95070034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
95170034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
95270034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
95370034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
95470034214SMatthew G. Knepley   PetscFunctionReturn(0);
95570034214SMatthew G. Knepley }
95670034214SMatthew G. Knepley 
95770034214SMatthew G. Knepley #undef __FUNCT__
958cc71bff1SMichael Lange #define __FUNCT__ "DMPlexDistributeCones"
959cc71bff1SMichael Lange PetscErrorCode DMPlexDistributeCones(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
960cc71bff1SMichael Lange {
961cc71bff1SMichael Lange   DM_Plex               *pmesh   = (DM_Plex*) (dmParallel)->data;
962cc71bff1SMichael Lange   MPI_Comm               comm;
963cc71bff1SMichael Lange   PetscSF                coneSF;
964cc71bff1SMichael Lange   PetscSection           originalConeSection, newConeSection;
965cc71bff1SMichael Lange   PetscInt              *remoteOffsets, *cones, *newCones, newConesSize;
966cc71bff1SMichael Lange   PetscBool              flg;
967cc71bff1SMichael Lange   PetscErrorCode         ierr;
968cc71bff1SMichael Lange 
969cc71bff1SMichael Lange   PetscFunctionBegin;
970cc71bff1SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
971cc71bff1SMichael Lange   PetscValidPointer(dmParallel,4);
972cc71bff1SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
973cc71bff1SMichael Lange 
974cc71bff1SMichael Lange   /* Distribute cone section */
975cc71bff1SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
976cc71bff1SMichael Lange   ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr);
977cc71bff1SMichael Lange   ierr = DMPlexGetConeSection(dmParallel, &newConeSection);CHKERRQ(ierr);
978cc71bff1SMichael Lange   ierr = PetscSFDistributeSection(migrationSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr);
979cc71bff1SMichael Lange   ierr = DMSetUp(dmParallel);CHKERRQ(ierr);
980cc71bff1SMichael Lange   {
981cc71bff1SMichael Lange     PetscInt pStart, pEnd, p;
982cc71bff1SMichael Lange 
983cc71bff1SMichael Lange     ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr);
984cc71bff1SMichael Lange     for (p = pStart; p < pEnd; ++p) {
985cc71bff1SMichael Lange       PetscInt coneSize;
986cc71bff1SMichael Lange       ierr               = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr);
987cc71bff1SMichael Lange       pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize);
988cc71bff1SMichael Lange     }
989cc71bff1SMichael Lange   }
990cc71bff1SMichael Lange   /* Communicate and renumber cones */
991cc71bff1SMichael Lange   ierr = PetscSFCreateSectionSF(migrationSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr);
992cc71bff1SMichael Lange   ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
993cc71bff1SMichael Lange   ierr = DMPlexGetCones(dmParallel, &newCones);CHKERRQ(ierr);
994cc71bff1SMichael Lange   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
995cc71bff1SMichael Lange   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
996cc71bff1SMichael Lange   ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr);
997cc71bff1SMichael Lange   ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr);
998*3533c52bSMatthew G. Knepley #if PETSC_USE_DEBUG
999*3533c52bSMatthew G. Knepley   {
1000*3533c52bSMatthew G. Knepley     PetscInt  p;
1001*3533c52bSMatthew G. Knepley     PetscBool valid = PETSC_TRUE;
1002*3533c52bSMatthew G. Knepley     for (p = 0; p < newConesSize; ++p) {
1003*3533c52bSMatthew G. Knepley       if (newCones[p] < 0) {valid = PETSC_FALSE; ierr = PetscPrintf(PETSC_COMM_SELF, "Point %d not in overlap SF\n", p);CHKERRQ(ierr);}
1004*3533c52bSMatthew G. Knepley     }
1005*3533c52bSMatthew G. Knepley     if (!valid) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map");
1006*3533c52bSMatthew G. Knepley   }
1007*3533c52bSMatthew G. Knepley #endif
1008cc71bff1SMichael Lange   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr);
1009cc71bff1SMichael Lange   if (flg) {
1010cc71bff1SMichael Lange     ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr);
1011cc71bff1SMichael Lange     ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1012cc71bff1SMichael Lange     ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr);
1013cc71bff1SMichael Lange     ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1014cc71bff1SMichael Lange     ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr);
1015cc71bff1SMichael Lange   }
1016cc71bff1SMichael Lange   ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr);
1017cc71bff1SMichael Lange   ierr = DMPlexGetConeOrientations(dmParallel, &newCones);CHKERRQ(ierr);
1018cc71bff1SMichael Lange   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
1019cc71bff1SMichael Lange   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
1020cc71bff1SMichael Lange   ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr);
1021cc71bff1SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
1022cc71bff1SMichael Lange   /* Create supports and stratify sieve */
1023cc71bff1SMichael Lange   {
1024cc71bff1SMichael Lange     PetscInt pStart, pEnd;
1025cc71bff1SMichael Lange 
1026cc71bff1SMichael Lange     ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr);
1027cc71bff1SMichael Lange     ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr);
1028cc71bff1SMichael Lange   }
1029cc71bff1SMichael Lange   ierr = DMPlexSymmetrize(dmParallel);CHKERRQ(ierr);
1030cc71bff1SMichael Lange   ierr = DMPlexStratify(dmParallel);CHKERRQ(ierr);
1031cc71bff1SMichael Lange   PetscFunctionReturn(0);
1032cc71bff1SMichael Lange }
1033cc71bff1SMichael Lange 
1034cc71bff1SMichael Lange #undef __FUNCT__
10350df0e737SMichael Lange #define __FUNCT__ "DMPlexDistributeCoordinates"
10360df0e737SMichael Lange PetscErrorCode DMPlexDistributeCoordinates(DM dm, PetscSF migrationSF, DM dmParallel)
10370df0e737SMichael Lange {
10380df0e737SMichael Lange   MPI_Comm         comm;
10390df0e737SMichael Lange   PetscSection     originalCoordSection, newCoordSection;
10400df0e737SMichael Lange   Vec              originalCoordinates, newCoordinates;
10410df0e737SMichael Lange   PetscInt         bs;
10420df0e737SMichael Lange   const char      *name;
10430df0e737SMichael Lange   const PetscReal *maxCell, *L;
10440df0e737SMichael Lange   PetscErrorCode   ierr;
10450df0e737SMichael Lange 
10460df0e737SMichael Lange   PetscFunctionBegin;
10470df0e737SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10480df0e737SMichael Lange   PetscValidPointer(dmParallel, 3);
10490df0e737SMichael Lange 
10500df0e737SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
10510df0e737SMichael Lange   ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr);
10520df0e737SMichael Lange   ierr = DMGetCoordinateSection(dmParallel, &newCoordSection);CHKERRQ(ierr);
10530df0e737SMichael Lange   ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr);
10540df0e737SMichael Lange   if (originalCoordinates) {
10550df0e737SMichael Lange     ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr);
10560df0e737SMichael Lange     ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr);
10570df0e737SMichael Lange     ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr);
10580df0e737SMichael Lange 
10590df0e737SMichael Lange     ierr = DMPlexDistributeField(dm, migrationSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr);
10600df0e737SMichael Lange     ierr = DMSetCoordinatesLocal(dmParallel, newCoordinates);CHKERRQ(ierr);
10610df0e737SMichael Lange     ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr);
10620df0e737SMichael Lange     ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr);
10630df0e737SMichael Lange     ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr);
10640df0e737SMichael Lange   }
10650df0e737SMichael Lange   ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr);
10660df0e737SMichael Lange   if (L) {ierr = DMSetPeriodicity(dmParallel, maxCell, L);CHKERRQ(ierr);}
10670df0e737SMichael Lange   PetscFunctionReturn(0);
10680df0e737SMichael Lange }
10690df0e737SMichael Lange 
10700df0e737SMichael Lange #undef __FUNCT__
10710df0e737SMichael Lange #define __FUNCT__ "DMPlexDistributeLabels"
10722eb1fa14SMichael Lange PetscErrorCode DMPlexDistributeLabels(DM dm, PetscSF migrationSF, DM dmParallel)
10730df0e737SMichael Lange {
10740df0e737SMichael Lange   MPI_Comm       comm;
10750df0e737SMichael Lange   PetscMPIInt    rank;
1076bdd2d751SMichael Lange   PetscInt       numLabels, l;
10770df0e737SMichael Lange   PetscErrorCode ierr;
10780df0e737SMichael Lange 
10790df0e737SMichael Lange   PetscFunctionBegin;
10800df0e737SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10812eb1fa14SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
10820df0e737SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
10830df0e737SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
10840df0e737SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
10850df0e737SMichael Lange 
10860df0e737SMichael Lange   /* Bcast number of labels */
1087bdd2d751SMichael Lange   ierr = DMPlexGetNumLabels(dm, &numLabels);CHKERRQ(ierr);
10880df0e737SMichael Lange   ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1089bdd2d751SMichael Lange   for (l = numLabels-1; l >= 0; --l) {
1090bdd2d751SMichael Lange     DMLabel     label = NULL, labelNew = NULL;
10910df0e737SMichael Lange     PetscBool   isdepth;
10920df0e737SMichael Lange 
1093bdd2d751SMichael Lange     if (!rank) {
1094bdd2d751SMichael Lange       ierr = DMPlexGetLabelByNum(dm, l, &label);CHKERRQ(ierr);
10950df0e737SMichael Lange       /* Skip "depth" because it is recreated */
1096bdd2d751SMichael Lange       ierr = PetscStrcmp(label->name, "depth", &isdepth);CHKERRQ(ierr);
1097bdd2d751SMichael Lange     }
10980df0e737SMichael Lange     ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
1099bdd2d751SMichael Lange     if (isdepth) continue;
1100bdd2d751SMichael Lange     ierr = DMLabelDistribute(label, migrationSF, &labelNew);CHKERRQ(ierr);
1101bdd2d751SMichael Lange     ierr = DMPlexAddLabel(dmParallel, labelNew);CHKERRQ(ierr);
11020df0e737SMichael Lange   }
11030df0e737SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
11040df0e737SMichael Lange   PetscFunctionReturn(0);
11050df0e737SMichael Lange }
11060df0e737SMichael Lange 
11079734c634SMichael Lange #undef __FUNCT__
11089734c634SMichael Lange #define __FUNCT__ "DMPlexDistributeSetupHybrid"
11099734c634SMichael Lange PetscErrorCode DMPlexDistributeSetupHybrid(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
11109734c634SMichael Lange {
11119734c634SMichael Lange   DM_Plex        *mesh  = (DM_Plex*) dm->data;
11129734c634SMichael Lange   DM_Plex        *pmesh = (DM_Plex*) (dmParallel)->data;
11139734c634SMichael Lange   MPI_Comm        comm;
11149734c634SMichael Lange   const PetscInt *gpoints;
11159734c634SMichael Lange   PetscInt        dim, depth, n, d;
11169734c634SMichael Lange   PetscErrorCode  ierr;
11179734c634SMichael Lange 
11189734c634SMichael Lange   PetscFunctionBegin;
11199734c634SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11209734c634SMichael Lange   PetscValidPointer(dmParallel, 4);
11219734c634SMichael Lange 
11229734c634SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
11239734c634SMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
11249734c634SMichael Lange 
11259734c634SMichael Lange   /* Setup hybrid structure */
11269734c634SMichael Lange   for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];}
11279734c634SMichael Lange   ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr);
11289734c634SMichael Lange   ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr);
11299734c634SMichael Lange   ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr);
11309734c634SMichael Lange   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
11319734c634SMichael Lange   for (d = 0; d <= dim; ++d) {
11329734c634SMichael Lange     PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p;
11339734c634SMichael Lange 
11349734c634SMichael Lange     if (pmax < 0) continue;
11359734c634SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr);
11369734c634SMichael Lange     ierr = DMPlexGetDepthStratum(dmParallel, d, NULL, &pEnd);CHKERRQ(ierr);
11379734c634SMichael Lange     ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr);
11389734c634SMichael Lange     for (p = 0; p < n; ++p) {
11399734c634SMichael Lange       const PetscInt point = gpoints[p];
11409734c634SMichael Lange 
11419734c634SMichael Lange       if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax;
11429734c634SMichael Lange     }
11439734c634SMichael Lange     if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax;
11449734c634SMichael Lange     else            pmesh->hybridPointMax[d] = -1;
11459734c634SMichael Lange   }
11469734c634SMichael Lange   ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr);
11479734c634SMichael Lange   PetscFunctionReturn(0);
11489734c634SMichael Lange }
11490df0e737SMichael Lange 
1150a6f36705SMichael Lange #undef __FUNCT__
1151a6f36705SMichael Lange #define __FUNCT__ "DMPlexDistributeSetupTree"
1152a6f36705SMichael Lange PetscErrorCode DMPlexDistributeSetupTree(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
1153a6f36705SMichael Lange {
1154a6f36705SMichael Lange   MPI_Comm        comm;
1155a6f36705SMichael Lange   DM              refTree;
1156a6f36705SMichael Lange   PetscSection    origParentSection, newParentSection;
1157a6f36705SMichael Lange   PetscInt        *origParents, *origChildIDs;
1158a6f36705SMichael Lange   PetscBool       flg;
1159a6f36705SMichael Lange   PetscErrorCode  ierr;
1160a6f36705SMichael Lange 
1161a6f36705SMichael Lange   PetscFunctionBegin;
1162a6f36705SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1163a6f36705SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 4);
1164a6f36705SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
1165a6f36705SMichael Lange 
1166a6f36705SMichael Lange   /* Set up tree */
1167a6f36705SMichael Lange   ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr);
1168a6f36705SMichael Lange   ierr = DMPlexSetReferenceTree(dmParallel,refTree);CHKERRQ(ierr);
1169a6f36705SMichael Lange   ierr = DMPlexGetTree(dm,&origParentSection,&origParents,&origChildIDs,NULL,NULL);CHKERRQ(ierr);
1170a6f36705SMichael Lange   if (origParentSection) {
1171a6f36705SMichael Lange     PetscInt        pStart, pEnd;
1172a6f36705SMichael Lange     PetscInt        *newParents, *newChildIDs;
1173a6f36705SMichael Lange     PetscInt        *remoteOffsetsParents, newParentSize;
1174a6f36705SMichael Lange     PetscSF         parentSF;
1175a6f36705SMichael Lange 
1176a6f36705SMichael Lange     ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
1177a6f36705SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dmParallel),&newParentSection);CHKERRQ(ierr);
1178a6f36705SMichael Lange     ierr = PetscSectionSetChart(newParentSection,pStart,pEnd);CHKERRQ(ierr);
1179a6f36705SMichael Lange     ierr = PetscSFDistributeSection(migrationSF, origParentSection, &remoteOffsetsParents, newParentSection);CHKERRQ(ierr);
1180a6f36705SMichael Lange     ierr = PetscSFCreateSectionSF(migrationSF, origParentSection, remoteOffsetsParents, newParentSection, &parentSF);CHKERRQ(ierr);
1181a6f36705SMichael Lange     ierr = PetscSectionGetStorageSize(newParentSection,&newParentSize);CHKERRQ(ierr);
1182a6f36705SMichael Lange     ierr = PetscMalloc2(newParentSize,&newParents,newParentSize,&newChildIDs);CHKERRQ(ierr);
1183a6f36705SMichael Lange     ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
1184a6f36705SMichael Lange     ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
1185a6f36705SMichael Lange     ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1186a6f36705SMichael Lange     ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1187a6f36705SMichael Lange     ierr = ISGlobalToLocalMappingApplyBlock(renumbering,IS_GTOLM_MASK, newParentSize, newParents, NULL, newParents);CHKERRQ(ierr);
1188a6f36705SMichael Lange     ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-parents_view", &flg);CHKERRQ(ierr);
1189a6f36705SMichael Lange     if (flg) {
1190a6f36705SMichael Lange       ierr = PetscPrintf(comm, "Serial Parent Section: \n");CHKERRQ(ierr);
1191a6f36705SMichael Lange       ierr = PetscSectionView(origParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1192a6f36705SMichael Lange       ierr = PetscPrintf(comm, "Parallel Parent Section: \n");CHKERRQ(ierr);
1193a6f36705SMichael Lange       ierr = PetscSectionView(newParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1194a6f36705SMichael Lange       ierr = PetscSFView(parentSF, NULL);CHKERRQ(ierr);
1195a6f36705SMichael Lange     }
1196a6f36705SMichael Lange     ierr = DMPlexSetTree(dmParallel,newParentSection,newParents,newChildIDs);CHKERRQ(ierr);
1197a6f36705SMichael Lange     ierr = PetscSectionDestroy(&newParentSection);CHKERRQ(ierr);
1198a6f36705SMichael Lange     ierr = PetscFree2(newParents,newChildIDs);CHKERRQ(ierr);
1199a6f36705SMichael Lange     ierr = PetscSFDestroy(&parentSF);CHKERRQ(ierr);
1200a6f36705SMichael Lange   }
1201a6f36705SMichael Lange   PetscFunctionReturn(0);
1202a6f36705SMichael Lange }
12030df0e737SMichael Lange 
12040df0e737SMichael Lange #undef __FUNCT__
12054eca1733SMichael Lange #define __FUNCT__ "DMPlexDistributeSF"
12064eca1733SMichael Lange PetscErrorCode DMPlexDistributeSF(DM dm, PetscSF migrationSF, PetscSection partSection, IS part, PetscSection origPartSection, IS origPart, DM dmParallel)
12074eca1733SMichael Lange {
12084eca1733SMichael Lange   DM_Plex               *mesh  = (DM_Plex*) dm->data;
12094eca1733SMichael Lange   DM_Plex               *pmesh = (DM_Plex*) (dmParallel)->data;
12104eca1733SMichael Lange   PetscMPIInt            rank, numProcs;
12114eca1733SMichael Lange   MPI_Comm               comm;
12124eca1733SMichael Lange   PetscErrorCode         ierr;
12134eca1733SMichael Lange 
12144eca1733SMichael Lange   PetscFunctionBegin;
12154eca1733SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12164eca1733SMichael Lange   PetscValidPointer(dmParallel,7);
12174eca1733SMichael Lange 
12184eca1733SMichael Lange   /* Create point SF for parallel mesh */
12194eca1733SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
12204eca1733SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
12214eca1733SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
12224eca1733SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
12234eca1733SMichael Lange   {
12244eca1733SMichael Lange     const PetscInt *leaves;
12254eca1733SMichael Lange     PetscSFNode    *remotePoints, *rowners, *lowners;
12264eca1733SMichael Lange     PetscInt        numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints;
12274eca1733SMichael Lange     PetscInt        pStart, pEnd;
12284eca1733SMichael Lange 
12294eca1733SMichael Lange     ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
12304eca1733SMichael Lange     ierr = PetscSFGetGraph(migrationSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr);
12314eca1733SMichael Lange     ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr);
12324eca1733SMichael Lange     for (p=0; p<numRoots; p++) {
12334eca1733SMichael Lange       rowners[p].rank  = -1;
12344eca1733SMichael Lange       rowners[p].index = -1;
12354eca1733SMichael Lange     }
12364eca1733SMichael Lange     if (origPart) {
12374eca1733SMichael Lange       /* Make sure points in the original partition are not assigned to other procs */
12384eca1733SMichael Lange       const PetscInt *origPoints;
12394eca1733SMichael Lange 
12404eca1733SMichael Lange       ierr = ISGetIndices(origPart, &origPoints);CHKERRQ(ierr);
12414eca1733SMichael Lange       for (p = 0; p < numProcs; ++p) {
12424eca1733SMichael Lange         PetscInt dof, off, d;
12434eca1733SMichael Lange 
12444eca1733SMichael Lange         ierr = PetscSectionGetDof(origPartSection, p, &dof);CHKERRQ(ierr);
12454eca1733SMichael Lange         ierr = PetscSectionGetOffset(origPartSection, p, &off);CHKERRQ(ierr);
12464eca1733SMichael Lange         for (d = off; d < off+dof; ++d) {
12474eca1733SMichael Lange           rowners[origPoints[d]].rank = p;
12484eca1733SMichael Lange         }
12494eca1733SMichael Lange       }
12504eca1733SMichael Lange       ierr = ISRestoreIndices(origPart, &origPoints);CHKERRQ(ierr);
12514eca1733SMichael Lange     }
12524eca1733SMichael Lange     ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12534eca1733SMichael Lange     ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12544eca1733SMichael Lange     for (p = 0; p < numLeaves; ++p) {
12554eca1733SMichael Lange       if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */
12564eca1733SMichael Lange         lowners[p].rank  = rank;
12574eca1733SMichael Lange         lowners[p].index = leaves ? leaves[p] : p;
12584eca1733SMichael Lange       } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */
12594eca1733SMichael Lange         lowners[p].rank  = -2;
12604eca1733SMichael Lange         lowners[p].index = -2;
12614eca1733SMichael Lange       }
12624eca1733SMichael Lange     }
12634eca1733SMichael Lange     for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */
12644eca1733SMichael Lange       rowners[p].rank  = -3;
12654eca1733SMichael Lange       rowners[p].index = -3;
12664eca1733SMichael Lange     }
12674eca1733SMichael Lange     ierr = PetscSFReduceBegin(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
12684eca1733SMichael Lange     ierr = PetscSFReduceEnd(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
12694eca1733SMichael Lange     ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12704eca1733SMichael Lange     ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12714eca1733SMichael Lange     for (p = 0; p < numLeaves; ++p) {
12724eca1733SMichael Lange       if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed");
12734eca1733SMichael Lange       if (lowners[p].rank != rank) ++numGhostPoints;
12744eca1733SMichael Lange     }
12754eca1733SMichael Lange     ierr = PetscMalloc1(numGhostPoints, &ghostPoints);CHKERRQ(ierr);
12764eca1733SMichael Lange     ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr);
12774eca1733SMichael Lange     for (p = 0, gp = 0; p < numLeaves; ++p) {
12784eca1733SMichael Lange       if (lowners[p].rank != rank) {
12794eca1733SMichael Lange         ghostPoints[gp]        = leaves ? leaves[p] : p;
12804eca1733SMichael Lange         remotePoints[gp].rank  = lowners[p].rank;
12814eca1733SMichael Lange         remotePoints[gp].index = lowners[p].index;
12824eca1733SMichael Lange         ++gp;
12834eca1733SMichael Lange       }
12844eca1733SMichael Lange     }
12854eca1733SMichael Lange     ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr);
12864eca1733SMichael Lange     ierr = PetscSFSetGraph((dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
12874eca1733SMichael Lange     ierr = PetscSFSetFromOptions((dmParallel)->sf);CHKERRQ(ierr);
12884eca1733SMichael Lange   }
12894eca1733SMichael Lange   pmesh->useCone    = mesh->useCone;
12904eca1733SMichael Lange   pmesh->useClosure = mesh->useClosure;
12914eca1733SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
12924eca1733SMichael Lange   PetscFunctionReturn(0);
12934eca1733SMichael Lange }
12944eca1733SMichael Lange 
12954eca1733SMichael Lange #undef __FUNCT__
129670034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistribute"
129770034214SMatthew G. Knepley /*@C
129870034214SMatthew G. Knepley   DMPlexDistribute - Distributes the mesh and any associated sections.
129970034214SMatthew G. Knepley 
130070034214SMatthew G. Knepley   Not Collective
130170034214SMatthew G. Knepley 
130270034214SMatthew G. Knepley   Input Parameter:
130370034214SMatthew G. Knepley + dm  - The original DMPlex object
130470034214SMatthew G. Knepley - overlap - The overlap of partitions, 0 is the default
130570034214SMatthew G. Knepley 
130670034214SMatthew G. Knepley   Output Parameter:
130770034214SMatthew G. Knepley + sf - The PetscSF used for point distribution
130870034214SMatthew G. Knepley - parallelMesh - The distributed DMPlex object, or NULL
130970034214SMatthew G. Knepley 
1310a9c22940SMatthew G. Knepley   Note: If the mesh was not distributed, the return value is NULL.
1311a9c22940SMatthew G. Knepley 
1312a9c22940SMatthew G. Knepley   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
1313a9c22940SMatthew G. Knepley   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
1314a9c22940SMatthew G. Knepley   representation on the mesh.
131570034214SMatthew G. Knepley 
131670034214SMatthew G. Knepley   Level: intermediate
131770034214SMatthew G. Knepley 
131870034214SMatthew G. Knepley .keywords: mesh, elements
1319a9c22940SMatthew G. Knepley .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
132070034214SMatthew G. Knepley @*/
132180cf41d5SMatthew G. Knepley PetscErrorCode DMPlexDistribute(DM dm, PetscInt overlap, PetscSF *sf, DM *dmParallel)
132270034214SMatthew G. Knepley {
132370034214SMatthew G. Knepley   DM_Plex               *mesh   = (DM_Plex*) dm->data, *pmesh;
132470034214SMatthew G. Knepley   MPI_Comm               comm;
132543331d4aSMichael Lange   PetscInt               dim, numRemoteRanks, nroots, nleaves;
1326a157612eSMichael Lange   DM                     dmOverlap;
1327a157612eSMichael Lange   IS                     cellPart,        part;
1328a157612eSMichael Lange   PetscSection           cellPartSection, partSection;
132943331d4aSMichael Lange   PetscSFNode           *remoteRanks, *newRemote;
133043331d4aSMichael Lange   const PetscSFNode     *oldRemote;
133143331d4aSMichael Lange   PetscSF                partSF, pointSF, overlapPointSF, overlapSF;
133270034214SMatthew G. Knepley   ISLocalToGlobalMapping renumbering;
133370034214SMatthew G. Knepley   PetscBool              flg;
133470034214SMatthew G. Knepley   PetscMPIInt            rank, numProcs, p;
133570034214SMatthew G. Knepley   PetscErrorCode         ierr;
133670034214SMatthew G. Knepley 
133770034214SMatthew G. Knepley   PetscFunctionBegin;
133870034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
133970034214SMatthew G. Knepley   if (sf) PetscValidPointer(sf,4);
134070034214SMatthew G. Knepley   PetscValidPointer(dmParallel,5);
134170034214SMatthew G. Knepley 
134270034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
134370034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
134470034214SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
134570034214SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
134670034214SMatthew G. Knepley 
134770034214SMatthew G. Knepley   *dmParallel = NULL;
134870034214SMatthew G. Knepley   if (numProcs == 1) PetscFunctionReturn(0);
134970034214SMatthew G. Knepley 
1350c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
135170034214SMatthew G. Knepley   /* Create cell partition - We need to rewrite to use IS, use the MatPartition stuff */
135277623264SMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
135377623264SMatthew G. Knepley   if (overlap > 1) SETERRQ(comm, PETSC_ERR_SUP, "Overlap > 1 not yet implemented");
135477623264SMatthew G. Knepley   ierr = PetscSectionCreate(comm, &cellPartSection);CHKERRQ(ierr);
1355a157612eSMichael Lange   ierr = PetscPartitionerPartition(mesh->partitioner, dm, PETSC_FALSE, cellPartSection, &cellPart, NULL, NULL);CHKERRQ(ierr);
135670034214SMatthew G. Knepley   /* Create SF assuming a serial partition for all processes: Could check for IS length here */
135770034214SMatthew G. Knepley   if (!rank) numRemoteRanks = numProcs;
135870034214SMatthew G. Knepley   else       numRemoteRanks = 0;
135970034214SMatthew G. Knepley   ierr = PetscMalloc1(numRemoteRanks, &remoteRanks);CHKERRQ(ierr);
136070034214SMatthew G. Knepley   for (p = 0; p < numRemoteRanks; ++p) {
136170034214SMatthew G. Knepley     remoteRanks[p].rank  = p;
136270034214SMatthew G. Knepley     remoteRanks[p].index = 0;
136370034214SMatthew G. Knepley   }
136470034214SMatthew G. Knepley   ierr = PetscSFCreate(comm, &partSF);CHKERRQ(ierr);
136570034214SMatthew G. Knepley   ierr = PetscSFSetGraph(partSF, 1, numRemoteRanks, NULL, PETSC_OWN_POINTER, remoteRanks, PETSC_OWN_POINTER);CHKERRQ(ierr);
136670034214SMatthew G. Knepley   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
136770034214SMatthew G. Knepley   if (flg) {
1368a157612eSMichael Lange     ierr = PetscPrintf(comm, "Original Cell Partition:\n");CHKERRQ(ierr);
136970034214SMatthew G. Knepley     ierr = PetscSectionView(cellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
137070034214SMatthew G. Knepley     ierr = ISView(cellPart, NULL);CHKERRQ(ierr);
137170034214SMatthew G. Knepley     ierr = PetscSFView(partSF, NULL);CHKERRQ(ierr);
137270034214SMatthew G. Knepley   }
137370034214SMatthew G. Knepley   /* Close the partition over the mesh */
137470034214SMatthew G. Knepley   ierr = DMPlexCreatePartitionClosure(dm, cellPartSection, cellPart, &partSection, &part);CHKERRQ(ierr);
137570034214SMatthew G. Knepley   /* Create new mesh */
137670034214SMatthew G. Knepley   ierr  = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr);
1377c73cfb54SMatthew G. Knepley   ierr  = DMSetDimension(*dmParallel, dim);CHKERRQ(ierr);
137870034214SMatthew G. Knepley   ierr  = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr);
137970034214SMatthew G. Knepley   pmesh = (DM_Plex*) (*dmParallel)->data;
13804eca1733SMichael Lange   pmesh->useAnchors = mesh->useAnchors;
13814eca1733SMichael Lange 
138270034214SMatthew G. Knepley   /* Distribute sieve points and the global point numbering (replaces creating remote bases) */
138370034214SMatthew G. Knepley   ierr = PetscSFConvertPartition(partSF, partSection, part, &renumbering, &pointSF);CHKERRQ(ierr);
138470034214SMatthew G. Knepley   if (flg) {
138570034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Point Partition:\n");CHKERRQ(ierr);
138670034214SMatthew G. Knepley     ierr = PetscSectionView(partSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
138770034214SMatthew G. Knepley     ierr = ISView(part, NULL);CHKERRQ(ierr);
138870034214SMatthew G. Knepley     ierr = PetscSFView(pointSF, NULL);CHKERRQ(ierr);
138970034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Point Renumbering after partition:\n");CHKERRQ(ierr);
139070034214SMatthew G. Knepley     ierr = ISLocalToGlobalMappingView(renumbering, NULL);CHKERRQ(ierr);
139170034214SMatthew G. Knepley   }
139277623264SMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
139370034214SMatthew G. Knepley 
1394a157612eSMichael Lange   /* Migrate data to a non-overlapping parallel DM */
1395cc71bff1SMichael Lange   ierr = DMPlexDistributeCones(dm, pointSF, renumbering, *dmParallel);CHKERRQ(ierr);
13960df0e737SMichael Lange   ierr = DMPlexDistributeCoordinates(dm, pointSF, *dmParallel);CHKERRQ(ierr);
13972eb1fa14SMichael Lange   ierr = DMPlexDistributeLabels(dm, pointSF, *dmParallel);CHKERRQ(ierr);
13989734c634SMichael Lange   ierr = DMPlexDistributeSetupHybrid(dm, pointSF, renumbering, *dmParallel);CHKERRQ(ierr);
1399a6f36705SMichael Lange   ierr = DMPlexDistributeSetupTree(dm, pointSF, renumbering, *dmParallel);CHKERRQ(ierr);
140070034214SMatthew G. Knepley 
1401a157612eSMichael Lange   /* Build the point SF without overlap */
1402a157612eSMichael Lange   ierr = DMPlexDistributeSF(dm, pointSF, partSection, part, NULL, NULL, *dmParallel);CHKERRQ(ierr);
140315fff7beSMatthew G. Knepley   if (flg) {
140415fff7beSMatthew G. Knepley     ierr = PetscPrintf(comm, "Point SF:\n");CHKERRQ(ierr);
140515fff7beSMatthew G. Knepley     ierr = PetscSFView((*dmParallel)->sf, NULL);CHKERRQ(ierr);
140615fff7beSMatthew G. Knepley   }
1407a157612eSMichael Lange 
1408a157612eSMichael Lange   if (overlap > 0) {
14093d822a50SMichael Lange     ierr = PetscLogEventBegin(DMPLEX_DistributeOverlap,dm,0,0,0);CHKERRQ(ierr);
1410a157612eSMichael Lange     /* Add the partition overlap to the distributed DM */
1411a157612eSMichael Lange     ierr = DMPlexDistributeOverlap(*dmParallel, overlap, renumbering, &overlapSF, &dmOverlap);CHKERRQ(ierr);
1412a157612eSMichael Lange     ierr = DMDestroy(dmParallel);CHKERRQ(ierr);
1413a157612eSMichael Lange     *dmParallel = dmOverlap;
14143d822a50SMichael Lange     if (flg) {
14153d822a50SMichael Lange       ierr = PetscPrintf(comm, "Overlap Migration SF:\n");CHKERRQ(ierr);
14163d822a50SMichael Lange       ierr = PetscSFView(overlapSF, NULL);CHKERRQ(ierr);
14173d822a50SMichael Lange     }
141843331d4aSMichael Lange 
141943331d4aSMichael Lange     /* Re-map the pointSF to establish the full migration pattern */
142043331d4aSMichael Lange     ierr = PetscSFGetGraph(pointSF, &nroots, NULL, NULL, &oldRemote);CHKERRQ(ierr);
142143331d4aSMichael Lange     ierr = PetscSFGetGraph(overlapSF, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr);
142243331d4aSMichael Lange     ierr = PetscMalloc1(nleaves, &newRemote);CHKERRQ(ierr);
142343331d4aSMichael Lange     ierr = PetscSFBcastBegin(overlapSF, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr);
142443331d4aSMichael Lange     ierr = PetscSFBcastEnd(overlapSF, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr);
142543331d4aSMichael Lange     ierr = PetscSFCreate(comm, &overlapPointSF);CHKERRQ(ierr);
142643331d4aSMichael Lange     ierr = PetscSFSetGraph(overlapPointSF, nroots, nleaves, NULL, PETSC_OWN_POINTER, newRemote, PETSC_OWN_POINTER);CHKERRQ(ierr);
142715fff7beSMatthew G. Knepley     ierr = PetscSFDestroy(&overlapSF);CHKERRQ(ierr);
142843331d4aSMichael Lange     ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr);
142943331d4aSMichael Lange     pointSF = overlapPointSF;
14303d822a50SMichael Lange     ierr = PetscLogEventEnd(DMPLEX_DistributeOverlap,dm,0,0,0);CHKERRQ(ierr);
1431a157612eSMichael Lange   }
1432bf5244c7SMatthew G. Knepley   /* Cleanup Partition */
1433bf5244c7SMatthew G. Knepley   ierr = ISLocalToGlobalMappingDestroy(&renumbering);CHKERRQ(ierr);
1434bf5244c7SMatthew G. Knepley   ierr = PetscSFDestroy(&partSF);CHKERRQ(ierr);
1435bf5244c7SMatthew G. Knepley   ierr = PetscSectionDestroy(&partSection);CHKERRQ(ierr);
1436bf5244c7SMatthew G. Knepley   ierr = ISDestroy(&part);CHKERRQ(ierr);
14374eca1733SMichael Lange   ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr);
14384eca1733SMichael Lange   ierr = ISDestroy(&cellPart);CHKERRQ(ierr);
143986719b60SMatthew G. Knepley   /* Copy BC */
144086719b60SMatthew G. Knepley   ierr = DMPlexCopyBoundary(dm, *dmParallel);CHKERRQ(ierr);
144170034214SMatthew G. Knepley   /* Cleanup */
144270034214SMatthew G. Knepley   if (sf) {*sf = pointSF;}
144370034214SMatthew G. Knepley   else    {ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr);}
144470034214SMatthew G. Knepley   ierr = DMSetFromOptions(*dmParallel);CHKERRQ(ierr);
144570034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
144670034214SMatthew G. Knepley   PetscFunctionReturn(0);
144770034214SMatthew G. Knepley }
1448a157612eSMichael Lange 
1449a157612eSMichael Lange #undef __FUNCT__
1450a157612eSMichael Lange #define __FUNCT__ "DMPlexDistributeOverlap"
1451a157612eSMichael Lange /*@C
1452a157612eSMichael Lange   DMPlexDistribute - Add partition overlap to a distributed non-overlapping DM.
1453a157612eSMichael Lange 
1454a157612eSMichael Lange   Not Collective
1455a157612eSMichael Lange 
1456a157612eSMichael Lange   Input Parameter:
1457a157612eSMichael Lange + dm  - The non-overlapping distrbuted DMPlex object
1458a157612eSMichael Lange - overlap - The overlap of partitions, 0 is the default
1459a157612eSMichael Lange 
1460a157612eSMichael Lange   Output Parameter:
1461a157612eSMichael Lange + sf - The PetscSF used for point distribution
1462a157612eSMichael Lange - dmOverlap - The overlapping distributed DMPlex object, or NULL
1463a157612eSMichael Lange 
1464a157612eSMichael Lange   Note: If the mesh was not distributed, the return value is NULL.
1465a157612eSMichael Lange 
1466a157612eSMichael Lange   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
1467a157612eSMichael Lange   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
1468a157612eSMichael Lange   representation on the mesh.
1469a157612eSMichael Lange 
1470a157612eSMichael Lange   Level: intermediate
1471a157612eSMichael Lange 
1472a157612eSMichael Lange .keywords: mesh, elements
1473a157612eSMichael Lange .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
1474a157612eSMichael Lange @*/
1475a157612eSMichael Lange PetscErrorCode DMPlexDistributeOverlap(DM dm, PetscInt overlap, ISLocalToGlobalMapping renumbering, PetscSF *sf, DM *dmOverlap)
1476a157612eSMichael Lange {
1477a157612eSMichael Lange   MPI_Comm               comm;
1478a157612eSMichael Lange   PetscMPIInt            rank;
1479a157612eSMichael Lange   PetscSection           rootSection, leafSection;
1480a157612eSMichael Lange   IS                     rootrank, leafrank;
1481a157612eSMichael Lange   PetscSection           coneSection;
1482a157612eSMichael Lange   PetscSF                overlapSF, migrationSF, pointSF, newPointSF;
1483a157612eSMichael Lange   PetscSFNode           *ghostRemote;
1484a157612eSMichael Lange   const PetscSFNode     *overlapRemote;
1485a157612eSMichael Lange   ISLocalToGlobalMapping overlapRenumbering;
1486a157612eSMichael Lange   const PetscInt        *renumberingArray, *overlapLocal;
1487a157612eSMichael Lange   PetscInt               dim, p, pStart, pEnd, conesSize, idx;
1488a157612eSMichael Lange   PetscInt               numGhostPoints, numOverlapPoints, numSharedPoints, overlapLeaves;
1489a157612eSMichael Lange   PetscInt              *cones, *ghostLocal, *overlapRenumberingArray, *pointIDs, *recvPointIDs;
1490a157612eSMichael Lange   PetscErrorCode         ierr;
1491a157612eSMichael Lange 
1492a157612eSMichael Lange   PetscFunctionBegin;
1493a157612eSMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1494a157612eSMichael Lange   if (sf) PetscValidPointer(sf, 3);
1495a157612eSMichael Lange   PetscValidPointer(dmOverlap, 4);
1496a157612eSMichael Lange 
1497a157612eSMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1498a157612eSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1499a157612eSMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1500a157612eSMichael Lange 
1501a157612eSMichael Lange   /* Compute point overlap with neighbouring processes on the distributed DM */
1502a157612eSMichael Lange   ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1503a157612eSMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
1504a157612eSMichael Lange   ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr);
1505a157612eSMichael Lange   ierr = DMPlexDistributeOwnership(dm, rootSection, &rootrank, leafSection, &leafrank);CHKERRQ(ierr);
1506a157612eSMichael Lange   ierr = DMPlexCreateOverlap(dm, rootSection, rootrank, leafSection, leafrank, &overlapSF);CHKERRQ(ierr);
150715fff7beSMatthew G. Knepley   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
150815fff7beSMatthew G. Knepley   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
150915fff7beSMatthew G. Knepley   ierr = ISDestroy(&rootrank);CHKERRQ(ierr);
151015fff7beSMatthew G. Knepley   ierr = ISDestroy(&leafrank);CHKERRQ(ierr);
1511a157612eSMichael Lange   ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1512a157612eSMichael Lange 
1513a157612eSMichael Lange   /* Build dense migration SF that maps the non-overlapping partition to the overlapping one */
1514a157612eSMichael Lange   ierr = DMPlexCreateOverlapMigrationSF(dm, overlapSF, &migrationSF);CHKERRQ(ierr);
1515a157612eSMichael Lange 
1516a157612eSMichael Lange   /* Convert cones to global numbering before migrating them */
1517a157612eSMichael Lange   ierr = DMPlexGetConeSection(dm, &coneSection);CHKERRQ(ierr);
1518a157612eSMichael Lange   ierr = PetscSectionGetStorageSize(coneSection, &conesSize);CHKERRQ(ierr);
1519a157612eSMichael Lange   ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
1520a157612eSMichael Lange   ierr = ISLocalToGlobalMappingApplyBlock(renumbering, conesSize, cones, cones);CHKERRQ(ierr);
1521a157612eSMichael Lange 
1522a157612eSMichael Lange   /* Derive the new local-to-global mapping from the old one */
1523a157612eSMichael Lange   ierr = PetscSFGetGraph(migrationSF, NULL, &overlapLeaves, &overlapLocal, &overlapRemote);CHKERRQ(ierr);
1524a157612eSMichael Lange   ierr = PetscMalloc1(overlapLeaves, &overlapRenumberingArray);CHKERRQ(ierr);
1525a157612eSMichael Lange   ierr = ISLocalToGlobalMappingGetBlockIndices(renumbering, &renumberingArray);CHKERRQ(ierr);
1526729b3788SMatthew G. Knepley   ierr = PetscSFBcastBegin(migrationSF, MPIU_INT, (PetscInt *) renumberingArray, overlapRenumberingArray);CHKERRQ(ierr);
1527729b3788SMatthew G. Knepley   ierr = PetscSFBcastEnd(migrationSF, MPIU_INT, (PetscInt *) renumberingArray, overlapRenumberingArray);CHKERRQ(ierr);
1528a157612eSMichael Lange   ierr = ISLocalToGlobalMappingCreate(comm, 1, overlapLeaves, (const PetscInt*) overlapRenumberingArray, PETSC_OWN_POINTER, &overlapRenumbering);CHKERRQ(ierr);
1529a157612eSMichael Lange 
1530a157612eSMichael Lange   /* Build the overlapping DM */
1531a157612eSMichael Lange   ierr = DMPlexCreate(comm, dmOverlap);CHKERRQ(ierr);
1532a157612eSMichael Lange   ierr = DMSetDimension(*dmOverlap, dim);CHKERRQ(ierr);
1533a157612eSMichael Lange   ierr = PetscObjectSetName((PetscObject) *dmOverlap, "Parallel Mesh");CHKERRQ(ierr);
1534a157612eSMichael Lange   ierr = DMPlexDistributeCones(dm, migrationSF, overlapRenumbering, *dmOverlap);CHKERRQ(ierr);
1535a157612eSMichael Lange   ierr = DMPlexDistributeCoordinates(dm, migrationSF, *dmOverlap);CHKERRQ(ierr);
1536a157612eSMichael Lange   ierr = DMPlexDistributeLabels(dm, migrationSF, *dmOverlap);CHKERRQ(ierr);
1537a157612eSMichael Lange   ierr = DMPlexDistributeSetupHybrid(dm, migrationSF, overlapRenumbering, *dmOverlap);CHKERRQ(ierr);
1538a157612eSMichael Lange 
1539a157612eSMichael Lange   /* Build the new point SF by propagating the depthShift generate remote root indices */
1540a157612eSMichael Lange   ierr = DMGetPointSF(dm, &pointSF);CHKERRQ(ierr);
1541a157612eSMichael Lange   ierr = PetscSFGetGraph(pointSF, NULL, &numSharedPoints, NULL, NULL);CHKERRQ(ierr);
1542a157612eSMichael Lange   ierr = PetscSFGetGraph(overlapSF, NULL, &numOverlapPoints, NULL, NULL);CHKERRQ(ierr);
1543a157612eSMichael Lange   numGhostPoints = numSharedPoints + numOverlapPoints;
154409b7985cSMatthew G. Knepley   ierr = PetscMalloc1(numGhostPoints, &ghostLocal);CHKERRQ(ierr);
154509b7985cSMatthew G. Knepley   ierr = PetscMalloc1(numGhostPoints, &ghostRemote);CHKERRQ(ierr);
1546a157612eSMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
1547a157612eSMichael Lange   ierr = PetscMalloc2(pEnd-pStart, &pointIDs, overlapLeaves, &recvPointIDs);CHKERRQ(ierr);
1548a157612eSMichael Lange   for (p=0; p<overlapLeaves; p++) {
1549a157612eSMichael Lange     if (overlapRemote[p].rank == rank) pointIDs[overlapRemote[p].index] = overlapLocal[p];
1550a157612eSMichael Lange   }
1551a157612eSMichael Lange   ierr = PetscSFBcastBegin(migrationSF, MPIU_INT, pointIDs, recvPointIDs);CHKERRQ(ierr);
1552a157612eSMichael Lange   ierr = PetscSFBcastEnd(migrationSF, MPIU_INT, pointIDs, recvPointIDs);CHKERRQ(ierr);
1553a157612eSMichael Lange   for (idx=0, p=0; p<overlapLeaves; p++) {
1554a157612eSMichael Lange     if (overlapRemote[p].rank != rank) {
1555a157612eSMichael Lange       ghostLocal[idx] = overlapLocal[p];
1556a157612eSMichael Lange       ghostRemote[idx].index = recvPointIDs[p];
1557a157612eSMichael Lange       ghostRemote[idx].rank = overlapRemote[p].rank;
1558a157612eSMichael Lange       idx++;
1559a157612eSMichael Lange     }
1560a157612eSMichael Lange   }
1561a157612eSMichael Lange   ierr = DMPlexGetChart(*dmOverlap, &pStart, &pEnd);CHKERRQ(ierr);
1562a157612eSMichael Lange   ierr = PetscSFCreate(comm, &newPointSF);;CHKERRQ(ierr);
1563a157612eSMichael Lange   ierr = PetscSFSetGraph(newPointSF, pEnd - pStart, numGhostPoints, ghostLocal, PETSC_OWN_POINTER, ghostRemote, PETSC_OWN_POINTER);CHKERRQ(ierr);
1564a157612eSMichael Lange   ierr = DMSetPointSF(*dmOverlap, newPointSF);CHKERRQ(ierr);
156515fff7beSMatthew G. Knepley   ierr = PetscSFDestroy(&newPointSF);CHKERRQ(ierr);
1566a157612eSMichael Lange   /* Cleanup overlap partition */
1567a157612eSMichael Lange   ierr = ISLocalToGlobalMappingDestroy(&overlapRenumbering);CHKERRQ(ierr);
1568a157612eSMichael Lange   ierr = PetscSFDestroy(&overlapSF);CHKERRQ(ierr);
1569a157612eSMichael Lange   ierr = PetscFree2(pointIDs, recvPointIDs);CHKERRQ(ierr);
1570a157612eSMichael Lange   if (sf) *sf = migrationSF;
1571a157612eSMichael Lange   else    {ierr = PetscSFDestroy(&migrationSF);CHKERRQ(ierr);}
1572a157612eSMichael Lange   ierr = DMSetFromOptions(*dmOverlap);CHKERRQ(ierr);
1573a157612eSMichael Lange   PetscFunctionReturn(0);
1574a157612eSMichael Lange }
1575