xref: /petsc/src/dm/impls/plex/plexdistribute.c (revision 8a713f3bec5a20e7d0b4a32fd244f7bc2dd1ea43)
1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
270034214SMatthew G. Knepley 
370034214SMatthew G. Knepley #undef __FUNCT__
470034214SMatthew G. Knepley #define __FUNCT__ "DMPlexSetAdjacencyUseCone"
570034214SMatthew G. Knepley /*@
670034214SMatthew G. Knepley   DMPlexSetAdjacencyUseCone - Define adjacency in the mesh using either the cone or the support first
770034214SMatthew G. Knepley 
870034214SMatthew G. Knepley   Input Parameters:
970034214SMatthew G. Knepley + dm      - The DM object
1070034214SMatthew G. Knepley - useCone - Flag to use the cone first
1170034214SMatthew G. Knepley 
1270034214SMatthew G. Knepley   Level: intermediate
1370034214SMatthew G. Knepley 
1470034214SMatthew G. Knepley   Notes:
1570034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
164b6b44bdSMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
1770034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
1870034214SMatthew G. Knepley 
1970034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
2070034214SMatthew G. Knepley @*/
2170034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseCone(DM dm, PetscBool useCone)
2270034214SMatthew G. Knepley {
2370034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
2470034214SMatthew G. Knepley 
2570034214SMatthew G. Knepley   PetscFunctionBegin;
2670034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2770034214SMatthew G. Knepley   mesh->useCone = useCone;
2870034214SMatthew G. Knepley   PetscFunctionReturn(0);
2970034214SMatthew G. Knepley }
3070034214SMatthew G. Knepley 
3170034214SMatthew G. Knepley #undef __FUNCT__
3270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacencyUseCone"
3370034214SMatthew G. Knepley /*@
3470034214SMatthew G. Knepley   DMPlexGetAdjacencyUseCone - Query whether adjacency in the mesh uses the cone or the support first
3570034214SMatthew G. Knepley 
3670034214SMatthew G. Knepley   Input Parameter:
3770034214SMatthew G. Knepley . dm      - The DM object
3870034214SMatthew G. Knepley 
3970034214SMatthew G. Knepley   Output Parameter:
4070034214SMatthew G. Knepley . useCone - Flag to use the cone first
4170034214SMatthew G. Knepley 
4270034214SMatthew G. Knepley   Level: intermediate
4370034214SMatthew G. Knepley 
4470034214SMatthew G. Knepley   Notes:
4570034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
464b6b44bdSMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
4770034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
4870034214SMatthew G. Knepley 
4970034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
5070034214SMatthew G. Knepley @*/
5170034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseCone(DM dm, PetscBool *useCone)
5270034214SMatthew G. Knepley {
5370034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
5470034214SMatthew G. Knepley 
5570034214SMatthew G. Knepley   PetscFunctionBegin;
5670034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5770034214SMatthew G. Knepley   PetscValidIntPointer(useCone, 2);
5870034214SMatthew G. Knepley   *useCone = mesh->useCone;
5970034214SMatthew G. Knepley   PetscFunctionReturn(0);
6070034214SMatthew G. Knepley }
6170034214SMatthew G. Knepley 
6270034214SMatthew G. Knepley #undef __FUNCT__
6370034214SMatthew G. Knepley #define __FUNCT__ "DMPlexSetAdjacencyUseClosure"
6470034214SMatthew G. Knepley /*@
6570034214SMatthew G. Knepley   DMPlexSetAdjacencyUseClosure - Define adjacency in the mesh using the transitive closure
6670034214SMatthew G. Knepley 
6770034214SMatthew G. Knepley   Input Parameters:
6870034214SMatthew G. Knepley + dm      - The DM object
6970034214SMatthew G. Knepley - useClosure - Flag to use the closure
7070034214SMatthew G. Knepley 
7170034214SMatthew G. Knepley   Level: intermediate
7270034214SMatthew G. Knepley 
7370034214SMatthew G. Knepley   Notes:
7470034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
754b6b44bdSMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
7670034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
7770034214SMatthew G. Knepley 
7870034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
7970034214SMatthew G. Knepley @*/
8070034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseClosure(DM dm, PetscBool useClosure)
8170034214SMatthew G. Knepley {
8270034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
8370034214SMatthew G. Knepley 
8470034214SMatthew G. Knepley   PetscFunctionBegin;
8570034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8670034214SMatthew G. Knepley   mesh->useClosure = useClosure;
8770034214SMatthew G. Knepley   PetscFunctionReturn(0);
8870034214SMatthew G. Knepley }
8970034214SMatthew G. Knepley 
9070034214SMatthew G. Knepley #undef __FUNCT__
9170034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacencyUseClosure"
9270034214SMatthew G. Knepley /*@
9370034214SMatthew G. Knepley   DMPlexGetAdjacencyUseClosure - Query whether adjacency in the mesh uses the transitive closure
9470034214SMatthew G. Knepley 
9570034214SMatthew G. Knepley   Input Parameter:
9670034214SMatthew G. Knepley . dm      - The DM object
9770034214SMatthew G. Knepley 
9870034214SMatthew G. Knepley   Output Parameter:
9970034214SMatthew G. Knepley . useClosure - Flag to use the closure
10070034214SMatthew G. Knepley 
10170034214SMatthew G. Knepley   Level: intermediate
10270034214SMatthew G. Knepley 
10370034214SMatthew G. Knepley   Notes:
10470034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
1054b6b44bdSMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
10670034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
10770034214SMatthew G. Knepley 
10870034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
10970034214SMatthew G. Knepley @*/
11070034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseClosure(DM dm, PetscBool *useClosure)
11170034214SMatthew G. Knepley {
11270034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
11370034214SMatthew G. Knepley 
11470034214SMatthew G. Knepley   PetscFunctionBegin;
11570034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11670034214SMatthew G. Knepley   PetscValidIntPointer(useClosure, 2);
11770034214SMatthew G. Knepley   *useClosure = mesh->useClosure;
11870034214SMatthew G. Knepley   PetscFunctionReturn(0);
11970034214SMatthew G. Knepley }
12070034214SMatthew G. Knepley 
12170034214SMatthew G. Knepley #undef __FUNCT__
122a17985deSToby Isaac #define __FUNCT__ "DMPlexSetAdjacencyUseAnchors"
1238b0b4c70SToby Isaac /*@
124a17985deSToby Isaac   DMPlexSetAdjacencyUseAnchors - Define adjacency in the mesh using the point-to-point constraints.
1258b0b4c70SToby Isaac 
1268b0b4c70SToby Isaac   Input Parameters:
1278b0b4c70SToby Isaac + dm      - The DM object
1285b317d89SToby 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.
1298b0b4c70SToby Isaac 
1308b0b4c70SToby Isaac   Level: intermediate
1318b0b4c70SToby Isaac 
132a17985deSToby Isaac .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
1338b0b4c70SToby Isaac @*/
1345b317d89SToby Isaac PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM dm, PetscBool useAnchors)
1358b0b4c70SToby Isaac {
1368b0b4c70SToby Isaac   DM_Plex *mesh = (DM_Plex *) dm->data;
1378b0b4c70SToby Isaac 
1388b0b4c70SToby Isaac   PetscFunctionBegin;
1398b0b4c70SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1405b317d89SToby Isaac   mesh->useAnchors = useAnchors;
1418b0b4c70SToby Isaac   PetscFunctionReturn(0);
1428b0b4c70SToby Isaac }
1438b0b4c70SToby Isaac 
1448b0b4c70SToby Isaac #undef __FUNCT__
145a17985deSToby Isaac #define __FUNCT__ "DMPlexGetAdjacencyUseAnchors"
1468b0b4c70SToby Isaac /*@
147a17985deSToby Isaac   DMPlexGetAdjacencyUseAnchors - Query whether adjacency in the mesh uses the point-to-point constraints.
1488b0b4c70SToby Isaac 
1498b0b4c70SToby Isaac   Input Parameter:
1508b0b4c70SToby Isaac . dm      - The DM object
1518b0b4c70SToby Isaac 
1528b0b4c70SToby Isaac   Output Parameter:
1535b317d89SToby 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.
1548b0b4c70SToby Isaac 
1558b0b4c70SToby Isaac   Level: intermediate
1568b0b4c70SToby Isaac 
157a17985deSToby Isaac .seealso: DMPlexSetAdjacencyUseAnchors(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
1588b0b4c70SToby Isaac @*/
1595b317d89SToby Isaac PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM dm, PetscBool *useAnchors)
1608b0b4c70SToby Isaac {
1618b0b4c70SToby Isaac   DM_Plex *mesh = (DM_Plex *) dm->data;
1628b0b4c70SToby Isaac 
1638b0b4c70SToby Isaac   PetscFunctionBegin;
1648b0b4c70SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1655b317d89SToby Isaac   PetscValidIntPointer(useAnchors, 2);
1665b317d89SToby Isaac   *useAnchors = mesh->useAnchors;
1678b0b4c70SToby Isaac   PetscFunctionReturn(0);
1688b0b4c70SToby Isaac }
1698b0b4c70SToby Isaac 
1708b0b4c70SToby Isaac #undef __FUNCT__
17170034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Cone_Internal"
17270034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Cone_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
17370034214SMatthew G. Knepley {
17470034214SMatthew G. Knepley   const PetscInt *cone = NULL;
17570034214SMatthew G. Knepley   PetscInt        numAdj = 0, maxAdjSize = *adjSize, coneSize, c;
17670034214SMatthew G. Knepley   PetscErrorCode  ierr;
17770034214SMatthew G. Knepley 
17870034214SMatthew G. Knepley   PetscFunctionBeginHot;
17970034214SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
18070034214SMatthew G. Knepley   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
1814b6b44bdSMatthew G. Knepley   for (c = 0; c <= coneSize; ++c) {
1824b6b44bdSMatthew G. Knepley     const PetscInt  point   = !c ? p : cone[c-1];
18370034214SMatthew G. Knepley     const PetscInt *support = NULL;
18470034214SMatthew G. Knepley     PetscInt        supportSize, s, q;
18570034214SMatthew G. Knepley 
1864b6b44bdSMatthew G. Knepley     ierr = DMPlexGetSupportSize(dm, point, &supportSize);CHKERRQ(ierr);
1874b6b44bdSMatthew G. Knepley     ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr);
18870034214SMatthew G. Knepley     for (s = 0; s < supportSize; ++s) {
18970034214SMatthew G. Knepley       for (q = 0; q < numAdj || (adj[numAdj++] = support[s],0); ++q) {
19070034214SMatthew G. Knepley         if (support[s] == adj[q]) break;
19170034214SMatthew G. Knepley       }
19270034214SMatthew G. Knepley       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
19370034214SMatthew G. Knepley     }
19470034214SMatthew G. Knepley   }
19570034214SMatthew G. Knepley   *adjSize = numAdj;
19670034214SMatthew G. Knepley   PetscFunctionReturn(0);
19770034214SMatthew G. Knepley }
19870034214SMatthew G. Knepley 
19970034214SMatthew G. Knepley #undef __FUNCT__
20070034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Support_Internal"
20170034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Support_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
20270034214SMatthew G. Knepley {
20370034214SMatthew G. Knepley   const PetscInt *support = NULL;
20470034214SMatthew G. Knepley   PetscInt        numAdj   = 0, maxAdjSize = *adjSize, supportSize, s;
20570034214SMatthew G. Knepley   PetscErrorCode  ierr;
20670034214SMatthew G. Knepley 
20770034214SMatthew G. Knepley   PetscFunctionBeginHot;
20870034214SMatthew G. Knepley   ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr);
20970034214SMatthew G. Knepley   ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr);
2104b6b44bdSMatthew G. Knepley   for (s = 0; s <= supportSize; ++s) {
2114b6b44bdSMatthew G. Knepley     const PetscInt  point = !s ? p : support[s-1];
21270034214SMatthew G. Knepley     const PetscInt *cone  = NULL;
21370034214SMatthew G. Knepley     PetscInt        coneSize, c, q;
21470034214SMatthew G. Knepley 
2154b6b44bdSMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
2164b6b44bdSMatthew G. Knepley     ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
21770034214SMatthew G. Knepley     for (c = 0; c < coneSize; ++c) {
21870034214SMatthew G. Knepley       for (q = 0; q < numAdj || (adj[numAdj++] = cone[c],0); ++q) {
21970034214SMatthew G. Knepley         if (cone[c] == adj[q]) break;
22070034214SMatthew G. Knepley       }
22170034214SMatthew G. Knepley       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
22270034214SMatthew G. Knepley     }
22370034214SMatthew G. Knepley   }
22470034214SMatthew G. Knepley   *adjSize = numAdj;
22570034214SMatthew G. Knepley   PetscFunctionReturn(0);
22670034214SMatthew G. Knepley }
22770034214SMatthew G. Knepley 
22870034214SMatthew G. Knepley #undef __FUNCT__
22970034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Transitive_Internal"
23070034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Transitive_Internal(DM dm, PetscInt p, PetscBool useClosure, PetscInt *adjSize, PetscInt adj[])
23170034214SMatthew G. Knepley {
23270034214SMatthew G. Knepley   PetscInt      *star = NULL;
23370034214SMatthew G. Knepley   PetscInt       numAdj = 0, maxAdjSize = *adjSize, starSize, s;
23470034214SMatthew G. Knepley   PetscErrorCode ierr;
23570034214SMatthew G. Knepley 
23670034214SMatthew G. Knepley   PetscFunctionBeginHot;
23770034214SMatthew G. Knepley   ierr = DMPlexGetTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr);
23870034214SMatthew G. Knepley   for (s = 0; s < starSize*2; s += 2) {
23970034214SMatthew G. Knepley     const PetscInt *closure = NULL;
24070034214SMatthew G. Knepley     PetscInt        closureSize, c, q;
24170034214SMatthew G. Knepley 
24270034214SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
24370034214SMatthew G. Knepley     for (c = 0; c < closureSize*2; c += 2) {
24470034214SMatthew G. Knepley       for (q = 0; q < numAdj || (adj[numAdj++] = closure[c],0); ++q) {
24570034214SMatthew G. Knepley         if (closure[c] == adj[q]) break;
24670034214SMatthew G. Knepley       }
24770034214SMatthew G. Knepley       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
24870034214SMatthew G. Knepley     }
24970034214SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
25070034214SMatthew G. Knepley   }
25170034214SMatthew G. Knepley   ierr = DMPlexRestoreTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr);
25270034214SMatthew G. Knepley   *adjSize = numAdj;
25370034214SMatthew G. Knepley   PetscFunctionReturn(0);
25470034214SMatthew G. Knepley }
25570034214SMatthew G. Knepley 
25670034214SMatthew G. Knepley #undef __FUNCT__
25770034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Internal"
2585b317d89SToby Isaac PetscErrorCode DMPlexGetAdjacency_Internal(DM dm, PetscInt p, PetscBool useCone, PetscBool useTransitiveClosure, PetscBool useAnchors, PetscInt *adjSize, PetscInt *adj[])
25970034214SMatthew G. Knepley {
26079bad331SMatthew G. Knepley   static PetscInt asiz = 0;
2618b0b4c70SToby Isaac   PetscInt maxAnchors = 1;
2628b0b4c70SToby Isaac   PetscInt aStart = -1, aEnd = -1;
2638b0b4c70SToby Isaac   PetscInt maxAdjSize;
2648b0b4c70SToby Isaac   PetscSection aSec = NULL;
2658b0b4c70SToby Isaac   IS aIS = NULL;
2668b0b4c70SToby Isaac   const PetscInt *anchors;
26770034214SMatthew G. Knepley   PetscErrorCode  ierr;
26870034214SMatthew G. Knepley 
26970034214SMatthew G. Knepley   PetscFunctionBeginHot;
2705b317d89SToby Isaac   if (useAnchors) {
271a17985deSToby Isaac     ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr);
2728b0b4c70SToby Isaac     if (aSec) {
2738b0b4c70SToby Isaac       ierr = PetscSectionGetMaxDof(aSec,&maxAnchors);CHKERRQ(ierr);
27424c766afSToby Isaac       maxAnchors = PetscMax(1,maxAnchors);
2758b0b4c70SToby Isaac       ierr = PetscSectionGetChart(aSec,&aStart,&aEnd);CHKERRQ(ierr);
2768b0b4c70SToby Isaac       ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr);
2778b0b4c70SToby Isaac     }
2788b0b4c70SToby Isaac   }
27979bad331SMatthew G. Knepley   if (!*adj) {
28024c766afSToby Isaac     PetscInt depth, coneSeries, supportSeries, maxC, maxS, pStart, pEnd;
28179bad331SMatthew G. Knepley 
28224c766afSToby Isaac     ierr  = DMPlexGetChart(dm, &pStart,&pEnd);CHKERRQ(ierr);
28379bad331SMatthew G. Knepley     ierr  = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
28424c766afSToby Isaac     ierr  = DMPlexGetMaxSizes(dm, &maxC, &maxS);CHKERRQ(ierr);
28524c766afSToby Isaac     coneSeries    = (maxC > 1) ? ((PetscPowInt(maxC,depth+1)-1)/(maxC-1)) : depth+1;
28624c766afSToby Isaac     supportSeries = (maxS > 1) ? ((PetscPowInt(maxS,depth+1)-1)/(maxS-1)) : depth+1;
28724c766afSToby Isaac     asiz  = PetscMax(PetscPowInt(maxS,depth)*coneSeries,PetscPowInt(maxC,depth)*supportSeries);
2888b0b4c70SToby Isaac     asiz *= maxAnchors;
28924c766afSToby Isaac     asiz  = PetscMin(asiz,pEnd-pStart);
29079bad331SMatthew G. Knepley     ierr  = PetscMalloc1(asiz,adj);CHKERRQ(ierr);
29179bad331SMatthew G. Knepley   }
29279bad331SMatthew G. Knepley   if (*adjSize < 0) *adjSize = asiz;
2938b0b4c70SToby Isaac   maxAdjSize = *adjSize;
29470034214SMatthew G. Knepley   if (useTransitiveClosure) {
29579bad331SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Transitive_Internal(dm, p, useCone, adjSize, *adj);CHKERRQ(ierr);
29670034214SMatthew G. Knepley   } else if (useCone) {
29779bad331SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Cone_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
29870034214SMatthew G. Knepley   } else {
29979bad331SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Support_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
30070034214SMatthew G. Knepley   }
3015b317d89SToby Isaac   if (useAnchors && aSec) {
3028b0b4c70SToby Isaac     PetscInt origSize = *adjSize;
3038b0b4c70SToby Isaac     PetscInt numAdj = origSize;
3048b0b4c70SToby Isaac     PetscInt i = 0, j;
3058b0b4c70SToby Isaac     PetscInt *orig = *adj;
3068b0b4c70SToby Isaac 
3078b0b4c70SToby Isaac     while (i < origSize) {
3088b0b4c70SToby Isaac       PetscInt p = orig[i];
3098b0b4c70SToby Isaac       PetscInt aDof = 0;
3108b0b4c70SToby Isaac 
3118b0b4c70SToby Isaac       if (p >= aStart && p < aEnd) {
3128b0b4c70SToby Isaac         ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr);
3138b0b4c70SToby Isaac       }
3148b0b4c70SToby Isaac       if (aDof) {
3158b0b4c70SToby Isaac         PetscInt aOff;
3168b0b4c70SToby Isaac         PetscInt s, q;
3178b0b4c70SToby Isaac 
3188b0b4c70SToby Isaac         for (j = i + 1; j < numAdj; j++) {
3198b0b4c70SToby Isaac           orig[j - 1] = orig[j];
3208b0b4c70SToby Isaac         }
3218b0b4c70SToby Isaac         origSize--;
3228b0b4c70SToby Isaac         numAdj--;
3238b0b4c70SToby Isaac         ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr);
3248b0b4c70SToby Isaac         for (s = 0; s < aDof; ++s) {
3258b0b4c70SToby Isaac           for (q = 0; q < numAdj || (orig[numAdj++] = anchors[aOff+s],0); ++q) {
3268b0b4c70SToby Isaac             if (anchors[aOff+s] == orig[q]) break;
3278b0b4c70SToby Isaac           }
3288b0b4c70SToby Isaac           if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
3298b0b4c70SToby Isaac         }
3308b0b4c70SToby Isaac       }
3318b0b4c70SToby Isaac       else {
3328b0b4c70SToby Isaac         i++;
3338b0b4c70SToby Isaac       }
3348b0b4c70SToby Isaac     }
3358b0b4c70SToby Isaac     *adjSize = numAdj;
3368b0b4c70SToby Isaac     ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr);
3378b0b4c70SToby Isaac   }
33870034214SMatthew G. Knepley   PetscFunctionReturn(0);
33970034214SMatthew G. Knepley }
34070034214SMatthew G. Knepley 
34170034214SMatthew G. Knepley #undef __FUNCT__
34270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency"
34370034214SMatthew G. Knepley /*@
34470034214SMatthew G. Knepley   DMPlexGetAdjacency - Return all points adjacent to the given point
34570034214SMatthew G. Knepley 
34670034214SMatthew G. Knepley   Input Parameters:
34770034214SMatthew G. Knepley + dm - The DM object
34870034214SMatthew G. Knepley . p  - The point
34970034214SMatthew G. Knepley . adjSize - The maximum size of adj if it is non-NULL, or PETSC_DETERMINE
35070034214SMatthew G. Knepley - adj - Either NULL so that the array is allocated, or an existing array with size adjSize
35170034214SMatthew G. Knepley 
35270034214SMatthew G. Knepley   Output Parameters:
35370034214SMatthew G. Knepley + adjSize - The number of adjacent points
35470034214SMatthew G. Knepley - adj - The adjacent points
35570034214SMatthew G. Knepley 
35670034214SMatthew G. Knepley   Level: advanced
35770034214SMatthew G. Knepley 
35870034214SMatthew G. Knepley   Notes: The user must PetscFree the adj array if it was not passed in.
35970034214SMatthew G. Knepley 
36070034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexDistribute(), DMCreateMatrix(), DMPlexPreallocateOperator()
36170034214SMatthew G. Knepley @*/
36270034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacency(DM dm, PetscInt p, PetscInt *adjSize, PetscInt *adj[])
36370034214SMatthew G. Knepley {
36470034214SMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
36570034214SMatthew G. Knepley   PetscErrorCode ierr;
36670034214SMatthew G. Knepley 
36770034214SMatthew G. Knepley   PetscFunctionBeginHot;
36870034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
36970034214SMatthew G. Knepley   PetscValidPointer(adjSize,3);
37070034214SMatthew G. Knepley   PetscValidPointer(adj,4);
3715b317d89SToby Isaac   ierr = DMPlexGetAdjacency_Internal(dm, p, mesh->useCone, mesh->useClosure, mesh->useAnchors, adjSize, adj);CHKERRQ(ierr);
37270034214SMatthew G. Knepley   PetscFunctionReturn(0);
37370034214SMatthew G. Knepley }
374b0a623aaSMatthew G. Knepley #undef __FUNCT__
375b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexCreateTwoSidedProcessSF"
376b0a623aaSMatthew G. Knepley /*@
377b0a623aaSMatthew G. Knepley   DMPlexCreateTwoSidedProcessSF - Create an SF which just has process connectivity
378b0a623aaSMatthew G. Knepley 
379b0a623aaSMatthew G. Knepley   Collective on DM
380b0a623aaSMatthew G. Knepley 
381b0a623aaSMatthew G. Knepley   Input Parameters:
382b0a623aaSMatthew G. Knepley + dm      - The DM
383b0a623aaSMatthew G. Knepley - sfPoint - The PetscSF which encodes point connectivity
384b0a623aaSMatthew G. Knepley 
385b0a623aaSMatthew G. Knepley   Output Parameters:
386b0a623aaSMatthew G. Knepley + processRanks - A list of process neighbors, or NULL
387b0a623aaSMatthew G. Knepley - sfProcess    - An SF encoding the two-sided process connectivity, or NULL
388b0a623aaSMatthew G. Knepley 
389b0a623aaSMatthew G. Knepley   Level: developer
390b0a623aaSMatthew G. Knepley 
391b0a623aaSMatthew G. Knepley .seealso: PetscSFCreate(), DMPlexCreateProcessSF()
392b0a623aaSMatthew G. Knepley @*/
393b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM dm, PetscSF sfPoint, PetscSection rootRankSection, IS rootRanks, PetscSection leafRankSection, IS leafRanks, IS *processRanks, PetscSF *sfProcess)
394b0a623aaSMatthew G. Knepley {
395b0a623aaSMatthew G. Knepley   const PetscSFNode *remotePoints;
396b0a623aaSMatthew G. Knepley   PetscInt          *localPointsNew;
397b0a623aaSMatthew G. Knepley   PetscSFNode       *remotePointsNew;
398b0a623aaSMatthew G. Knepley   const PetscInt    *nranks;
399b0a623aaSMatthew G. Knepley   PetscInt          *ranksNew;
400b0a623aaSMatthew G. Knepley   PetscBT            neighbors;
401b0a623aaSMatthew G. Knepley   PetscInt           pStart, pEnd, p, numLeaves, l, numNeighbors, n;
402b0a623aaSMatthew G. Knepley   PetscMPIInt        numProcs, proc, rank;
403b0a623aaSMatthew G. Knepley   PetscErrorCode     ierr;
404b0a623aaSMatthew G. Knepley 
405b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
406b0a623aaSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
407b0a623aaSMatthew G. Knepley   PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2);
408b0a623aaSMatthew G. Knepley   if (processRanks) {PetscValidPointer(processRanks, 3);}
409b0a623aaSMatthew G. Knepley   if (sfProcess)    {PetscValidPointer(sfProcess, 4);}
410b0a623aaSMatthew G. Knepley   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr);
411b0a623aaSMatthew G. Knepley   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
412b0a623aaSMatthew G. Knepley   ierr = PetscSFGetGraph(sfPoint, NULL, &numLeaves, NULL, &remotePoints);CHKERRQ(ierr);
413b0a623aaSMatthew G. Knepley   ierr = PetscBTCreate(numProcs, &neighbors);CHKERRQ(ierr);
414b0a623aaSMatthew G. Knepley   ierr = PetscBTMemzero(numProcs, neighbors);CHKERRQ(ierr);
415b0a623aaSMatthew G. Knepley   /* Compute root-to-leaf process connectivity */
416b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetChart(rootRankSection, &pStart, &pEnd);CHKERRQ(ierr);
417b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(rootRanks, &nranks);CHKERRQ(ierr);
418b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
419b0a623aaSMatthew G. Knepley     PetscInt ndof, noff, n;
420b0a623aaSMatthew G. Knepley 
421b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetDof(rootRankSection, p, &ndof);CHKERRQ(ierr);
422b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(rootRankSection, p, &noff);CHKERRQ(ierr);
423302440fdSBarry Smith     for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);CHKERRQ(ierr);}
424b0a623aaSMatthew G. Knepley   }
425b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(rootRanks, &nranks);CHKERRQ(ierr);
426b0a623aaSMatthew G. Knepley   /* Compute leaf-to-neighbor process connectivity */
427b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetChart(leafRankSection, &pStart, &pEnd);CHKERRQ(ierr);
428b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(leafRanks, &nranks);CHKERRQ(ierr);
429b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
430b0a623aaSMatthew G. Knepley     PetscInt ndof, noff, n;
431b0a623aaSMatthew G. Knepley 
432b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetDof(leafRankSection, p, &ndof);CHKERRQ(ierr);
433b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(leafRankSection, p, &noff);CHKERRQ(ierr);
434302440fdSBarry Smith     for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);CHKERRQ(ierr);}
435b0a623aaSMatthew G. Knepley   }
436b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(leafRanks, &nranks);CHKERRQ(ierr);
437b0a623aaSMatthew G. Knepley   /* Compute leaf-to-root process connectivity */
438b0a623aaSMatthew G. Knepley   for (l = 0; l < numLeaves; ++l) {PetscBTSet(neighbors, remotePoints[l].rank);}
439b0a623aaSMatthew G. Knepley   /* Calculate edges */
440b0a623aaSMatthew G. Knepley   PetscBTClear(neighbors, rank);
441b0a623aaSMatthew G. Knepley   for(proc = 0, numNeighbors = 0; proc < numProcs; ++proc) {if (PetscBTLookup(neighbors, proc)) ++numNeighbors;}
442b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(numNeighbors, &ranksNew);CHKERRQ(ierr);
443b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(numNeighbors, &localPointsNew);CHKERRQ(ierr);
444b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(numNeighbors, &remotePointsNew);CHKERRQ(ierr);
445b0a623aaSMatthew G. Knepley   for(proc = 0, n = 0; proc < numProcs; ++proc) {
446b0a623aaSMatthew G. Knepley     if (PetscBTLookup(neighbors, proc)) {
447b0a623aaSMatthew G. Knepley       ranksNew[n]              = proc;
448b0a623aaSMatthew G. Knepley       localPointsNew[n]        = proc;
449b0a623aaSMatthew G. Knepley       remotePointsNew[n].index = rank;
450b0a623aaSMatthew G. Knepley       remotePointsNew[n].rank  = proc;
451b0a623aaSMatthew G. Knepley       ++n;
452b0a623aaSMatthew G. Knepley     }
453b0a623aaSMatthew G. Knepley   }
454b0a623aaSMatthew G. Knepley   ierr = PetscBTDestroy(&neighbors);CHKERRQ(ierr);
455b0a623aaSMatthew G. Knepley   if (processRanks) {ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), numNeighbors, ranksNew, PETSC_OWN_POINTER, processRanks);CHKERRQ(ierr);}
456b0a623aaSMatthew G. Knepley   else              {ierr = PetscFree(ranksNew);CHKERRQ(ierr);}
457b0a623aaSMatthew G. Knepley   if (sfProcess) {
458b0a623aaSMatthew G. Knepley     ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess);CHKERRQ(ierr);
459b0a623aaSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) *sfProcess, "Two-Sided Process SF");CHKERRQ(ierr);
460b0a623aaSMatthew G. Knepley     ierr = PetscSFSetFromOptions(*sfProcess);CHKERRQ(ierr);
461b0a623aaSMatthew G. Knepley     ierr = PetscSFSetGraph(*sfProcess, numProcs, numNeighbors, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr);
462b0a623aaSMatthew G. Knepley   }
463b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
464b0a623aaSMatthew G. Knepley }
465b0a623aaSMatthew G. Knepley 
466b0a623aaSMatthew G. Knepley #undef __FUNCT__
467b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeOwnership"
468b0a623aaSMatthew G. Knepley /*@
469b0a623aaSMatthew G. Knepley   DMPlexDistributeOwnership - Compute owner information for shared points. This basically gets two-sided for an SF.
470b0a623aaSMatthew G. Knepley 
471b0a623aaSMatthew G. Knepley   Collective on DM
472b0a623aaSMatthew G. Knepley 
473b0a623aaSMatthew G. Knepley   Input Parameter:
474b0a623aaSMatthew G. Knepley . dm - The DM
475b0a623aaSMatthew G. Knepley 
476b0a623aaSMatthew G. Knepley   Output Parameters:
477b0a623aaSMatthew G. Knepley + rootSection - The number of leaves for a given root point
478b0a623aaSMatthew G. Knepley . rootrank    - The rank of each edge into the root point
479b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point
480b0a623aaSMatthew G. Knepley - leafrank    - The rank of each process sharing a leaf point
481b0a623aaSMatthew G. Knepley 
482b0a623aaSMatthew G. Knepley   Level: developer
483b0a623aaSMatthew G. Knepley 
484b0a623aaSMatthew G. Knepley .seealso: DMPlexCreateOverlap()
485b0a623aaSMatthew G. Knepley @*/
486b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexDistributeOwnership(DM dm, PetscSection rootSection, IS *rootrank, PetscSection leafSection, IS *leafrank)
487b0a623aaSMatthew G. Knepley {
488b0a623aaSMatthew G. Knepley   MPI_Comm        comm;
489b0a623aaSMatthew G. Knepley   PetscSF         sfPoint;
490b0a623aaSMatthew G. Knepley   const PetscInt *rootdegree;
491b0a623aaSMatthew G. Knepley   PetscInt       *myrank, *remoterank;
492b0a623aaSMatthew G. Knepley   PetscInt        pStart, pEnd, p, nedges;
493b0a623aaSMatthew G. Knepley   PetscMPIInt     rank;
494b0a623aaSMatthew G. Knepley   PetscErrorCode  ierr;
495b0a623aaSMatthew G. Knepley 
496b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
497b0a623aaSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
498b0a623aaSMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
499b0a623aaSMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
500b0a623aaSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
501b0a623aaSMatthew G. Knepley   /* Compute number of leaves for each root */
502b0a623aaSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) rootSection, "Root Section");CHKERRQ(ierr);
503b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetChart(rootSection, pStart, pEnd);CHKERRQ(ierr);
504b0a623aaSMatthew G. Knepley   ierr = PetscSFComputeDegreeBegin(sfPoint, &rootdegree);CHKERRQ(ierr);
505b0a623aaSMatthew G. Knepley   ierr = PetscSFComputeDegreeEnd(sfPoint, &rootdegree);CHKERRQ(ierr);
506b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {ierr = PetscSectionSetDof(rootSection, p, rootdegree[p-pStart]);CHKERRQ(ierr);}
507b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
508b0a623aaSMatthew G. Knepley   /* Gather rank of each leaf to root */
509b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(rootSection, &nedges);CHKERRQ(ierr);
510b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(pEnd-pStart, &myrank);CHKERRQ(ierr);
511b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(nedges,  &remoterank);CHKERRQ(ierr);
512b0a623aaSMatthew G. Knepley   for (p = 0; p < pEnd-pStart; ++p) myrank[p] = rank;
513b0a623aaSMatthew G. Knepley   ierr = PetscSFGatherBegin(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr);
514b0a623aaSMatthew G. Knepley   ierr = PetscSFGatherEnd(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr);
515b0a623aaSMatthew G. Knepley   ierr = PetscFree(myrank);CHKERRQ(ierr);
516b0a623aaSMatthew G. Knepley   ierr = ISCreateGeneral(comm, nedges, remoterank, PETSC_OWN_POINTER, rootrank);CHKERRQ(ierr);
517b0a623aaSMatthew G. Knepley   /* Distribute remote ranks to leaves */
518b0a623aaSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) leafSection, "Leaf Section");CHKERRQ(ierr);
519b0a623aaSMatthew G. Knepley   ierr = DMPlexDistributeFieldIS(dm, sfPoint, rootSection, *rootrank, leafSection, leafrank);CHKERRQ(ierr);
520b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
521b0a623aaSMatthew G. Knepley }
522b0a623aaSMatthew G. Knepley 
523b0a623aaSMatthew G. Knepley #undef __FUNCT__
524b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexCreateOverlap"
525278397a0SMatthew G. Knepley /*@C
526b0a623aaSMatthew G. Knepley   DMPlexCreateOverlap - Compute owner information for shared points. This basically gets two-sided for an SF.
527b0a623aaSMatthew G. Knepley 
528b0a623aaSMatthew G. Knepley   Collective on DM
529b0a623aaSMatthew G. Knepley 
530b0a623aaSMatthew G. Knepley   Input Parameters:
531b0a623aaSMatthew G. Knepley + dm          - The DM
53224d039d7SMichael Lange . levels      - Number of overlap levels
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:
539a9f1d5b2SMichael Lange + ovLabel     - DMLabel containing remote overlap contributions as point/rank pairings
540b0a623aaSMatthew G. Knepley 
541b0a623aaSMatthew G. Knepley   Level: developer
542b0a623aaSMatthew G. Knepley 
5431fd9873aSMichael Lange .seealso: DMPlexDistributeOwnership(), DMPlexDistribute()
544b0a623aaSMatthew G. Knepley @*/
545a9f1d5b2SMichael Lange PetscErrorCode DMPlexCreateOverlap(DM dm, PetscInt levels, PetscSection rootSection, IS rootrank, PetscSection leafSection, IS leafrank, DMLabel *ovLabel)
546b0a623aaSMatthew G. Knepley {
547e540f424SMichael Lange   MPI_Comm           comm;
548b0a623aaSMatthew G. Knepley   DMLabel            ovAdjByRank; /* A DMLabel containing all points adjacent to shared points, separated by rank (value in label) */
549b0a623aaSMatthew G. Knepley   PetscSF            sfPoint, sfProc;
550b0a623aaSMatthew G. Knepley   const PetscSFNode *remote;
551b0a623aaSMatthew G. Knepley   const PetscInt    *local;
5521fd9873aSMichael Lange   const PetscInt    *nrank, *rrank;
553b0a623aaSMatthew G. Knepley   PetscInt          *adj = NULL;
5541fd9873aSMichael Lange   PetscInt           pStart, pEnd, p, sStart, sEnd, nleaves, l;
555b0a623aaSMatthew G. Knepley   PetscMPIInt        rank, numProcs;
55626a7d390SMatthew G. Knepley   PetscBool          useCone, useClosure, flg;
557b0a623aaSMatthew G. Knepley   PetscErrorCode     ierr;
558b0a623aaSMatthew G. Knepley 
559b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
560e540f424SMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
561e540f424SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
562e540f424SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
563b0a623aaSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
564b0a623aaSMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
565b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetChart(leafSection, &sStart, &sEnd);CHKERRQ(ierr);
566b0a623aaSMatthew G. Knepley   ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr);
567b0a623aaSMatthew G. Knepley   ierr = DMLabelCreate("Overlap adjacency", &ovAdjByRank);CHKERRQ(ierr);
568b0a623aaSMatthew G. Knepley   /* Handle leaves: shared with the root point */
569b0a623aaSMatthew G. Knepley   for (l = 0; l < nleaves; ++l) {
570b0a623aaSMatthew G. Knepley     PetscInt adjSize = PETSC_DETERMINE, a;
571b0a623aaSMatthew G. Knepley 
572b0a623aaSMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, local[l], &adjSize, &adj);CHKERRQ(ierr);
573b0a623aaSMatthew G. Knepley     for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remote[l].rank);CHKERRQ(ierr);}
574b0a623aaSMatthew G. Knepley   }
575b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(rootrank, &rrank);CHKERRQ(ierr);
576b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(leafrank, &nrank);CHKERRQ(ierr);
577b0a623aaSMatthew G. Knepley   /* Handle roots */
578b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
579b0a623aaSMatthew G. Knepley     PetscInt adjSize = PETSC_DETERMINE, neighbors = 0, noff, n, a;
580b0a623aaSMatthew G. Knepley 
581b0a623aaSMatthew G. Knepley     if ((p >= sStart) && (p < sEnd)) {
582b0a623aaSMatthew G. Knepley       /* Some leaves share a root with other leaves on different processes */
583b0a623aaSMatthew G. Knepley       ierr = PetscSectionGetDof(leafSection, p, &neighbors);CHKERRQ(ierr);
584b0a623aaSMatthew G. Knepley       if (neighbors) {
585b0a623aaSMatthew G. Knepley         ierr = PetscSectionGetOffset(leafSection, p, &noff);CHKERRQ(ierr);
586b0a623aaSMatthew G. Knepley         ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
587b0a623aaSMatthew G. Knepley         for (n = 0; n < neighbors; ++n) {
588b0a623aaSMatthew G. Knepley           const PetscInt remoteRank = nrank[noff+n];
589b0a623aaSMatthew G. Knepley 
590b0a623aaSMatthew G. Knepley           if (remoteRank == rank) continue;
591b0a623aaSMatthew G. Knepley           for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);}
592b0a623aaSMatthew G. Knepley         }
593b0a623aaSMatthew G. Knepley       }
594b0a623aaSMatthew G. Knepley     }
595b0a623aaSMatthew G. Knepley     /* Roots are shared with leaves */
596b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetDof(rootSection, p, &neighbors);CHKERRQ(ierr);
597b0a623aaSMatthew G. Knepley     if (!neighbors) continue;
598b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(rootSection, p, &noff);CHKERRQ(ierr);
599b0a623aaSMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
600b0a623aaSMatthew G. Knepley     for (n = 0; n < neighbors; ++n) {
601b0a623aaSMatthew G. Knepley       const PetscInt remoteRank = rrank[noff+n];
602b0a623aaSMatthew G. Knepley 
603b0a623aaSMatthew G. Knepley       if (remoteRank == rank) continue;
604b0a623aaSMatthew G. Knepley       for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);}
605b0a623aaSMatthew G. Knepley     }
606b0a623aaSMatthew G. Knepley   }
607b0a623aaSMatthew G. Knepley   ierr = PetscFree(adj);CHKERRQ(ierr);
608b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(rootrank, &rrank);CHKERRQ(ierr);
609b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(leafrank, &nrank);CHKERRQ(ierr);
61024d039d7SMichael Lange   /* Add additional overlap levels */
61124d039d7SMichael Lange   for (l = 1; l < levels; l++) {ierr = DMPlexPartitionLabelAdjacency(dm, ovAdjByRank);CHKERRQ(ierr);}
61226a7d390SMatthew G. Knepley   /* We require the closure in the overlap */
61326a7d390SMatthew G. Knepley   ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
61426a7d390SMatthew G. Knepley   ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr);
61526a7d390SMatthew G. Knepley   if (useCone || !useClosure) {
6165abbe4feSMichael Lange     ierr = DMPlexPartitionLabelClosure(dm, ovAdjByRank);CHKERRQ(ierr);
61726a7d390SMatthew G. Knepley   }
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   }
6231fd9873aSMichael Lange   /* Make process SF and invert sender to receiver label */
624b0a623aaSMatthew G. Knepley   ierr = DMPlexCreateTwoSidedProcessSF(dm, sfPoint, rootSection, rootrank, leafSection, leafrank, NULL, &sfProc);CHKERRQ(ierr);
625a9f1d5b2SMichael Lange   ierr = DMLabelCreate("Overlap label", ovLabel);CHKERRQ(ierr);
626a9f1d5b2SMichael Lange   ierr = DMPlexPartitionLabelInvert(dm, ovAdjByRank, sfProc, *ovLabel);CHKERRQ(ierr);
627a9f1d5b2SMichael Lange   /* Add owned points, except for shared local points */
628a9f1d5b2SMichael Lange   for (p = pStart; p < pEnd; ++p) {ierr = DMLabelSetValue(*ovLabel, p, rank);CHKERRQ(ierr);}
629e540f424SMichael Lange   for (l = 0; l < nleaves; ++l) {
630a9f1d5b2SMichael Lange     ierr = DMLabelClearValue(*ovLabel, local[l], rank);CHKERRQ(ierr);
631a9f1d5b2SMichael Lange     ierr = DMLabelSetValue(*ovLabel, remote[l].index, remote[l].rank);CHKERRQ(ierr);
632e540f424SMichael Lange   }
633e540f424SMichael Lange   /* Clean up */
6341fd9873aSMichael Lange   ierr = DMLabelDestroy(&ovAdjByRank);CHKERRQ(ierr);
635b0a623aaSMatthew G. Knepley   ierr = PetscSFDestroy(&sfProc);CHKERRQ(ierr);
636b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
637b0a623aaSMatthew G. Knepley }
63870034214SMatthew G. Knepley 
63970034214SMatthew G. Knepley #undef __FUNCT__
64046f9b1c3SMichael Lange #define __FUNCT__ "DMPlexCreateOverlapMigrationSF"
64146f9b1c3SMichael Lange PetscErrorCode DMPlexCreateOverlapMigrationSF(DM dm, PetscSF overlapSF, PetscSF *migrationSF)
64246f9b1c3SMichael Lange {
64346f9b1c3SMichael Lange   MPI_Comm           comm;
64446f9b1c3SMichael Lange   PetscMPIInt        rank, numProcs;
64546f9b1c3SMichael Lange   PetscInt           d, dim, p, pStart, pEnd, nroots, nleaves, newLeaves, point, numSharedPoints;
64646f9b1c3SMichael Lange   PetscInt          *pointDepths, *remoteDepths, *ilocal;
64746f9b1c3SMichael Lange   PetscInt          *depthRecv, *depthShift, *depthIdx;
64846f9b1c3SMichael Lange   PetscSFNode       *iremote;
64946f9b1c3SMichael Lange   PetscSF            pointSF;
65046f9b1c3SMichael Lange   const PetscInt    *sharedLocal;
65146f9b1c3SMichael Lange   const PetscSFNode *overlapRemote, *sharedRemote;
65246f9b1c3SMichael Lange   PetscErrorCode     ierr;
65346f9b1c3SMichael Lange 
65446f9b1c3SMichael Lange   PetscFunctionBegin;
65546f9b1c3SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
65646f9b1c3SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
65746f9b1c3SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
65846f9b1c3SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
65946f9b1c3SMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
66046f9b1c3SMichael Lange 
66146f9b1c3SMichael Lange   /* Before building the migration SF we need to know the new stratum offsets */
66246f9b1c3SMichael Lange   ierr = PetscSFGetGraph(overlapSF, &nroots, &nleaves, NULL, &overlapRemote);CHKERRQ(ierr);
66346f9b1c3SMichael Lange   ierr = PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);CHKERRQ(ierr);
66446f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
66546f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
66646f9b1c3SMichael Lange     for (p=pStart; p<pEnd; p++) pointDepths[p] = d;
66746f9b1c3SMichael Lange   }
66846f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) remoteDepths[p] = -1;
66946f9b1c3SMichael Lange   ierr = PetscSFBcastBegin(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
67046f9b1c3SMichael Lange   ierr = PetscSFBcastEnd(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
67146f9b1c3SMichael Lange 
67246f9b1c3SMichael Lange   /* Count recevied points in each stratum and compute the internal strata shift */
67346f9b1c3SMichael Lange   ierr = PetscMalloc3(dim+1, &depthRecv, dim+1, &depthShift, dim+1, &depthIdx);CHKERRQ(ierr);
67446f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) depthRecv[d]=0;
67546f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) depthRecv[remoteDepths[p]]++;
67646f9b1c3SMichael Lange   depthShift[dim] = 0;
67746f9b1c3SMichael Lange   for (d=0; d<dim; d++) depthShift[d] = depthRecv[dim];
67846f9b1c3SMichael Lange   for (d=1; d<dim; d++) depthShift[d] += depthRecv[0];
67946f9b1c3SMichael Lange   for (d=dim-2; d>0; d--) depthShift[d] += depthRecv[d+1];
68046f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
68146f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
68246f9b1c3SMichael Lange     depthIdx[d] = pStart + depthShift[d];
68346f9b1c3SMichael Lange   }
68446f9b1c3SMichael Lange 
68546f9b1c3SMichael Lange   /* Form the overlap SF build an SF that describes the full overlap migration SF */
68646f9b1c3SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
68746f9b1c3SMichael Lange   newLeaves = pEnd - pStart + nleaves;
68809b7985cSMatthew G. Knepley   ierr = PetscMalloc1(newLeaves, &ilocal);CHKERRQ(ierr);
68909b7985cSMatthew G. Knepley   ierr = PetscMalloc1(newLeaves, &iremote);CHKERRQ(ierr);
69046f9b1c3SMichael Lange   /* First map local points to themselves */
69146f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
69246f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
69346f9b1c3SMichael Lange     for (p=pStart; p<pEnd; p++) {
69446f9b1c3SMichael Lange       point = p + depthShift[d];
69546f9b1c3SMichael Lange       ilocal[point] = point;
69646f9b1c3SMichael Lange       iremote[point].index = p;
69746f9b1c3SMichael Lange       iremote[point].rank = rank;
69846f9b1c3SMichael Lange       depthIdx[d]++;
69946f9b1c3SMichael Lange     }
70046f9b1c3SMichael Lange   }
70146f9b1c3SMichael Lange 
70246f9b1c3SMichael Lange   /* Add in the remote roots for currently shared points */
70346f9b1c3SMichael Lange   ierr = DMGetPointSF(dm, &pointSF);CHKERRQ(ierr);
70446f9b1c3SMichael Lange   ierr = PetscSFGetGraph(pointSF, NULL, &numSharedPoints, &sharedLocal, &sharedRemote);CHKERRQ(ierr);
70546f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
70646f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
70746f9b1c3SMichael Lange     for (p=0; p<numSharedPoints; p++) {
70846f9b1c3SMichael Lange       if (pStart <= sharedLocal[p] && sharedLocal[p] < pEnd) {
70946f9b1c3SMichael Lange         point = sharedLocal[p] + depthShift[d];
71046f9b1c3SMichael Lange         iremote[point].index = sharedRemote[p].index;
71146f9b1c3SMichael Lange         iremote[point].rank = sharedRemote[p].rank;
71246f9b1c3SMichael Lange       }
71346f9b1c3SMichael Lange     }
71446f9b1c3SMichael Lange   }
71546f9b1c3SMichael Lange 
71646f9b1c3SMichael Lange   /* Now add the incoming overlap points */
71746f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) {
71846f9b1c3SMichael Lange     point = depthIdx[remoteDepths[p]];
71946f9b1c3SMichael Lange     ilocal[point] = point;
72046f9b1c3SMichael Lange     iremote[point].index = overlapRemote[p].index;
72146f9b1c3SMichael Lange     iremote[point].rank = overlapRemote[p].rank;
72246f9b1c3SMichael Lange     depthIdx[remoteDepths[p]]++;
72346f9b1c3SMichael Lange   }
72415fff7beSMatthew G. Knepley   ierr = PetscFree2(pointDepths,remoteDepths);CHKERRQ(ierr);
72546f9b1c3SMichael Lange 
72646f9b1c3SMichael Lange   ierr = PetscSFCreate(comm, migrationSF);CHKERRQ(ierr);
72746f9b1c3SMichael Lange   ierr = PetscObjectSetName((PetscObject) *migrationSF, "Overlap Migration SF");CHKERRQ(ierr);
72846f9b1c3SMichael Lange   ierr = PetscSFSetFromOptions(*migrationSF);CHKERRQ(ierr);
72946f9b1c3SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
73046f9b1c3SMichael Lange   ierr = PetscSFSetGraph(*migrationSF, pEnd-pStart, newLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
73146f9b1c3SMichael Lange 
73246f9b1c3SMichael Lange   ierr = PetscFree3(depthRecv, depthShift, depthIdx);CHKERRQ(ierr);
73370034214SMatthew G. Knepley   PetscFunctionReturn(0);
73470034214SMatthew G. Knepley }
73570034214SMatthew G. Knepley 
73670034214SMatthew G. Knepley #undef __FUNCT__
737a9f1d5b2SMichael Lange #define __FUNCT__ "DMPlexStratifyMigrationSF"
738a9f1d5b2SMichael Lange /*@
739a9f1d5b2SMichael Lange   DMPlexStratifyMigrationSF - Add partition overlap to a distributed non-overlapping DM.
740a9f1d5b2SMichael Lange 
741a9f1d5b2SMichael Lange   Input Parameter:
742a9f1d5b2SMichael Lange + dm          - The DM
743a9f1d5b2SMichael Lange - sf          - A star forest with non-ordered leaves, usually defining a DM point migration
744a9f1d5b2SMichael Lange 
745a9f1d5b2SMichael Lange   Output Parameter:
746a9f1d5b2SMichael Lange . migrationSF - A star forest with added leaf indirection that ensures the resulting DM is stratified
747a9f1d5b2SMichael Lange 
748a9f1d5b2SMichael Lange   Level: developer
749a9f1d5b2SMichael Lange 
750a9f1d5b2SMichael Lange .seealso: DMPlexPartitionLabelCreateSF(), DMPlexDistribute(), DMPlexDistributeOverlap()
751a9f1d5b2SMichael Lange @*/
752a9f1d5b2SMichael Lange PetscErrorCode DMPlexStratifyMigrationSF(DM dm, PetscSF sf, PetscSF *migrationSF)
753a9f1d5b2SMichael Lange {
754a9f1d5b2SMichael Lange   MPI_Comm           comm;
755a9f1d5b2SMichael Lange   PetscMPIInt        rank, numProcs;
7567fab53ddSMatthew G. Knepley   PetscInt           d, ldepth, depth, p, pStart, pEnd, nroots, nleaves;
757a9f1d5b2SMichael Lange   PetscInt          *pointDepths, *remoteDepths, *ilocal;
758a9f1d5b2SMichael Lange   PetscInt          *depthRecv, *depthShift, *depthIdx;
759a9f1d5b2SMichael Lange   const PetscSFNode *iremote;
760a9f1d5b2SMichael Lange   PetscErrorCode     ierr;
761a9f1d5b2SMichael Lange 
762a9f1d5b2SMichael Lange   PetscFunctionBegin;
763a9f1d5b2SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
764a9f1d5b2SMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
765a9f1d5b2SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
766a9f1d5b2SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
7677fab53ddSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &ldepth);CHKERRQ(ierr);
7687fab53ddSMatthew G. Knepley   ierr = MPI_Allreduce(&ldepth, &depth, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
7697fab53ddSMatthew G. Knepley   if ((ldepth >= 0) && (depth != ldepth)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent Plex depth %d != %d", ldepth, depth);
770a9f1d5b2SMichael Lange 
771a9f1d5b2SMichael Lange   /* Before building the migration SF we need to know the new stratum offsets */
772a9f1d5b2SMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, &iremote);CHKERRQ(ierr);
773a9f1d5b2SMichael Lange   ierr = PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);CHKERRQ(ierr);
7747fab53ddSMatthew G. Knepley   for (d = 0; d < depth+1; ++d) {
775a9f1d5b2SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
7767fab53ddSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) pointDepths[p] = d;
777a9f1d5b2SMichael Lange   }
7787fab53ddSMatthew G. Knepley   for (p = 0; p < nleaves; ++p) remoteDepths[p] = -1;
779a9f1d5b2SMichael Lange   ierr = PetscSFBcastBegin(sf, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
780a9f1d5b2SMichael Lange   ierr = PetscSFBcastEnd(sf, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
781a9f1d5b2SMichael Lange   /* Count recevied points in each stratum and compute the internal strata shift */
7827fab53ddSMatthew G. Knepley   ierr = PetscMalloc3(depth+1, &depthRecv, depth+1, &depthShift, depth+1, &depthIdx);CHKERRQ(ierr);
7837fab53ddSMatthew G. Knepley   for (d = 0; d < depth+1; ++d) depthRecv[d] = 0;
7847fab53ddSMatthew G. Knepley   for (p = 0; p < nleaves; ++p) depthRecv[remoteDepths[p]]++;
7857fab53ddSMatthew G. Knepley   depthShift[depth] = 0;
7867fab53ddSMatthew G. Knepley   for (d = 0; d < depth; ++d) depthShift[d] = depthRecv[depth];
7877fab53ddSMatthew G. Knepley   for (d = 1; d < depth; ++d) depthShift[d] += depthRecv[0];
7887fab53ddSMatthew G. Knepley   for (d = depth-2; d > 0; --d) depthShift[d] += depthRecv[d+1];
7897fab53ddSMatthew G. Knepley   for (d = 0; d < depth+1; ++d) {depthIdx[d] = 0;}
790a9f1d5b2SMichael Lange   /* Derive a new local permutation based on stratified indices */
791a9f1d5b2SMichael Lange   ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr);
7927fab53ddSMatthew G. Knepley   for (p = 0; p < nleaves; ++p) {
7937fab53ddSMatthew G. Knepley     const PetscInt dep = remoteDepths[p];
7947fab53ddSMatthew G. Knepley 
7957fab53ddSMatthew G. Knepley     ilocal[p] = depthShift[dep] + depthIdx[dep];
7967fab53ddSMatthew G. Knepley     depthIdx[dep]++;
797a9f1d5b2SMichael Lange   }
798a9f1d5b2SMichael Lange   ierr = PetscSFCreate(comm, migrationSF);CHKERRQ(ierr);
799a9f1d5b2SMichael Lange   ierr = PetscObjectSetName((PetscObject) *migrationSF, "Migration SF");CHKERRQ(ierr);
800a9f1d5b2SMichael Lange   ierr = PetscSFSetGraph(*migrationSF, nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_COPY_VALUES);CHKERRQ(ierr);
801a9f1d5b2SMichael Lange   ierr = PetscFree2(pointDepths,remoteDepths);CHKERRQ(ierr);
802a9f1d5b2SMichael Lange   ierr = PetscFree3(depthRecv, depthShift, depthIdx);CHKERRQ(ierr);
803a9f1d5b2SMichael Lange   PetscFunctionReturn(0);
804a9f1d5b2SMichael Lange }
805a9f1d5b2SMichael Lange 
806a9f1d5b2SMichael Lange #undef __FUNCT__
80770034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeField"
80870034214SMatthew G. Knepley /*@
80970034214SMatthew G. Knepley   DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
81070034214SMatthew G. Knepley 
81170034214SMatthew G. Knepley   Collective on DM
81270034214SMatthew G. Knepley 
81370034214SMatthew G. Knepley   Input Parameters:
81470034214SMatthew G. Knepley + dm - The DMPlex object
81570034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
81670034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
81770034214SMatthew G. Knepley - originalVec - The existing data
81870034214SMatthew G. Knepley 
81970034214SMatthew G. Knepley   Output Parameters:
82070034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout
82170034214SMatthew G. Knepley - newVec - The new data
82270034214SMatthew G. Knepley 
82370034214SMatthew G. Knepley   Level: developer
82470034214SMatthew G. Knepley 
82570034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeFieldIS(), DMPlexDistributeData()
82670034214SMatthew G. Knepley @*/
82770034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec)
82870034214SMatthew G. Knepley {
82970034214SMatthew G. Knepley   PetscSF        fieldSF;
83070034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
83170034214SMatthew G. Knepley   PetscScalar   *originalValues, *newValues;
83270034214SMatthew G. Knepley   PetscErrorCode ierr;
83370034214SMatthew G. Knepley 
83470034214SMatthew G. Knepley   PetscFunctionBegin;
83570034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
83670034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
83770034214SMatthew G. Knepley 
83870034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
83970034214SMatthew G. Knepley   ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr);
84070034214SMatthew G. Knepley   ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr);
84170034214SMatthew G. Knepley 
84270034214SMatthew G. Knepley   ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr);
84370034214SMatthew G. Knepley   ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr);
84470034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
845*8a713f3bSLawrence Mitchell   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
84670034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
84770034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
84870034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
84970034214SMatthew G. Knepley   ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr);
85070034214SMatthew G. Knepley   ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr);
85170034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
85270034214SMatthew G. Knepley   PetscFunctionReturn(0);
85370034214SMatthew G. Knepley }
85470034214SMatthew G. Knepley 
85570034214SMatthew G. Knepley #undef __FUNCT__
85670034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeFieldIS"
85770034214SMatthew G. Knepley /*@
85870034214SMatthew G. Knepley   DMPlexDistributeFieldIS - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
85970034214SMatthew G. Knepley 
86070034214SMatthew G. Knepley   Collective on DM
86170034214SMatthew G. Knepley 
86270034214SMatthew G. Knepley   Input Parameters:
86370034214SMatthew G. Knepley + dm - The DMPlex object
86470034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
86570034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
86670034214SMatthew G. Knepley - originalIS - The existing data
86770034214SMatthew G. Knepley 
86870034214SMatthew G. Knepley   Output Parameters:
86970034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout
87070034214SMatthew G. Knepley - newIS - The new data
87170034214SMatthew G. Knepley 
87270034214SMatthew G. Knepley   Level: developer
87370034214SMatthew G. Knepley 
87470034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexDistributeData()
87570034214SMatthew G. Knepley @*/
87670034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeFieldIS(DM dm, PetscSF pointSF, PetscSection originalSection, IS originalIS, PetscSection newSection, IS *newIS)
87770034214SMatthew G. Knepley {
87870034214SMatthew G. Knepley   PetscSF         fieldSF;
87970034214SMatthew G. Knepley   PetscInt       *newValues, *remoteOffsets, fieldSize;
88070034214SMatthew G. Knepley   const PetscInt *originalValues;
88170034214SMatthew G. Knepley   PetscErrorCode  ierr;
88270034214SMatthew G. Knepley 
88370034214SMatthew G. Knepley   PetscFunctionBegin;
88470034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
88570034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
88670034214SMatthew G. Knepley 
88770034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
888854ce69bSBarry Smith   ierr = PetscMalloc1(fieldSize, &newValues);CHKERRQ(ierr);
88970034214SMatthew G. Knepley 
89070034214SMatthew G. Knepley   ierr = ISGetIndices(originalIS, &originalValues);CHKERRQ(ierr);
89170034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
892*8a713f3bSLawrence Mitchell   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
893aaf8c182SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr);
894aaf8c182SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr);
89570034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
89670034214SMatthew G. Knepley   ierr = ISRestoreIndices(originalIS, &originalValues);CHKERRQ(ierr);
89770034214SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) pointSF), fieldSize, newValues, PETSC_OWN_POINTER, newIS);CHKERRQ(ierr);
89870034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
89970034214SMatthew G. Knepley   PetscFunctionReturn(0);
90070034214SMatthew G. Knepley }
90170034214SMatthew G. Knepley 
90270034214SMatthew G. Knepley #undef __FUNCT__
90370034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeData"
90470034214SMatthew G. Knepley /*@
90570034214SMatthew G. Knepley   DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
90670034214SMatthew G. Knepley 
90770034214SMatthew G. Knepley   Collective on DM
90870034214SMatthew G. Knepley 
90970034214SMatthew G. Knepley   Input Parameters:
91070034214SMatthew G. Knepley + dm - The DMPlex object
91170034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
91270034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
91370034214SMatthew G. Knepley . datatype - The type of data
91470034214SMatthew G. Knepley - originalData - The existing data
91570034214SMatthew G. Knepley 
91670034214SMatthew G. Knepley   Output Parameters:
91770034214SMatthew G. Knepley + newSection - The PetscSection describing the new data layout
91870034214SMatthew G. Knepley - newData - The new data
91970034214SMatthew G. Knepley 
92070034214SMatthew G. Knepley   Level: developer
92170034214SMatthew G. Knepley 
92270034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField()
92370034214SMatthew G. Knepley @*/
92470034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData)
92570034214SMatthew G. Knepley {
92670034214SMatthew G. Knepley   PetscSF        fieldSF;
92770034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
92870034214SMatthew G. Knepley   PetscMPIInt    dataSize;
92970034214SMatthew G. Knepley   PetscErrorCode ierr;
93070034214SMatthew G. Knepley 
93170034214SMatthew G. Knepley   PetscFunctionBegin;
93270034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
93370034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
93470034214SMatthew G. Knepley 
93570034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
93670034214SMatthew G. Knepley   ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr);
93770034214SMatthew G. Knepley   ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr);
93870034214SMatthew G. Knepley 
93970034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
940*8a713f3bSLawrence Mitchell   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
94170034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
94270034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
94370034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
94470034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
94570034214SMatthew G. Knepley   PetscFunctionReturn(0);
94670034214SMatthew G. Knepley }
94770034214SMatthew G. Knepley 
94870034214SMatthew G. Knepley #undef __FUNCT__
949cc71bff1SMichael Lange #define __FUNCT__ "DMPlexDistributeCones"
950cc71bff1SMichael Lange PetscErrorCode DMPlexDistributeCones(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
951cc71bff1SMichael Lange {
952f5bf2dbfSMichael Lange   DM_Plex               *mesh  = (DM_Plex*) dm->data;
953cc71bff1SMichael Lange   DM_Plex               *pmesh = (DM_Plex*) (dmParallel)->data;
954cc71bff1SMichael Lange   MPI_Comm               comm;
955cc71bff1SMichael Lange   PetscSF                coneSF;
956cc71bff1SMichael Lange   PetscSection           originalConeSection, newConeSection;
957cc71bff1SMichael Lange   PetscInt              *remoteOffsets, *cones, *newCones, newConesSize;
958cc71bff1SMichael Lange   PetscBool              flg;
959cc71bff1SMichael Lange   PetscErrorCode         ierr;
960cc71bff1SMichael Lange 
961cc71bff1SMichael Lange   PetscFunctionBegin;
962cc71bff1SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
963cc71bff1SMichael Lange   PetscValidPointer(dmParallel,4);
964cc71bff1SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
965cc71bff1SMichael Lange 
966cc71bff1SMichael Lange   /* Distribute cone section */
967cc71bff1SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
968cc71bff1SMichael Lange   ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr);
969cc71bff1SMichael Lange   ierr = DMPlexGetConeSection(dmParallel, &newConeSection);CHKERRQ(ierr);
970cc71bff1SMichael Lange   ierr = PetscSFDistributeSection(migrationSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr);
971cc71bff1SMichael Lange   ierr = DMSetUp(dmParallel);CHKERRQ(ierr);
972cc71bff1SMichael Lange   {
973cc71bff1SMichael Lange     PetscInt pStart, pEnd, p;
974cc71bff1SMichael Lange 
975cc71bff1SMichael Lange     ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr);
976cc71bff1SMichael Lange     for (p = pStart; p < pEnd; ++p) {
977cc71bff1SMichael Lange       PetscInt coneSize;
978cc71bff1SMichael Lange       ierr               = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr);
979cc71bff1SMichael Lange       pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize);
980cc71bff1SMichael Lange     }
981cc71bff1SMichael Lange   }
982cc71bff1SMichael Lange   /* Communicate and renumber cones */
983cc71bff1SMichael Lange   ierr = PetscSFCreateSectionSF(migrationSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr);
984*8a713f3bSLawrence Mitchell   ierr = PetscFree(remoteOffsets);CHKERRQ(ierr);
985cc71bff1SMichael Lange   ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
986cc71bff1SMichael Lange   ierr = DMPlexGetCones(dmParallel, &newCones);CHKERRQ(ierr);
987cc71bff1SMichael Lange   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
988cc71bff1SMichael Lange   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
989cc71bff1SMichael Lange   ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr);
990cc71bff1SMichael Lange   ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr);
9913533c52bSMatthew G. Knepley #if PETSC_USE_DEBUG
9923533c52bSMatthew G. Knepley   {
9933533c52bSMatthew G. Knepley     PetscInt  p;
9943533c52bSMatthew G. Knepley     PetscBool valid = PETSC_TRUE;
9953533c52bSMatthew G. Knepley     for (p = 0; p < newConesSize; ++p) {
9963533c52bSMatthew G. Knepley       if (newCones[p] < 0) {valid = PETSC_FALSE; ierr = PetscPrintf(PETSC_COMM_SELF, "Point %d not in overlap SF\n", p);CHKERRQ(ierr);}
9973533c52bSMatthew G. Knepley     }
9983533c52bSMatthew G. Knepley     if (!valid) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map");
9993533c52bSMatthew G. Knepley   }
10003533c52bSMatthew G. Knepley #endif
1001cc71bff1SMichael Lange   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr);
1002cc71bff1SMichael Lange   if (flg) {
1003cc71bff1SMichael Lange     ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr);
1004cc71bff1SMichael Lange     ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1005cc71bff1SMichael Lange     ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr);
1006cc71bff1SMichael Lange     ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1007cc71bff1SMichael Lange     ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr);
1008cc71bff1SMichael Lange   }
1009cc71bff1SMichael Lange   ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr);
1010cc71bff1SMichael Lange   ierr = DMPlexGetConeOrientations(dmParallel, &newCones);CHKERRQ(ierr);
1011cc71bff1SMichael Lange   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
1012cc71bff1SMichael Lange   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
1013cc71bff1SMichael Lange   ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr);
1014cc71bff1SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
1015cc71bff1SMichael Lange   /* Create supports and stratify sieve */
1016cc71bff1SMichael Lange   {
1017cc71bff1SMichael Lange     PetscInt pStart, pEnd;
1018cc71bff1SMichael Lange 
1019cc71bff1SMichael Lange     ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr);
1020cc71bff1SMichael Lange     ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr);
1021cc71bff1SMichael Lange   }
1022cc71bff1SMichael Lange   ierr = DMPlexSymmetrize(dmParallel);CHKERRQ(ierr);
1023cc71bff1SMichael Lange   ierr = DMPlexStratify(dmParallel);CHKERRQ(ierr);
1024f5bf2dbfSMichael Lange   pmesh->useCone    = mesh->useCone;
1025f5bf2dbfSMichael Lange   pmesh->useClosure = mesh->useClosure;
1026cc71bff1SMichael Lange   PetscFunctionReturn(0);
1027cc71bff1SMichael Lange }
1028cc71bff1SMichael Lange 
1029cc71bff1SMichael Lange #undef __FUNCT__
10300df0e737SMichael Lange #define __FUNCT__ "DMPlexDistributeCoordinates"
10310df0e737SMichael Lange PetscErrorCode DMPlexDistributeCoordinates(DM dm, PetscSF migrationSF, DM dmParallel)
10320df0e737SMichael Lange {
10330df0e737SMichael Lange   MPI_Comm         comm;
10340df0e737SMichael Lange   PetscSection     originalCoordSection, newCoordSection;
10350df0e737SMichael Lange   Vec              originalCoordinates, newCoordinates;
10360df0e737SMichael Lange   PetscInt         bs;
10370df0e737SMichael Lange   const char      *name;
10380df0e737SMichael Lange   const PetscReal *maxCell, *L;
10395dc8c3f7SMatthew G. Knepley   const DMBoundaryType *bd;
10400df0e737SMichael Lange   PetscErrorCode   ierr;
10410df0e737SMichael Lange 
10420df0e737SMichael Lange   PetscFunctionBegin;
10430df0e737SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10440df0e737SMichael Lange   PetscValidPointer(dmParallel, 3);
10450df0e737SMichael Lange 
10460df0e737SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
10470df0e737SMichael Lange   ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr);
10480df0e737SMichael Lange   ierr = DMGetCoordinateSection(dmParallel, &newCoordSection);CHKERRQ(ierr);
10490df0e737SMichael Lange   ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr);
10500df0e737SMichael Lange   if (originalCoordinates) {
10510df0e737SMichael Lange     ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr);
10520df0e737SMichael Lange     ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr);
10530df0e737SMichael Lange     ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr);
10540df0e737SMichael Lange 
10550df0e737SMichael Lange     ierr = DMPlexDistributeField(dm, migrationSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr);
10560df0e737SMichael Lange     ierr = DMSetCoordinatesLocal(dmParallel, newCoordinates);CHKERRQ(ierr);
10570df0e737SMichael Lange     ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr);
10580df0e737SMichael Lange     ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr);
10590df0e737SMichael Lange     ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr);
10600df0e737SMichael Lange   }
10615dc8c3f7SMatthew G. Knepley   ierr = DMGetPeriodicity(dm, &maxCell, &L, &bd);CHKERRQ(ierr);
10625dc8c3f7SMatthew G. Knepley   if (L) {ierr = DMSetPeriodicity(dmParallel, maxCell, L, bd);CHKERRQ(ierr);}
10630df0e737SMichael Lange   PetscFunctionReturn(0);
10640df0e737SMichael Lange }
10650df0e737SMichael Lange 
10660df0e737SMichael Lange #undef __FUNCT__
10670df0e737SMichael Lange #define __FUNCT__ "DMPlexDistributeLabels"
1068d995df53SMatthew G. Knepley /* Here we are assuming that process 0 always has everything */
10692eb1fa14SMichael Lange PetscErrorCode DMPlexDistributeLabels(DM dm, PetscSF migrationSF, DM dmParallel)
10700df0e737SMichael Lange {
10710df0e737SMichael Lange   MPI_Comm       comm;
10720df0e737SMichael Lange   PetscMPIInt    rank;
1073d995df53SMatthew G. Knepley   PetscInt       numLabels, numLocalLabels, l;
1074d995df53SMatthew G. Knepley   PetscBool      hasLabels = PETSC_FALSE;
10750df0e737SMichael Lange   PetscErrorCode ierr;
10760df0e737SMichael Lange 
10770df0e737SMichael Lange   PetscFunctionBegin;
10780df0e737SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10792eb1fa14SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
10800df0e737SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
10810df0e737SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
10820df0e737SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
10830df0e737SMichael Lange 
1084d995df53SMatthew G. Knepley   /* Everyone must have either the same number of labels, or none */
1085d995df53SMatthew G. Knepley   ierr = DMPlexGetNumLabels(dm, &numLocalLabels);CHKERRQ(ierr);
1086d995df53SMatthew G. Knepley   numLabels = numLocalLabels;
10870df0e737SMichael Lange   ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1088627847f0SMatthew G. Knepley   if (numLabels == numLocalLabels) hasLabels = PETSC_TRUE;
1089bdd2d751SMichael Lange   for (l = numLabels-1; l >= 0; --l) {
1090bdd2d751SMichael Lange     DMLabel     label = NULL, labelNew = NULL;
10910df0e737SMichael Lange     PetscBool   isdepth;
10920df0e737SMichael Lange 
1093d995df53SMatthew G. Knepley     if (hasLabels) {
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 {
115415078cd4SMichael Lange   DM_Plex        *mesh  = (DM_Plex*) dm->data;
115515078cd4SMichael Lange   DM_Plex        *pmesh = (DM_Plex*) (dmParallel)->data;
1156a6f36705SMichael Lange   MPI_Comm        comm;
1157a6f36705SMichael Lange   DM              refTree;
1158a6f36705SMichael Lange   PetscSection    origParentSection, newParentSection;
1159a6f36705SMichael Lange   PetscInt        *origParents, *origChildIDs;
1160a6f36705SMichael Lange   PetscBool       flg;
1161a6f36705SMichael Lange   PetscErrorCode  ierr;
1162a6f36705SMichael Lange 
1163a6f36705SMichael Lange   PetscFunctionBegin;
1164a6f36705SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1165a6f36705SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 4);
1166a6f36705SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
1167a6f36705SMichael Lange 
1168a6f36705SMichael Lange   /* Set up tree */
1169a6f36705SMichael Lange   ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr);
1170a6f36705SMichael Lange   ierr = DMPlexSetReferenceTree(dmParallel,refTree);CHKERRQ(ierr);
1171a6f36705SMichael Lange   ierr = DMPlexGetTree(dm,&origParentSection,&origParents,&origChildIDs,NULL,NULL);CHKERRQ(ierr);
1172a6f36705SMichael Lange   if (origParentSection) {
1173a6f36705SMichael Lange     PetscInt        pStart, pEnd;
1174a6f36705SMichael Lange     PetscInt        *newParents, *newChildIDs;
1175a6f36705SMichael Lange     PetscInt        *remoteOffsetsParents, newParentSize;
1176a6f36705SMichael Lange     PetscSF         parentSF;
1177a6f36705SMichael Lange 
1178a6f36705SMichael Lange     ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
1179a6f36705SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dmParallel),&newParentSection);CHKERRQ(ierr);
1180a6f36705SMichael Lange     ierr = PetscSectionSetChart(newParentSection,pStart,pEnd);CHKERRQ(ierr);
1181a6f36705SMichael Lange     ierr = PetscSFDistributeSection(migrationSF, origParentSection, &remoteOffsetsParents, newParentSection);CHKERRQ(ierr);
1182a6f36705SMichael Lange     ierr = PetscSFCreateSectionSF(migrationSF, origParentSection, remoteOffsetsParents, newParentSection, &parentSF);CHKERRQ(ierr);
1183*8a713f3bSLawrence Mitchell     ierr = PetscFree(remoteOffsetsParents);CHKERRQ(ierr);
1184a6f36705SMichael Lange     ierr = PetscSectionGetStorageSize(newParentSection,&newParentSize);CHKERRQ(ierr);
1185a6f36705SMichael Lange     ierr = PetscMalloc2(newParentSize,&newParents,newParentSize,&newChildIDs);CHKERRQ(ierr);
1186a6f36705SMichael Lange     ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
1187a6f36705SMichael Lange     ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
1188a6f36705SMichael Lange     ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1189a6f36705SMichael Lange     ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1190a6f36705SMichael Lange     ierr = ISGlobalToLocalMappingApplyBlock(renumbering,IS_GTOLM_MASK, newParentSize, newParents, NULL, newParents);CHKERRQ(ierr);
11914a54e071SToby Isaac #if PETSC_USE_DEBUG
11924a54e071SToby Isaac     {
11934a54e071SToby Isaac       PetscInt  p;
11944a54e071SToby Isaac       PetscBool valid = PETSC_TRUE;
11954a54e071SToby Isaac       for (p = 0; p < newParentSize; ++p) {
11964a54e071SToby Isaac         if (newParents[p] < 0) {valid = PETSC_FALSE; ierr = PetscPrintf(PETSC_COMM_SELF, "Point %d not in overlap SF\n", p);CHKERRQ(ierr);}
11974a54e071SToby Isaac       }
11984a54e071SToby Isaac       if (!valid) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map");
11994a54e071SToby Isaac     }
12004a54e071SToby Isaac #endif
1201a6f36705SMichael Lange     ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-parents_view", &flg);CHKERRQ(ierr);
1202a6f36705SMichael Lange     if (flg) {
1203a6f36705SMichael Lange       ierr = PetscPrintf(comm, "Serial Parent Section: \n");CHKERRQ(ierr);
1204a6f36705SMichael Lange       ierr = PetscSectionView(origParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1205a6f36705SMichael Lange       ierr = PetscPrintf(comm, "Parallel Parent Section: \n");CHKERRQ(ierr);
1206a6f36705SMichael Lange       ierr = PetscSectionView(newParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1207a6f36705SMichael Lange       ierr = PetscSFView(parentSF, NULL);CHKERRQ(ierr);
1208a6f36705SMichael Lange     }
1209a6f36705SMichael Lange     ierr = DMPlexSetTree(dmParallel,newParentSection,newParents,newChildIDs);CHKERRQ(ierr);
1210a6f36705SMichael Lange     ierr = PetscSectionDestroy(&newParentSection);CHKERRQ(ierr);
1211a6f36705SMichael Lange     ierr = PetscFree2(newParents,newChildIDs);CHKERRQ(ierr);
1212a6f36705SMichael Lange     ierr = PetscSFDestroy(&parentSF);CHKERRQ(ierr);
1213a6f36705SMichael Lange   }
121415078cd4SMichael Lange   pmesh->useAnchors = mesh->useAnchors;
1215a6f36705SMichael Lange   PetscFunctionReturn(0);
1216a6f36705SMichael Lange }
12170df0e737SMichael Lange 
12180df0e737SMichael Lange #undef __FUNCT__
12194eca1733SMichael Lange #define __FUNCT__ "DMPlexDistributeSF"
1220f8987ae8SMichael Lange PetscErrorCode DMPlexDistributeSF(DM dm, PetscSF migrationSF, DM dmParallel)
12214eca1733SMichael Lange {
12224eca1733SMichael Lange   DM_Plex               *mesh  = (DM_Plex*) dm->data;
12234eca1733SMichael Lange   DM_Plex               *pmesh = (DM_Plex*) (dmParallel)->data;
12244eca1733SMichael Lange   PetscMPIInt            rank, numProcs;
12254eca1733SMichael Lange   MPI_Comm               comm;
12264eca1733SMichael Lange   PetscErrorCode         ierr;
12274eca1733SMichael Lange 
12284eca1733SMichael Lange   PetscFunctionBegin;
12294eca1733SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12304eca1733SMichael Lange   PetscValidPointer(dmParallel,7);
12314eca1733SMichael Lange 
12324eca1733SMichael Lange   /* Create point SF for parallel mesh */
12334eca1733SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
12344eca1733SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
12354eca1733SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
12364eca1733SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
12374eca1733SMichael Lange   {
12384eca1733SMichael Lange     const PetscInt *leaves;
12394eca1733SMichael Lange     PetscSFNode    *remotePoints, *rowners, *lowners;
12404eca1733SMichael Lange     PetscInt        numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints;
12414eca1733SMichael Lange     PetscInt        pStart, pEnd;
12424eca1733SMichael Lange 
12434eca1733SMichael Lange     ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
12444eca1733SMichael Lange     ierr = PetscSFGetGraph(migrationSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr);
12454eca1733SMichael Lange     ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr);
12464eca1733SMichael Lange     for (p=0; p<numRoots; p++) {
12474eca1733SMichael Lange       rowners[p].rank  = -1;
12484eca1733SMichael Lange       rowners[p].index = -1;
12494eca1733SMichael Lange     }
12504eca1733SMichael Lange     ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12514eca1733SMichael Lange     ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12524eca1733SMichael Lange     for (p = 0; p < numLeaves; ++p) {
12534eca1733SMichael Lange       if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */
12544eca1733SMichael Lange         lowners[p].rank  = rank;
12554eca1733SMichael Lange         lowners[p].index = leaves ? leaves[p] : p;
12564eca1733SMichael Lange       } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */
12574eca1733SMichael Lange         lowners[p].rank  = -2;
12584eca1733SMichael Lange         lowners[p].index = -2;
12594eca1733SMichael Lange       }
12604eca1733SMichael Lange     }
12614eca1733SMichael Lange     for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */
12624eca1733SMichael Lange       rowners[p].rank  = -3;
12634eca1733SMichael Lange       rowners[p].index = -3;
12644eca1733SMichael Lange     }
12654eca1733SMichael Lange     ierr = PetscSFReduceBegin(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
12664eca1733SMichael Lange     ierr = PetscSFReduceEnd(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
12674eca1733SMichael Lange     ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12684eca1733SMichael Lange     ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12694eca1733SMichael Lange     for (p = 0; p < numLeaves; ++p) {
12704eca1733SMichael Lange       if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed");
12714eca1733SMichael Lange       if (lowners[p].rank != rank) ++numGhostPoints;
12724eca1733SMichael Lange     }
12734eca1733SMichael Lange     ierr = PetscMalloc1(numGhostPoints, &ghostPoints);CHKERRQ(ierr);
12744eca1733SMichael Lange     ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr);
12754eca1733SMichael Lange     for (p = 0, gp = 0; p < numLeaves; ++p) {
12764eca1733SMichael Lange       if (lowners[p].rank != rank) {
12774eca1733SMichael Lange         ghostPoints[gp]        = leaves ? leaves[p] : p;
12784eca1733SMichael Lange         remotePoints[gp].rank  = lowners[p].rank;
12794eca1733SMichael Lange         remotePoints[gp].index = lowners[p].index;
12804eca1733SMichael Lange         ++gp;
12814eca1733SMichael Lange       }
12824eca1733SMichael Lange     }
12834eca1733SMichael Lange     ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr);
12844eca1733SMichael Lange     ierr = PetscSFSetGraph((dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
12854eca1733SMichael Lange     ierr = PetscSFSetFromOptions((dmParallel)->sf);CHKERRQ(ierr);
12864eca1733SMichael Lange   }
12874eca1733SMichael Lange   pmesh->useCone    = mesh->useCone;
12884eca1733SMichael Lange   pmesh->useClosure = mesh->useClosure;
12894eca1733SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
12904eca1733SMichael Lange   PetscFunctionReturn(0);
12914eca1733SMichael Lange }
12924eca1733SMichael Lange 
12934eca1733SMichael Lange #undef __FUNCT__
1294f5bf2dbfSMichael Lange #define __FUNCT__ "DMPlexCreatePointSF"
1295f5bf2dbfSMichael Lange /*@C
1296f5bf2dbfSMichael Lange   DMPlexDerivePointSF - Build a point SF from an SF describing a point migration
1297f5bf2dbfSMichael Lange 
1298f5bf2dbfSMichael Lange   Input Parameter:
1299f5bf2dbfSMichael Lange + dm          - The source DMPlex object
1300f5bf2dbfSMichael Lange . migrationSF - The star forest that describes the parallel point remapping
13011627f6ccSMichael Lange . ownership   - Flag causing a vote to determine point ownership
1302f5bf2dbfSMichael Lange 
1303f5bf2dbfSMichael Lange   Output Parameter:
1304f5bf2dbfSMichael Lange - pointSF     - The star forest describing the point overlap in the remapped DM
1305f5bf2dbfSMichael Lange 
1306f5bf2dbfSMichael Lange   Level: developer
1307f5bf2dbfSMichael Lange 
1308f5bf2dbfSMichael Lange .seealso: DMPlexDistribute(), DMPlexDistributeOverlap()
1309f5bf2dbfSMichael Lange @*/
1310f5bf2dbfSMichael Lange PetscErrorCode DMPlexCreatePointSF(DM dm, PetscSF migrationSF, PetscBool ownership, PetscSF *pointSF)
1311f5bf2dbfSMichael Lange {
1312f5bf2dbfSMichael Lange   PetscMPIInt        rank;
13131627f6ccSMichael Lange   PetscInt           p, nroots, nleaves, idx, npointLeaves;
1314f5bf2dbfSMichael Lange   PetscInt          *pointLocal;
1315f5bf2dbfSMichael Lange   const PetscInt    *leaves;
1316f5bf2dbfSMichael Lange   const PetscSFNode *roots;
1317f5bf2dbfSMichael Lange   PetscSFNode       *rootNodes, *leafNodes, *pointRemote;
1318f5bf2dbfSMichael Lange   PetscErrorCode     ierr;
1319f5bf2dbfSMichael Lange 
1320f5bf2dbfSMichael Lange   PetscFunctionBegin;
1321f5bf2dbfSMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1322f5bf2dbfSMichael Lange   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
1323f5bf2dbfSMichael Lange 
1324f5bf2dbfSMichael Lange   ierr = PetscSFGetGraph(migrationSF, &nroots, &nleaves, &leaves, &roots);CHKERRQ(ierr);
1325f5bf2dbfSMichael Lange   ierr = PetscMalloc2(nroots, &rootNodes, nleaves, &leafNodes);CHKERRQ(ierr);
1326f5bf2dbfSMichael Lange   if (ownership) {
13271627f6ccSMichael Lange     /* Point ownership vote: Process with highest rank ownes shared points */
1328f5bf2dbfSMichael Lange     for (p = 0; p < nleaves; ++p) {
1329f5bf2dbfSMichael Lange       /* Either put in a bid or we know we own it */
1330f5bf2dbfSMichael Lange       leafNodes[p].rank  = rank;
133143f7d02bSMichael Lange       leafNodes[p].index = p;
1332f5bf2dbfSMichael Lange     }
1333f5bf2dbfSMichael Lange     for (p = 0; p < nroots; p++) {
13341627f6ccSMichael Lange       /* Root must not participate in the reduction, flag so that MAXLOC does not use */
1335f5bf2dbfSMichael Lange       rootNodes[p].rank  = -3;
1336f5bf2dbfSMichael Lange       rootNodes[p].index = -3;
1337f5bf2dbfSMichael Lange     }
1338f5bf2dbfSMichael Lange     ierr = PetscSFReduceBegin(migrationSF, MPIU_2INT, leafNodes, rootNodes, MPI_MAXLOC);CHKERRQ(ierr);
1339f5bf2dbfSMichael Lange     ierr = PetscSFReduceEnd(migrationSF, MPIU_2INT, leafNodes, rootNodes, MPI_MAXLOC);CHKERRQ(ierr);
1340f5bf2dbfSMichael Lange   } else {
1341f5bf2dbfSMichael Lange     for (p = 0; p < nroots; p++) {
1342f5bf2dbfSMichael Lange       rootNodes[p].index = -1;
1343f5bf2dbfSMichael Lange       rootNodes[p].rank = rank;
1344f5bf2dbfSMichael Lange     };
1345f5bf2dbfSMichael Lange     for (p = 0; p < nleaves; p++) {
1346f5bf2dbfSMichael Lange       /* Write new local id into old location */
1347f5bf2dbfSMichael Lange       if (roots[p].rank == rank) {
1348f5bf2dbfSMichael Lange         rootNodes[roots[p].index].index = leaves[p];
1349f5bf2dbfSMichael Lange       }
1350f5bf2dbfSMichael Lange     }
1351f5bf2dbfSMichael Lange   }
1352f5bf2dbfSMichael Lange   ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rootNodes, leafNodes);CHKERRQ(ierr);
1353f5bf2dbfSMichael Lange   ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rootNodes, leafNodes);CHKERRQ(ierr);
1354f5bf2dbfSMichael Lange 
1355f5bf2dbfSMichael Lange   for (npointLeaves = 0, p = 0; p < nleaves; p++) {if (leafNodes[p].rank != rank) npointLeaves++;}
13561627f6ccSMichael Lange   ierr = PetscMalloc1(npointLeaves, &pointLocal);CHKERRQ(ierr);
13571627f6ccSMichael Lange   ierr = PetscMalloc1(npointLeaves, &pointRemote);CHKERRQ(ierr);
1358f5bf2dbfSMichael Lange   for (idx = 0, p = 0; p < nleaves; p++) {
1359f5bf2dbfSMichael Lange     if (leafNodes[p].rank != rank) {
1360f5bf2dbfSMichael Lange       pointLocal[idx] = p;
1361f5bf2dbfSMichael Lange       pointRemote[idx] = leafNodes[p];
1362f5bf2dbfSMichael Lange       idx++;
1363f5bf2dbfSMichael Lange     }
1364f5bf2dbfSMichael Lange   }
1365f5bf2dbfSMichael Lange   ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), pointSF);CHKERRQ(ierr);
13661627f6ccSMichael Lange   ierr = PetscSFSetFromOptions(*pointSF);CHKERRQ(ierr);
1367f5bf2dbfSMichael Lange   ierr = PetscSFSetGraph(*pointSF, nleaves, npointLeaves, pointLocal, PETSC_OWN_POINTER, pointRemote, PETSC_OWN_POINTER);CHKERRQ(ierr);
1368f5bf2dbfSMichael Lange   ierr = PetscFree2(rootNodes, leafNodes);CHKERRQ(ierr);
1369f5bf2dbfSMichael Lange   PetscFunctionReturn(0);
1370f5bf2dbfSMichael Lange }
1371f5bf2dbfSMichael Lange 
1372f5bf2dbfSMichael Lange #undef __FUNCT__
137315078cd4SMichael Lange #define __FUNCT__ "DMPlexMigrate"
137415078cd4SMichael Lange /*@C
137515078cd4SMichael Lange   DMPlexMigrate  - Migrates internal DM data over the supplied star forest
137615078cd4SMichael Lange 
137715078cd4SMichael Lange   Input Parameter:
137815078cd4SMichael Lange + dm       - The source DMPlex object
13791627f6ccSMichael Lange . sf       - The star forest communication context describing the migration pattern
138015078cd4SMichael Lange 
138115078cd4SMichael Lange   Output Parameter:
138215078cd4SMichael Lange - targetDM - The target DMPlex object
138315078cd4SMichael Lange 
1384b9f40539SMichael Lange   Level: intermediate
138515078cd4SMichael Lange 
138615078cd4SMichael Lange .seealso: DMPlexDistribute(), DMPlexDistributeOverlap()
138715078cd4SMichael Lange @*/
1388b9f40539SMichael Lange PetscErrorCode DMPlexMigrate(DM dm, PetscSF sf, DM targetDM)
138915078cd4SMichael Lange {
1390b9f40539SMichael Lange   MPI_Comm               comm;
1391b9f40539SMichael Lange   PetscInt               dim, nroots;
1392b9f40539SMichael Lange   PetscSF                sfPoint;
139315078cd4SMichael Lange   ISLocalToGlobalMapping ltogMigration;
1394b9f40539SMichael Lange   PetscBool              flg;
139515078cd4SMichael Lange   PetscErrorCode         ierr;
139615078cd4SMichael Lange 
139715078cd4SMichael Lange   PetscFunctionBegin;
139815078cd4SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
13991b858b30SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_Migrate, dm, 0, 0, 0);CHKERRQ(ierr);
1400b9f40539SMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
140115078cd4SMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
140215078cd4SMichael Lange   ierr = DMSetDimension(targetDM, dim);CHKERRQ(ierr);
140315078cd4SMichael Lange 
1404bfb0467fSMichael Lange   /* Check for a one-to-all distribution pattern */
1405b9f40539SMichael Lange   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
1406b9f40539SMichael Lange   ierr = PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1407bfb0467fSMichael Lange   if (nroots >= 0) {
1408b9f40539SMichael Lange     IS                     isOriginal;
1409b9f40539SMichael Lange     ISLocalToGlobalMapping ltogOriginal;
1410b9f40539SMichael Lange     PetscInt               n, size, nleaves, conesSize;
1411b9f40539SMichael Lange     PetscInt              *numbering_orig, *numbering_new, *cones;
141215078cd4SMichael Lange     PetscSection           coneSection;
1413b9f40539SMichael Lange     /* Get the original point numbering */
1414b9f40539SMichael Lange     ierr = DMPlexCreatePointNumbering(dm, &isOriginal);CHKERRQ(ierr);
1415b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingCreateIS(isOriginal, &ltogOriginal);CHKERRQ(ierr);
1416b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingGetSize(ltogOriginal, &size);CHKERRQ(ierr);
1417b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingGetBlockIndices(ltogOriginal, (const PetscInt**)&numbering_orig);CHKERRQ(ierr);
1418b9f40539SMichael Lange     /* Convert to positive global numbers */
1419b9f40539SMichael Lange     for (n=0; n<size; n++) {if (numbering_orig[n] < 0) numbering_orig[n] = -(numbering_orig[n]+1);}
1420b9f40539SMichael Lange     /* Derive the new local-to-global mapping from the old one */
1421b9f40539SMichael Lange     ierr = PetscSFGetGraph(sf, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr);
1422b9f40539SMichael Lange     ierr = PetscMalloc1(nleaves, &numbering_new);CHKERRQ(ierr);
1423b9f40539SMichael Lange     ierr = PetscSFBcastBegin(sf, MPIU_INT, (PetscInt *) numbering_orig, numbering_new);CHKERRQ(ierr);
1424b9f40539SMichael Lange     ierr = PetscSFBcastEnd(sf, MPIU_INT, (PetscInt *) numbering_orig, numbering_new);CHKERRQ(ierr);
1425b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingCreate(comm, 1, nleaves, (const PetscInt*) numbering_new, PETSC_OWN_POINTER, &ltogMigration);CHKERRQ(ierr);
1426b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingRestoreIndices(ltogOriginal, (const PetscInt**)&numbering_orig);CHKERRQ(ierr);
142715078cd4SMichael Lange     /* Convert cones to global numbering before migrating them */
142815078cd4SMichael Lange     ierr = DMPlexGetConeSection(dm, &coneSection);CHKERRQ(ierr);
142915078cd4SMichael Lange     ierr = PetscSectionGetStorageSize(coneSection, &conesSize);CHKERRQ(ierr);
143015078cd4SMichael Lange     ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
1431b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingApplyBlock(ltogOriginal, conesSize, cones, cones);CHKERRQ(ierr);
1432b9f40539SMichael Lange     ierr = ISDestroy(&isOriginal);CHKERRQ(ierr);
1433302440fdSBarry Smith     ierr = ISLocalToGlobalMappingDestroy(&ltogOriginal);CHKERRQ(ierr);
143415078cd4SMichael Lange   } else {
1435bfb0467fSMichael Lange     /* One-to-all distribution pattern: We can derive LToG from SF */
143615078cd4SMichael Lange     ierr = ISLocalToGlobalMappingCreateSF(sf, 0, &ltogMigration);CHKERRQ(ierr);
143715078cd4SMichael Lange   }
1438b9f40539SMichael Lange   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
1439b9f40539SMichael Lange   if (flg) {
1440b9f40539SMichael Lange     ierr = PetscPrintf(comm, "Point renumbering for DM migration:\n");CHKERRQ(ierr);
1441b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingView(ltogMigration, NULL);CHKERRQ(ierr);
1442b9f40539SMichael Lange   }
144315078cd4SMichael Lange   /* Migrate DM data to target DM */
144415078cd4SMichael Lange   ierr = DMPlexDistributeCones(dm, sf, ltogMigration, targetDM);CHKERRQ(ierr);
144515078cd4SMichael Lange   ierr = DMPlexDistributeCoordinates(dm, sf, targetDM);CHKERRQ(ierr);
144615078cd4SMichael Lange   ierr = DMPlexDistributeLabels(dm, sf, targetDM);CHKERRQ(ierr);
144715078cd4SMichael Lange   ierr = DMPlexDistributeSetupHybrid(dm, sf, ltogMigration, targetDM);CHKERRQ(ierr);
144815078cd4SMichael Lange   ierr = DMPlexDistributeSetupTree(dm, sf, ltogMigration, targetDM);CHKERRQ(ierr);
1449302440fdSBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&ltogMigration);CHKERRQ(ierr);
14501b858b30SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_Migrate, dm, 0, 0, 0);CHKERRQ(ierr);
145115078cd4SMichael Lange   PetscFunctionReturn(0);
145215078cd4SMichael Lange }
145315078cd4SMichael Lange 
145415078cd4SMichael Lange #undef __FUNCT__
145570034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistribute"
14563b8f15a2SMatthew G. Knepley /*@C
145770034214SMatthew G. Knepley   DMPlexDistribute - Distributes the mesh and any associated sections.
145870034214SMatthew G. Knepley 
145970034214SMatthew G. Knepley   Not Collective
146070034214SMatthew G. Knepley 
146170034214SMatthew G. Knepley   Input Parameter:
146270034214SMatthew G. Knepley + dm  - The original DMPlex object
146370034214SMatthew G. Knepley - overlap - The overlap of partitions, 0 is the default
146470034214SMatthew G. Knepley 
146570034214SMatthew G. Knepley   Output Parameter:
146670034214SMatthew G. Knepley + sf - The PetscSF used for point distribution
146770034214SMatthew G. Knepley - parallelMesh - The distributed DMPlex object, or NULL
146870034214SMatthew G. Knepley 
146970034214SMatthew G. Knepley   Note: If the mesh was not distributed, the return value is NULL.
147070034214SMatthew G. Knepley 
147170034214SMatthew G. Knepley   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
147270034214SMatthew G. Knepley   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
147370034214SMatthew G. Knepley   representation on the mesh.
147470034214SMatthew G. Knepley 
147570034214SMatthew G. Knepley   Level: intermediate
147670034214SMatthew G. Knepley 
147770034214SMatthew G. Knepley .keywords: mesh, elements
147870034214SMatthew G. Knepley .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
147970034214SMatthew G. Knepley @*/
148080cf41d5SMatthew G. Knepley PetscErrorCode DMPlexDistribute(DM dm, PetscInt overlap, PetscSF *sf, DM *dmParallel)
148170034214SMatthew G. Knepley {
148270034214SMatthew G. Knepley   MPI_Comm               comm;
148315078cd4SMichael Lange   PetscPartitioner       partitioner;
1484f8987ae8SMichael Lange   IS                     cellPart;
1485f8987ae8SMichael Lange   PetscSection           cellPartSection;
1486cf86098cSMatthew G. Knepley   DM                     dmCoord;
1487f8987ae8SMichael Lange   DMLabel                lblPartition, lblMigration;
148843f7d02bSMichael Lange   PetscSF                sfProcess, sfMigration, sfStratified, sfPoint;
148970034214SMatthew G. Knepley   PetscBool              flg;
149070034214SMatthew G. Knepley   PetscMPIInt            rank, numProcs, p;
149170034214SMatthew G. Knepley   PetscErrorCode         ierr;
149270034214SMatthew G. Knepley 
149370034214SMatthew G. Knepley   PetscFunctionBegin;
149470034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
149570034214SMatthew G. Knepley   if (sf) PetscValidPointer(sf,4);
149670034214SMatthew G. Knepley   PetscValidPointer(dmParallel,5);
149770034214SMatthew G. Knepley 
149870034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
149970034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
150070034214SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
150170034214SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
150270034214SMatthew G. Knepley 
150370034214SMatthew G. Knepley   *dmParallel = NULL;
150470034214SMatthew G. Knepley   if (numProcs == 1) PetscFunctionReturn(0);
150570034214SMatthew G. Knepley 
150615078cd4SMichael Lange   /* Create cell partition */
150777623264SMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
150877623264SMatthew G. Knepley   ierr = PetscSectionCreate(comm, &cellPartSection);CHKERRQ(ierr);
150915078cd4SMichael Lange   ierr = DMPlexGetPartitioner(dm, &partitioner);CHKERRQ(ierr);
151015078cd4SMichael Lange   ierr = PetscPartitionerPartition(partitioner, dm, cellPartSection, &cellPart);CHKERRQ(ierr);
1511f8987ae8SMichael Lange   {
1512f8987ae8SMichael Lange     /* Convert partition to DMLabel */
1513f8987ae8SMichael Lange     PetscInt proc, pStart, pEnd, npoints, poffset;
1514f8987ae8SMichael Lange     const PetscInt *points;
1515f8987ae8SMichael Lange     ierr = DMLabelCreate("Point Partition", &lblPartition);CHKERRQ(ierr);
1516f8987ae8SMichael Lange     ierr = ISGetIndices(cellPart, &points);CHKERRQ(ierr);
1517f8987ae8SMichael Lange     ierr = PetscSectionGetChart(cellPartSection, &pStart, &pEnd);CHKERRQ(ierr);
1518f8987ae8SMichael Lange     for (proc = pStart; proc < pEnd; proc++) {
1519f8987ae8SMichael Lange       ierr = PetscSectionGetDof(cellPartSection, proc, &npoints);CHKERRQ(ierr);
1520f8987ae8SMichael Lange       ierr = PetscSectionGetOffset(cellPartSection, proc, &poffset);CHKERRQ(ierr);
1521f8987ae8SMichael Lange       for (p = poffset; p < poffset+npoints; p++) {
1522f8987ae8SMichael Lange         ierr = DMLabelSetValue(lblPartition, points[p], proc);CHKERRQ(ierr);
152370034214SMatthew G. Knepley       }
1524f8987ae8SMichael Lange     }
1525f8987ae8SMichael Lange     ierr = ISRestoreIndices(cellPart, &points);CHKERRQ(ierr);
1526f8987ae8SMichael Lange   }
1527f8987ae8SMichael Lange   ierr = DMPlexPartitionLabelClosure(dm, lblPartition);CHKERRQ(ierr);
1528f8987ae8SMichael Lange   {
1529f8987ae8SMichael Lange     /* Build a global process SF */
1530f8987ae8SMichael Lange     PetscSFNode *remoteProc;
1531f8987ae8SMichael Lange     ierr = PetscMalloc1(numProcs, &remoteProc);CHKERRQ(ierr);
1532f8987ae8SMichael Lange     for (p = 0; p < numProcs; ++p) {
1533f8987ae8SMichael Lange       remoteProc[p].rank  = p;
1534f8987ae8SMichael Lange       remoteProc[p].index = rank;
1535f8987ae8SMichael Lange     }
1536f8987ae8SMichael Lange     ierr = PetscSFCreate(comm, &sfProcess);CHKERRQ(ierr);
1537f8987ae8SMichael Lange     ierr = PetscObjectSetName((PetscObject) sfProcess, "Process SF");CHKERRQ(ierr);
1538f8987ae8SMichael Lange     ierr = PetscSFSetGraph(sfProcess, numProcs, numProcs, NULL, PETSC_OWN_POINTER, remoteProc, PETSC_OWN_POINTER);CHKERRQ(ierr);
1539f8987ae8SMichael Lange   }
1540f8987ae8SMichael Lange   ierr = DMLabelCreate("Point migration", &lblMigration);CHKERRQ(ierr);
1541f8987ae8SMichael Lange   ierr = DMPlexPartitionLabelInvert(dm, lblPartition, sfProcess, lblMigration);CHKERRQ(ierr);
1542302440fdSBarry Smith   ierr = DMPlexPartitionLabelCreateSF(dm, lblMigration, &sfMigration);CHKERRQ(ierr);
154343f7d02bSMichael Lange   /* Stratify the SF in case we are migrating an already parallel plex */
154443f7d02bSMichael Lange   ierr = DMPlexStratifyMigrationSF(dm, sfMigration, &sfStratified);CHKERRQ(ierr);
154543f7d02bSMichael Lange   ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr);
154643f7d02bSMichael Lange   sfMigration = sfStratified;
1547f8987ae8SMichael Lange   ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
154870034214SMatthew G. Knepley   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
154970034214SMatthew G. Knepley   if (flg) {
1550f8987ae8SMichael Lange     ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr);
1551f8987ae8SMichael Lange     ierr = DMLabelView(lblPartition, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1552f8987ae8SMichael Lange     ierr = PetscSFView(sfMigration, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
155370034214SMatthew G. Knepley   }
1554f8987ae8SMichael Lange 
155515078cd4SMichael Lange   /* Create non-overlapping parallel DM and migrate internal data */
155670034214SMatthew G. Knepley   ierr = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr);
155770034214SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr);
1558b9f40539SMichael Lange   ierr = DMPlexMigrate(dm, sfMigration, *dmParallel);CHKERRQ(ierr);
155970034214SMatthew G. Knepley 
1560a157612eSMichael Lange   /* Build the point SF without overlap */
15611627f6ccSMichael Lange   ierr = DMPlexCreatePointSF(*dmParallel, sfMigration, PETSC_TRUE, &sfPoint);CHKERRQ(ierr);
1562f5bf2dbfSMichael Lange   ierr = DMSetPointSF(*dmParallel, sfPoint);CHKERRQ(ierr);
1563cf86098cSMatthew G. Knepley   ierr = DMGetCoordinateDM(*dmParallel, &dmCoord);CHKERRQ(ierr);
1564cf86098cSMatthew G. Knepley   if (dmCoord) {ierr = DMSetPointSF(dmCoord, sfPoint);CHKERRQ(ierr);}
1565f5bf2dbfSMichael Lange   if (flg) {ierr = PetscSFView(sfPoint, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
156670034214SMatthew G. Knepley 
1567a157612eSMichael Lange   if (overlap > 0) {
156815078cd4SMichael Lange     DM                 dmOverlap;
156915078cd4SMichael Lange     PetscInt           nroots, nleaves;
157015078cd4SMichael Lange     PetscSFNode       *newRemote;
157115078cd4SMichael Lange     const PetscSFNode *oldRemote;
157215078cd4SMichael Lange     PetscSF            sfOverlap, sfOverlapPoint;
1573a157612eSMichael Lange     /* Add the partition overlap to the distributed DM */
1574b9f40539SMichael Lange     ierr = DMPlexDistributeOverlap(*dmParallel, overlap, &sfOverlap, &dmOverlap);CHKERRQ(ierr);
1575a157612eSMichael Lange     ierr = DMDestroy(dmParallel);CHKERRQ(ierr);
1576a157612eSMichael Lange     *dmParallel = dmOverlap;
1577c389ea9fSToby Isaac     if (flg) {
15783d822a50SMichael Lange       ierr = PetscPrintf(comm, "Overlap Migration SF:\n");CHKERRQ(ierr);
157915078cd4SMichael Lange       ierr = PetscSFView(sfOverlap, NULL);CHKERRQ(ierr);
1580c389ea9fSToby Isaac     }
158143331d4aSMichael Lange 
1582f8987ae8SMichael Lange     /* Re-map the migration SF to establish the full migration pattern */
1583f8987ae8SMichael Lange     ierr = PetscSFGetGraph(sfMigration, &nroots, NULL, NULL, &oldRemote);CHKERRQ(ierr);
158415078cd4SMichael Lange     ierr = PetscSFGetGraph(sfOverlap, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr);
158543331d4aSMichael Lange     ierr = PetscMalloc1(nleaves, &newRemote);CHKERRQ(ierr);
158615078cd4SMichael Lange     ierr = PetscSFBcastBegin(sfOverlap, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr);
158715078cd4SMichael Lange     ierr = PetscSFBcastEnd(sfOverlap, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr);
158815078cd4SMichael Lange     ierr = PetscSFCreate(comm, &sfOverlapPoint);CHKERRQ(ierr);
158915078cd4SMichael Lange     ierr = PetscSFSetGraph(sfOverlapPoint, nroots, nleaves, NULL, PETSC_OWN_POINTER, newRemote, PETSC_OWN_POINTER);CHKERRQ(ierr);
159015078cd4SMichael Lange     ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr);
1591f8987ae8SMichael Lange     ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr);
159215078cd4SMichael Lange     sfMigration = sfOverlapPoint;
1593c389ea9fSToby Isaac   }
1594bf5244c7SMatthew G. Knepley   /* Cleanup Partition */
1595f8987ae8SMichael Lange   ierr = PetscSFDestroy(&sfProcess);CHKERRQ(ierr);
1596f8987ae8SMichael Lange   ierr = DMLabelDestroy(&lblPartition);CHKERRQ(ierr);
1597f8987ae8SMichael Lange   ierr = DMLabelDestroy(&lblMigration);CHKERRQ(ierr);
15984eca1733SMichael Lange   ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr);
15994eca1733SMichael Lange   ierr = ISDestroy(&cellPart);CHKERRQ(ierr);
1600721cbd76SMatthew G. Knepley   /* Copy BC */
1601721cbd76SMatthew G. Knepley   ierr = DMPlexCopyBoundary(dm, *dmParallel);CHKERRQ(ierr);
1602721cbd76SMatthew G. Knepley   /* Cleanup */
1603f8987ae8SMichael Lange   if (sf) {*sf = sfMigration;}
1604f8987ae8SMichael Lange   else    {ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr);}
1605f5bf2dbfSMichael Lange   ierr = PetscSFDestroy(&sfPoint);CHKERRQ(ierr);
160670034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
160770034214SMatthew G. Knepley   PetscFunctionReturn(0);
160870034214SMatthew G. Knepley }
1609a157612eSMichael Lange 
1610a157612eSMichael Lange #undef __FUNCT__
1611a157612eSMichael Lange #define __FUNCT__ "DMPlexDistributeOverlap"
1612a157612eSMichael Lange /*@C
1613a157612eSMichael Lange   DMPlexDistribute - Add partition overlap to a distributed non-overlapping DM.
1614a157612eSMichael Lange 
1615a157612eSMichael Lange   Not Collective
1616a157612eSMichael Lange 
1617a157612eSMichael Lange   Input Parameter:
1618a157612eSMichael Lange + dm  - The non-overlapping distrbuted DMPlex object
1619a157612eSMichael Lange - overlap - The overlap of partitions, 0 is the default
1620a157612eSMichael Lange 
1621a157612eSMichael Lange   Output Parameter:
1622a157612eSMichael Lange + sf - The PetscSF used for point distribution
1623a157612eSMichael Lange - dmOverlap - The overlapping distributed DMPlex object, or NULL
1624a157612eSMichael Lange 
1625a157612eSMichael Lange   Note: If the mesh was not distributed, the return value is NULL.
1626a157612eSMichael Lange 
1627a157612eSMichael Lange   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
1628a157612eSMichael Lange   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
1629a157612eSMichael Lange   representation on the mesh.
1630a157612eSMichael Lange 
1631a157612eSMichael Lange   Level: intermediate
1632a157612eSMichael Lange 
1633a157612eSMichael Lange .keywords: mesh, elements
1634a157612eSMichael Lange .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
1635a157612eSMichael Lange @*/
1636b9f40539SMichael Lange PetscErrorCode DMPlexDistributeOverlap(DM dm, PetscInt overlap, PetscSF *sf, DM *dmOverlap)
1637a157612eSMichael Lange {
1638a157612eSMichael Lange   MPI_Comm               comm;
1639a157612eSMichael Lange   PetscMPIInt            rank;
1640a157612eSMichael Lange   PetscSection           rootSection, leafSection;
1641a157612eSMichael Lange   IS                     rootrank, leafrank;
1642cf86098cSMatthew G. Knepley   DM                     dmCoord;
1643a9f1d5b2SMichael Lange   DMLabel                lblOverlap;
1644f5bf2dbfSMichael Lange   PetscSF                sfOverlap, sfStratified, sfPoint;
1645a157612eSMichael Lange   PetscErrorCode         ierr;
1646a157612eSMichael Lange 
1647a157612eSMichael Lange   PetscFunctionBegin;
1648a157612eSMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1649a157612eSMichael Lange   if (sf) PetscValidPointer(sf, 3);
1650a157612eSMichael Lange   PetscValidPointer(dmOverlap, 4);
1651a157612eSMichael Lange 
16521b858b30SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeOverlap, dm, 0, 0, 0);CHKERRQ(ierr);
1653a157612eSMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1654a157612eSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1655a157612eSMichael Lange 
1656a157612eSMichael Lange   /* Compute point overlap with neighbouring processes on the distributed DM */
1657a157612eSMichael Lange   ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1658a157612eSMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
1659a157612eSMichael Lange   ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr);
1660a157612eSMichael Lange   ierr = DMPlexDistributeOwnership(dm, rootSection, &rootrank, leafSection, &leafrank);CHKERRQ(ierr);
1661a9f1d5b2SMichael Lange   ierr = DMPlexCreateOverlap(dm, overlap, rootSection, rootrank, leafSection, leafrank, &lblOverlap);CHKERRQ(ierr);
1662a9f1d5b2SMichael Lange   /* Convert overlap label to stratified migration SF */
1663a9f1d5b2SMichael Lange   ierr = DMPlexPartitionLabelCreateSF(dm, lblOverlap, &sfOverlap);CHKERRQ(ierr);
1664a9f1d5b2SMichael Lange   ierr = DMPlexStratifyMigrationSF(dm, sfOverlap, &sfStratified);CHKERRQ(ierr);
1665a9f1d5b2SMichael Lange   ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr);
1666a9f1d5b2SMichael Lange   sfOverlap = sfStratified;
1667a9f1d5b2SMichael Lange   ierr = PetscObjectSetName((PetscObject) sfOverlap, "Overlap SF");CHKERRQ(ierr);
1668a9f1d5b2SMichael Lange   ierr = PetscSFSetFromOptions(sfOverlap);CHKERRQ(ierr);
1669a9f1d5b2SMichael Lange 
167015fff7beSMatthew G. Knepley   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
167115fff7beSMatthew G. Knepley   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
167215fff7beSMatthew G. Knepley   ierr = ISDestroy(&rootrank);CHKERRQ(ierr);
167315fff7beSMatthew G. Knepley   ierr = ISDestroy(&leafrank);CHKERRQ(ierr);
1674a157612eSMichael Lange   ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1675a157612eSMichael Lange 
1676a157612eSMichael Lange   /* Build the overlapping DM */
1677a157612eSMichael Lange   ierr = DMPlexCreate(comm, dmOverlap);CHKERRQ(ierr);
1678a157612eSMichael Lange   ierr = PetscObjectSetName((PetscObject) *dmOverlap, "Parallel Mesh");CHKERRQ(ierr);
1679a9f1d5b2SMichael Lange   ierr = DMPlexMigrate(dm, sfOverlap, *dmOverlap);CHKERRQ(ierr);
1680f5bf2dbfSMichael Lange   /* Build the new point SF */
16811627f6ccSMichael Lange   ierr = DMPlexCreatePointSF(*dmOverlap, sfOverlap, PETSC_FALSE, &sfPoint);CHKERRQ(ierr);
1682f5bf2dbfSMichael Lange   ierr = DMSetPointSF(*dmOverlap, sfPoint);CHKERRQ(ierr);
1683cf86098cSMatthew G. Knepley   ierr = DMGetCoordinateDM(*dmOverlap, &dmCoord);CHKERRQ(ierr);
1684cf86098cSMatthew G. Knepley   if (dmCoord) {ierr = DMSetPointSF(dmCoord, sfPoint);CHKERRQ(ierr);}
1685f5bf2dbfSMichael Lange   ierr = PetscSFDestroy(&sfPoint);CHKERRQ(ierr);
1686a157612eSMichael Lange   /* Cleanup overlap partition */
1687a9f1d5b2SMichael Lange   ierr = DMLabelDestroy(&lblOverlap);CHKERRQ(ierr);
1688a9f1d5b2SMichael Lange   if (sf) *sf = sfOverlap;
1689a9f1d5b2SMichael Lange   else    {ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr);}
16901b858b30SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeOverlap, dm, 0, 0, 0);CHKERRQ(ierr);
1691a157612eSMichael Lange   PetscFunctionReturn(0);
1692a157612eSMichael Lange }
1693