xref: /petsc/src/dm/impls/plex/plexdistribute.c (revision 086331701050c551f5df120a53b850fb58faebb9)
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 }
374*08633170SToby Isaac 
375b0a623aaSMatthew G. Knepley #undef __FUNCT__
376b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexCreateTwoSidedProcessSF"
377b0a623aaSMatthew G. Knepley /*@
378b0a623aaSMatthew G. Knepley   DMPlexCreateTwoSidedProcessSF - Create an SF which just has process connectivity
379b0a623aaSMatthew G. Knepley 
380b0a623aaSMatthew G. Knepley   Collective on DM
381b0a623aaSMatthew G. Knepley 
382b0a623aaSMatthew G. Knepley   Input Parameters:
383b0a623aaSMatthew G. Knepley + dm      - The DM
384b0a623aaSMatthew G. Knepley - sfPoint - The PetscSF which encodes point connectivity
385b0a623aaSMatthew G. Knepley 
386b0a623aaSMatthew G. Knepley   Output Parameters:
387b0a623aaSMatthew G. Knepley + processRanks - A list of process neighbors, or NULL
388b0a623aaSMatthew G. Knepley - sfProcess    - An SF encoding the two-sided process connectivity, or NULL
389b0a623aaSMatthew G. Knepley 
390b0a623aaSMatthew G. Knepley   Level: developer
391b0a623aaSMatthew G. Knepley 
392b0a623aaSMatthew G. Knepley .seealso: PetscSFCreate(), DMPlexCreateProcessSF()
393b0a623aaSMatthew G. Knepley @*/
394b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM dm, PetscSF sfPoint, PetscSection rootRankSection, IS rootRanks, PetscSection leafRankSection, IS leafRanks, IS *processRanks, PetscSF *sfProcess)
395b0a623aaSMatthew G. Knepley {
396b0a623aaSMatthew G. Knepley   const PetscSFNode *remotePoints;
397b0a623aaSMatthew G. Knepley   PetscInt          *localPointsNew;
398b0a623aaSMatthew G. Knepley   PetscSFNode       *remotePointsNew;
399b0a623aaSMatthew G. Knepley   const PetscInt    *nranks;
400b0a623aaSMatthew G. Knepley   PetscInt          *ranksNew;
401b0a623aaSMatthew G. Knepley   PetscBT            neighbors;
402b0a623aaSMatthew G. Knepley   PetscInt           pStart, pEnd, p, numLeaves, l, numNeighbors, n;
403b0a623aaSMatthew G. Knepley   PetscMPIInt        numProcs, proc, rank;
404b0a623aaSMatthew G. Knepley   PetscErrorCode     ierr;
405b0a623aaSMatthew G. Knepley 
406b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
407b0a623aaSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
408b0a623aaSMatthew G. Knepley   PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2);
409b0a623aaSMatthew G. Knepley   if (processRanks) {PetscValidPointer(processRanks, 3);}
410b0a623aaSMatthew G. Knepley   if (sfProcess)    {PetscValidPointer(sfProcess, 4);}
411b0a623aaSMatthew G. Knepley   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr);
412b0a623aaSMatthew G. Knepley   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
413b0a623aaSMatthew G. Knepley   ierr = PetscSFGetGraph(sfPoint, NULL, &numLeaves, NULL, &remotePoints);CHKERRQ(ierr);
414b0a623aaSMatthew G. Knepley   ierr = PetscBTCreate(numProcs, &neighbors);CHKERRQ(ierr);
415b0a623aaSMatthew G. Knepley   ierr = PetscBTMemzero(numProcs, neighbors);CHKERRQ(ierr);
416b0a623aaSMatthew G. Knepley   /* Compute root-to-leaf process connectivity */
417b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetChart(rootRankSection, &pStart, &pEnd);CHKERRQ(ierr);
418b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(rootRanks, &nranks);CHKERRQ(ierr);
419b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
420b0a623aaSMatthew G. Knepley     PetscInt ndof, noff, n;
421b0a623aaSMatthew G. Knepley 
422b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetDof(rootRankSection, p, &ndof);CHKERRQ(ierr);
423b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(rootRankSection, p, &noff);CHKERRQ(ierr);
424302440fdSBarry Smith     for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);CHKERRQ(ierr);}
425b0a623aaSMatthew G. Knepley   }
426b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(rootRanks, &nranks);CHKERRQ(ierr);
427b0a623aaSMatthew G. Knepley   /* Compute leaf-to-neighbor process connectivity */
428b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetChart(leafRankSection, &pStart, &pEnd);CHKERRQ(ierr);
429b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(leafRanks, &nranks);CHKERRQ(ierr);
430b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
431b0a623aaSMatthew G. Knepley     PetscInt ndof, noff, n;
432b0a623aaSMatthew G. Knepley 
433b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetDof(leafRankSection, p, &ndof);CHKERRQ(ierr);
434b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(leafRankSection, p, &noff);CHKERRQ(ierr);
435302440fdSBarry Smith     for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);CHKERRQ(ierr);}
436b0a623aaSMatthew G. Knepley   }
437b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(leafRanks, &nranks);CHKERRQ(ierr);
438b0a623aaSMatthew G. Knepley   /* Compute leaf-to-root process connectivity */
439b0a623aaSMatthew G. Knepley   for (l = 0; l < numLeaves; ++l) {PetscBTSet(neighbors, remotePoints[l].rank);}
440b0a623aaSMatthew G. Knepley   /* Calculate edges */
441b0a623aaSMatthew G. Knepley   PetscBTClear(neighbors, rank);
442b0a623aaSMatthew G. Knepley   for(proc = 0, numNeighbors = 0; proc < numProcs; ++proc) {if (PetscBTLookup(neighbors, proc)) ++numNeighbors;}
443b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(numNeighbors, &ranksNew);CHKERRQ(ierr);
444b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(numNeighbors, &localPointsNew);CHKERRQ(ierr);
445b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(numNeighbors, &remotePointsNew);CHKERRQ(ierr);
446b0a623aaSMatthew G. Knepley   for(proc = 0, n = 0; proc < numProcs; ++proc) {
447b0a623aaSMatthew G. Knepley     if (PetscBTLookup(neighbors, proc)) {
448b0a623aaSMatthew G. Knepley       ranksNew[n]              = proc;
449b0a623aaSMatthew G. Knepley       localPointsNew[n]        = proc;
450b0a623aaSMatthew G. Knepley       remotePointsNew[n].index = rank;
451b0a623aaSMatthew G. Knepley       remotePointsNew[n].rank  = proc;
452b0a623aaSMatthew G. Knepley       ++n;
453b0a623aaSMatthew G. Knepley     }
454b0a623aaSMatthew G. Knepley   }
455b0a623aaSMatthew G. Knepley   ierr = PetscBTDestroy(&neighbors);CHKERRQ(ierr);
456b0a623aaSMatthew G. Knepley   if (processRanks) {ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), numNeighbors, ranksNew, PETSC_OWN_POINTER, processRanks);CHKERRQ(ierr);}
457b0a623aaSMatthew G. Knepley   else              {ierr = PetscFree(ranksNew);CHKERRQ(ierr);}
458b0a623aaSMatthew G. Knepley   if (sfProcess) {
459b0a623aaSMatthew G. Knepley     ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess);CHKERRQ(ierr);
460b0a623aaSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) *sfProcess, "Two-Sided Process SF");CHKERRQ(ierr);
461b0a623aaSMatthew G. Knepley     ierr = PetscSFSetFromOptions(*sfProcess);CHKERRQ(ierr);
462b0a623aaSMatthew G. Knepley     ierr = PetscSFSetGraph(*sfProcess, numProcs, numNeighbors, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr);
463b0a623aaSMatthew G. Knepley   }
464b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
465b0a623aaSMatthew G. Knepley }
466b0a623aaSMatthew G. Knepley 
467b0a623aaSMatthew G. Knepley #undef __FUNCT__
468b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeOwnership"
469b0a623aaSMatthew G. Knepley /*@
470b0a623aaSMatthew G. Knepley   DMPlexDistributeOwnership - Compute owner information for shared points. This basically gets two-sided for an SF.
471b0a623aaSMatthew G. Knepley 
472b0a623aaSMatthew G. Knepley   Collective on DM
473b0a623aaSMatthew G. Knepley 
474b0a623aaSMatthew G. Knepley   Input Parameter:
475b0a623aaSMatthew G. Knepley . dm - The DM
476b0a623aaSMatthew G. Knepley 
477b0a623aaSMatthew G. Knepley   Output Parameters:
478b0a623aaSMatthew G. Knepley + rootSection - The number of leaves for a given root point
479b0a623aaSMatthew G. Knepley . rootrank    - The rank of each edge into the root point
480b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point
481b0a623aaSMatthew G. Knepley - leafrank    - The rank of each process sharing a leaf point
482b0a623aaSMatthew G. Knepley 
483b0a623aaSMatthew G. Knepley   Level: developer
484b0a623aaSMatthew G. Knepley 
485b0a623aaSMatthew G. Knepley .seealso: DMPlexCreateOverlap()
486b0a623aaSMatthew G. Knepley @*/
487b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexDistributeOwnership(DM dm, PetscSection rootSection, IS *rootrank, PetscSection leafSection, IS *leafrank)
488b0a623aaSMatthew G. Knepley {
489b0a623aaSMatthew G. Knepley   MPI_Comm        comm;
490b0a623aaSMatthew G. Knepley   PetscSF         sfPoint;
491b0a623aaSMatthew G. Knepley   const PetscInt *rootdegree;
492b0a623aaSMatthew G. Knepley   PetscInt       *myrank, *remoterank;
493b0a623aaSMatthew G. Knepley   PetscInt        pStart, pEnd, p, nedges;
494b0a623aaSMatthew G. Knepley   PetscMPIInt     rank;
495b0a623aaSMatthew G. Knepley   PetscErrorCode  ierr;
496b0a623aaSMatthew G. Knepley 
497b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
498b0a623aaSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
499b0a623aaSMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
500b0a623aaSMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
501b0a623aaSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
502b0a623aaSMatthew G. Knepley   /* Compute number of leaves for each root */
503b0a623aaSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) rootSection, "Root Section");CHKERRQ(ierr);
504b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetChart(rootSection, pStart, pEnd);CHKERRQ(ierr);
505b0a623aaSMatthew G. Knepley   ierr = PetscSFComputeDegreeBegin(sfPoint, &rootdegree);CHKERRQ(ierr);
506b0a623aaSMatthew G. Knepley   ierr = PetscSFComputeDegreeEnd(sfPoint, &rootdegree);CHKERRQ(ierr);
507b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {ierr = PetscSectionSetDof(rootSection, p, rootdegree[p-pStart]);CHKERRQ(ierr);}
508b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
509b0a623aaSMatthew G. Knepley   /* Gather rank of each leaf to root */
510b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(rootSection, &nedges);CHKERRQ(ierr);
511b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(pEnd-pStart, &myrank);CHKERRQ(ierr);
512b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(nedges,  &remoterank);CHKERRQ(ierr);
513b0a623aaSMatthew G. Knepley   for (p = 0; p < pEnd-pStart; ++p) myrank[p] = rank;
514b0a623aaSMatthew G. Knepley   ierr = PetscSFGatherBegin(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr);
515b0a623aaSMatthew G. Knepley   ierr = PetscSFGatherEnd(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr);
516b0a623aaSMatthew G. Knepley   ierr = PetscFree(myrank);CHKERRQ(ierr);
517b0a623aaSMatthew G. Knepley   ierr = ISCreateGeneral(comm, nedges, remoterank, PETSC_OWN_POINTER, rootrank);CHKERRQ(ierr);
518b0a623aaSMatthew G. Knepley   /* Distribute remote ranks to leaves */
519b0a623aaSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) leafSection, "Leaf Section");CHKERRQ(ierr);
520b0a623aaSMatthew G. Knepley   ierr = DMPlexDistributeFieldIS(dm, sfPoint, rootSection, *rootrank, leafSection, leafrank);CHKERRQ(ierr);
521b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
522b0a623aaSMatthew G. Knepley }
523b0a623aaSMatthew G. Knepley 
524b0a623aaSMatthew G. Knepley #undef __FUNCT__
525b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexCreateOverlap"
526278397a0SMatthew G. Knepley /*@C
527b0a623aaSMatthew G. Knepley   DMPlexCreateOverlap - Compute owner information for shared points. This basically gets two-sided for an SF.
528b0a623aaSMatthew G. Knepley 
529b0a623aaSMatthew G. Knepley   Collective on DM
530b0a623aaSMatthew G. Knepley 
531b0a623aaSMatthew G. Knepley   Input Parameters:
532b0a623aaSMatthew G. Knepley + dm          - The DM
53324d039d7SMichael Lange . levels      - Number of overlap levels
534b0a623aaSMatthew G. Knepley . rootSection - The number of leaves for a given root point
535b0a623aaSMatthew G. Knepley . rootrank    - The rank of each edge into the root point
536b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point
537b0a623aaSMatthew G. Knepley - leafrank    - The rank of each process sharing a leaf point
538b0a623aaSMatthew G. Knepley 
539b0a623aaSMatthew G. Knepley   Output Parameters:
540a9f1d5b2SMichael Lange + ovLabel     - DMLabel containing remote overlap contributions as point/rank pairings
541b0a623aaSMatthew G. Knepley 
542b0a623aaSMatthew G. Knepley   Level: developer
543b0a623aaSMatthew G. Knepley 
5441fd9873aSMichael Lange .seealso: DMPlexDistributeOwnership(), DMPlexDistribute()
545b0a623aaSMatthew G. Knepley @*/
546a9f1d5b2SMichael Lange PetscErrorCode DMPlexCreateOverlap(DM dm, PetscInt levels, PetscSection rootSection, IS rootrank, PetscSection leafSection, IS leafrank, DMLabel *ovLabel)
547b0a623aaSMatthew G. Knepley {
548e540f424SMichael Lange   MPI_Comm           comm;
549b0a623aaSMatthew G. Knepley   DMLabel            ovAdjByRank; /* A DMLabel containing all points adjacent to shared points, separated by rank (value in label) */
550b0a623aaSMatthew G. Knepley   PetscSF            sfPoint, sfProc;
551b0a623aaSMatthew G. Knepley   const PetscSFNode *remote;
552b0a623aaSMatthew G. Knepley   const PetscInt    *local;
5531fd9873aSMichael Lange   const PetscInt    *nrank, *rrank;
554b0a623aaSMatthew G. Knepley   PetscInt          *adj = NULL;
5551fd9873aSMichael Lange   PetscInt           pStart, pEnd, p, sStart, sEnd, nleaves, l;
556b0a623aaSMatthew G. Knepley   PetscMPIInt        rank, numProcs;
55726a7d390SMatthew G. Knepley   PetscBool          useCone, useClosure, flg;
558b0a623aaSMatthew G. Knepley   PetscErrorCode     ierr;
559b0a623aaSMatthew G. Knepley 
560b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
561e540f424SMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
562e540f424SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
563e540f424SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
564b0a623aaSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
565b0a623aaSMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
566b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetChart(leafSection, &sStart, &sEnd);CHKERRQ(ierr);
567b0a623aaSMatthew G. Knepley   ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr);
568b0a623aaSMatthew G. Knepley   ierr = DMLabelCreate("Overlap adjacency", &ovAdjByRank);CHKERRQ(ierr);
569b0a623aaSMatthew G. Knepley   /* Handle leaves: shared with the root point */
570b0a623aaSMatthew G. Knepley   for (l = 0; l < nleaves; ++l) {
571b0a623aaSMatthew G. Knepley     PetscInt adjSize = PETSC_DETERMINE, a;
572b0a623aaSMatthew G. Knepley 
573b0a623aaSMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, local[l], &adjSize, &adj);CHKERRQ(ierr);
574b0a623aaSMatthew G. Knepley     for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remote[l].rank);CHKERRQ(ierr);}
575b0a623aaSMatthew G. Knepley   }
576b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(rootrank, &rrank);CHKERRQ(ierr);
577b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(leafrank, &nrank);CHKERRQ(ierr);
578b0a623aaSMatthew G. Knepley   /* Handle roots */
579b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
580b0a623aaSMatthew G. Knepley     PetscInt adjSize = PETSC_DETERMINE, neighbors = 0, noff, n, a;
581b0a623aaSMatthew G. Knepley 
582b0a623aaSMatthew G. Knepley     if ((p >= sStart) && (p < sEnd)) {
583b0a623aaSMatthew G. Knepley       /* Some leaves share a root with other leaves on different processes */
584b0a623aaSMatthew G. Knepley       ierr = PetscSectionGetDof(leafSection, p, &neighbors);CHKERRQ(ierr);
585b0a623aaSMatthew G. Knepley       if (neighbors) {
586b0a623aaSMatthew G. Knepley         ierr = PetscSectionGetOffset(leafSection, p, &noff);CHKERRQ(ierr);
587b0a623aaSMatthew G. Knepley         ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
588b0a623aaSMatthew G. Knepley         for (n = 0; n < neighbors; ++n) {
589b0a623aaSMatthew G. Knepley           const PetscInt remoteRank = nrank[noff+n];
590b0a623aaSMatthew G. Knepley 
591b0a623aaSMatthew G. Knepley           if (remoteRank == rank) continue;
592b0a623aaSMatthew G. Knepley           for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);}
593b0a623aaSMatthew G. Knepley         }
594b0a623aaSMatthew G. Knepley       }
595b0a623aaSMatthew G. Knepley     }
596b0a623aaSMatthew G. Knepley     /* Roots are shared with leaves */
597b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetDof(rootSection, p, &neighbors);CHKERRQ(ierr);
598b0a623aaSMatthew G. Knepley     if (!neighbors) continue;
599b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(rootSection, p, &noff);CHKERRQ(ierr);
600b0a623aaSMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
601b0a623aaSMatthew G. Knepley     for (n = 0; n < neighbors; ++n) {
602b0a623aaSMatthew G. Knepley       const PetscInt remoteRank = rrank[noff+n];
603b0a623aaSMatthew G. Knepley 
604b0a623aaSMatthew G. Knepley       if (remoteRank == rank) continue;
605b0a623aaSMatthew G. Knepley       for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);}
606b0a623aaSMatthew G. Knepley     }
607b0a623aaSMatthew G. Knepley   }
608b0a623aaSMatthew G. Knepley   ierr = PetscFree(adj);CHKERRQ(ierr);
609b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(rootrank, &rrank);CHKERRQ(ierr);
610b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(leafrank, &nrank);CHKERRQ(ierr);
61124d039d7SMichael Lange   /* Add additional overlap levels */
61224d039d7SMichael Lange   for (l = 1; l < levels; l++) {ierr = DMPlexPartitionLabelAdjacency(dm, ovAdjByRank);CHKERRQ(ierr);}
61326a7d390SMatthew G. Knepley   /* We require the closure in the overlap */
61426a7d390SMatthew G. Knepley   ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
61526a7d390SMatthew G. Knepley   ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr);
61626a7d390SMatthew G. Knepley   if (useCone || !useClosure) {
6175abbe4feSMichael Lange     ierr = DMPlexPartitionLabelClosure(dm, ovAdjByRank);CHKERRQ(ierr);
61826a7d390SMatthew G. Knepley   }
619e540f424SMichael Lange   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-overlap_view", &flg);CHKERRQ(ierr);
620e540f424SMichael Lange   if (flg) {
621b0a623aaSMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr);
622b0a623aaSMatthew G. Knepley     ierr = DMLabelView(ovAdjByRank, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
623b0a623aaSMatthew G. Knepley   }
6241fd9873aSMichael Lange   /* Make process SF and invert sender to receiver label */
625b0a623aaSMatthew G. Knepley   ierr = DMPlexCreateTwoSidedProcessSF(dm, sfPoint, rootSection, rootrank, leafSection, leafrank, NULL, &sfProc);CHKERRQ(ierr);
626a9f1d5b2SMichael Lange   ierr = DMLabelCreate("Overlap label", ovLabel);CHKERRQ(ierr);
627a9f1d5b2SMichael Lange   ierr = DMPlexPartitionLabelInvert(dm, ovAdjByRank, sfProc, *ovLabel);CHKERRQ(ierr);
628a9f1d5b2SMichael Lange   /* Add owned points, except for shared local points */
629a9f1d5b2SMichael Lange   for (p = pStart; p < pEnd; ++p) {ierr = DMLabelSetValue(*ovLabel, p, rank);CHKERRQ(ierr);}
630e540f424SMichael Lange   for (l = 0; l < nleaves; ++l) {
631a9f1d5b2SMichael Lange     ierr = DMLabelClearValue(*ovLabel, local[l], rank);CHKERRQ(ierr);
632a9f1d5b2SMichael Lange     ierr = DMLabelSetValue(*ovLabel, remote[l].index, remote[l].rank);CHKERRQ(ierr);
633e540f424SMichael Lange   }
634e540f424SMichael Lange   /* Clean up */
6351fd9873aSMichael Lange   ierr = DMLabelDestroy(&ovAdjByRank);CHKERRQ(ierr);
636b0a623aaSMatthew G. Knepley   ierr = PetscSFDestroy(&sfProc);CHKERRQ(ierr);
637b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
638b0a623aaSMatthew G. Knepley }
63970034214SMatthew G. Knepley 
64070034214SMatthew G. Knepley #undef __FUNCT__
64146f9b1c3SMichael Lange #define __FUNCT__ "DMPlexCreateOverlapMigrationSF"
64246f9b1c3SMichael Lange PetscErrorCode DMPlexCreateOverlapMigrationSF(DM dm, PetscSF overlapSF, PetscSF *migrationSF)
64346f9b1c3SMichael Lange {
64446f9b1c3SMichael Lange   MPI_Comm           comm;
64546f9b1c3SMichael Lange   PetscMPIInt        rank, numProcs;
64646f9b1c3SMichael Lange   PetscInt           d, dim, p, pStart, pEnd, nroots, nleaves, newLeaves, point, numSharedPoints;
64746f9b1c3SMichael Lange   PetscInt          *pointDepths, *remoteDepths, *ilocal;
64846f9b1c3SMichael Lange   PetscInt          *depthRecv, *depthShift, *depthIdx;
64946f9b1c3SMichael Lange   PetscSFNode       *iremote;
65046f9b1c3SMichael Lange   PetscSF            pointSF;
65146f9b1c3SMichael Lange   const PetscInt    *sharedLocal;
65246f9b1c3SMichael Lange   const PetscSFNode *overlapRemote, *sharedRemote;
65346f9b1c3SMichael Lange   PetscErrorCode     ierr;
65446f9b1c3SMichael Lange 
65546f9b1c3SMichael Lange   PetscFunctionBegin;
65646f9b1c3SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
65746f9b1c3SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
65846f9b1c3SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
65946f9b1c3SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
66046f9b1c3SMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
66146f9b1c3SMichael Lange 
66246f9b1c3SMichael Lange   /* Before building the migration SF we need to know the new stratum offsets */
66346f9b1c3SMichael Lange   ierr = PetscSFGetGraph(overlapSF, &nroots, &nleaves, NULL, &overlapRemote);CHKERRQ(ierr);
66446f9b1c3SMichael Lange   ierr = PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);CHKERRQ(ierr);
66546f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
66646f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
66746f9b1c3SMichael Lange     for (p=pStart; p<pEnd; p++) pointDepths[p] = d;
66846f9b1c3SMichael Lange   }
66946f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) remoteDepths[p] = -1;
67046f9b1c3SMichael Lange   ierr = PetscSFBcastBegin(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
67146f9b1c3SMichael Lange   ierr = PetscSFBcastEnd(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
67246f9b1c3SMichael Lange 
67346f9b1c3SMichael Lange   /* Count recevied points in each stratum and compute the internal strata shift */
67446f9b1c3SMichael Lange   ierr = PetscMalloc3(dim+1, &depthRecv, dim+1, &depthShift, dim+1, &depthIdx);CHKERRQ(ierr);
67546f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) depthRecv[d]=0;
67646f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) depthRecv[remoteDepths[p]]++;
67746f9b1c3SMichael Lange   depthShift[dim] = 0;
67846f9b1c3SMichael Lange   for (d=0; d<dim; d++) depthShift[d] = depthRecv[dim];
67946f9b1c3SMichael Lange   for (d=1; d<dim; d++) depthShift[d] += depthRecv[0];
68046f9b1c3SMichael Lange   for (d=dim-2; d>0; d--) depthShift[d] += depthRecv[d+1];
68146f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
68246f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
68346f9b1c3SMichael Lange     depthIdx[d] = pStart + depthShift[d];
68446f9b1c3SMichael Lange   }
68546f9b1c3SMichael Lange 
68646f9b1c3SMichael Lange   /* Form the overlap SF build an SF that describes the full overlap migration SF */
68746f9b1c3SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
68846f9b1c3SMichael Lange   newLeaves = pEnd - pStart + nleaves;
68909b7985cSMatthew G. Knepley   ierr = PetscMalloc1(newLeaves, &ilocal);CHKERRQ(ierr);
69009b7985cSMatthew G. Knepley   ierr = PetscMalloc1(newLeaves, &iremote);CHKERRQ(ierr);
69146f9b1c3SMichael Lange   /* First map local points to themselves */
69246f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
69346f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
69446f9b1c3SMichael Lange     for (p=pStart; p<pEnd; p++) {
69546f9b1c3SMichael Lange       point = p + depthShift[d];
69646f9b1c3SMichael Lange       ilocal[point] = point;
69746f9b1c3SMichael Lange       iremote[point].index = p;
69846f9b1c3SMichael Lange       iremote[point].rank = rank;
69946f9b1c3SMichael Lange       depthIdx[d]++;
70046f9b1c3SMichael Lange     }
70146f9b1c3SMichael Lange   }
70246f9b1c3SMichael Lange 
70346f9b1c3SMichael Lange   /* Add in the remote roots for currently shared points */
70446f9b1c3SMichael Lange   ierr = DMGetPointSF(dm, &pointSF);CHKERRQ(ierr);
70546f9b1c3SMichael Lange   ierr = PetscSFGetGraph(pointSF, NULL, &numSharedPoints, &sharedLocal, &sharedRemote);CHKERRQ(ierr);
70646f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
70746f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
70846f9b1c3SMichael Lange     for (p=0; p<numSharedPoints; p++) {
70946f9b1c3SMichael Lange       if (pStart <= sharedLocal[p] && sharedLocal[p] < pEnd) {
71046f9b1c3SMichael Lange         point = sharedLocal[p] + depthShift[d];
71146f9b1c3SMichael Lange         iremote[point].index = sharedRemote[p].index;
71246f9b1c3SMichael Lange         iremote[point].rank = sharedRemote[p].rank;
71346f9b1c3SMichael Lange       }
71446f9b1c3SMichael Lange     }
71546f9b1c3SMichael Lange   }
71646f9b1c3SMichael Lange 
71746f9b1c3SMichael Lange   /* Now add the incoming overlap points */
71846f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) {
71946f9b1c3SMichael Lange     point = depthIdx[remoteDepths[p]];
72046f9b1c3SMichael Lange     ilocal[point] = point;
72146f9b1c3SMichael Lange     iremote[point].index = overlapRemote[p].index;
72246f9b1c3SMichael Lange     iremote[point].rank = overlapRemote[p].rank;
72346f9b1c3SMichael Lange     depthIdx[remoteDepths[p]]++;
72446f9b1c3SMichael Lange   }
72515fff7beSMatthew G. Knepley   ierr = PetscFree2(pointDepths,remoteDepths);CHKERRQ(ierr);
72646f9b1c3SMichael Lange 
72746f9b1c3SMichael Lange   ierr = PetscSFCreate(comm, migrationSF);CHKERRQ(ierr);
72846f9b1c3SMichael Lange   ierr = PetscObjectSetName((PetscObject) *migrationSF, "Overlap Migration SF");CHKERRQ(ierr);
72946f9b1c3SMichael Lange   ierr = PetscSFSetFromOptions(*migrationSF);CHKERRQ(ierr);
73046f9b1c3SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
73146f9b1c3SMichael Lange   ierr = PetscSFSetGraph(*migrationSF, pEnd-pStart, newLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
73246f9b1c3SMichael Lange 
73346f9b1c3SMichael Lange   ierr = PetscFree3(depthRecv, depthShift, depthIdx);CHKERRQ(ierr);
73470034214SMatthew G. Knepley   PetscFunctionReturn(0);
73570034214SMatthew G. Knepley }
73670034214SMatthew G. Knepley 
73770034214SMatthew G. Knepley #undef __FUNCT__
738a9f1d5b2SMichael Lange #define __FUNCT__ "DMPlexStratifyMigrationSF"
739a9f1d5b2SMichael Lange /*@
740a9f1d5b2SMichael Lange   DMPlexStratifyMigrationSF - Add partition overlap to a distributed non-overlapping DM.
741a9f1d5b2SMichael Lange 
742a9f1d5b2SMichael Lange   Input Parameter:
743a9f1d5b2SMichael Lange + dm          - The DM
744a9f1d5b2SMichael Lange - sf          - A star forest with non-ordered leaves, usually defining a DM point migration
745a9f1d5b2SMichael Lange 
746a9f1d5b2SMichael Lange   Output Parameter:
747a9f1d5b2SMichael Lange . migrationSF - A star forest with added leaf indirection that ensures the resulting DM is stratified
748a9f1d5b2SMichael Lange 
749a9f1d5b2SMichael Lange   Level: developer
750a9f1d5b2SMichael Lange 
751a9f1d5b2SMichael Lange .seealso: DMPlexPartitionLabelCreateSF(), DMPlexDistribute(), DMPlexDistributeOverlap()
752a9f1d5b2SMichael Lange @*/
753a9f1d5b2SMichael Lange PetscErrorCode DMPlexStratifyMigrationSF(DM dm, PetscSF sf, PetscSF *migrationSF)
754a9f1d5b2SMichael Lange {
755a9f1d5b2SMichael Lange   MPI_Comm           comm;
756a9f1d5b2SMichael Lange   PetscMPIInt        rank, numProcs;
7577fab53ddSMatthew G. Knepley   PetscInt           d, ldepth, depth, p, pStart, pEnd, nroots, nleaves;
758a9f1d5b2SMichael Lange   PetscInt          *pointDepths, *remoteDepths, *ilocal;
759a9f1d5b2SMichael Lange   PetscInt          *depthRecv, *depthShift, *depthIdx;
760a9f1d5b2SMichael Lange   const PetscSFNode *iremote;
761a9f1d5b2SMichael Lange   PetscErrorCode     ierr;
762a9f1d5b2SMichael Lange 
763a9f1d5b2SMichael Lange   PetscFunctionBegin;
764a9f1d5b2SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
765a9f1d5b2SMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
766a9f1d5b2SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
767a9f1d5b2SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
7687fab53ddSMatthew G. Knepley   ierr = DMPlexGetDepth(dm, &ldepth);CHKERRQ(ierr);
7697fab53ddSMatthew G. Knepley   ierr = MPI_Allreduce(&ldepth, &depth, 1, MPIU_INT, MPI_MAX, comm);CHKERRQ(ierr);
7707fab53ddSMatthew G. Knepley   if ((ldepth >= 0) && (depth != ldepth)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Inconsistent Plex depth %d != %d", ldepth, depth);
771a9f1d5b2SMichael Lange 
772a9f1d5b2SMichael Lange   /* Before building the migration SF we need to know the new stratum offsets */
773a9f1d5b2SMichael Lange   ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, &iremote);CHKERRQ(ierr);
774a9f1d5b2SMichael Lange   ierr = PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);CHKERRQ(ierr);
7757fab53ddSMatthew G. Knepley   for (d = 0; d < depth+1; ++d) {
776a9f1d5b2SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
7777fab53ddSMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) pointDepths[p] = d;
778a9f1d5b2SMichael Lange   }
7797fab53ddSMatthew G. Knepley   for (p = 0; p < nleaves; ++p) remoteDepths[p] = -1;
780a9f1d5b2SMichael Lange   ierr = PetscSFBcastBegin(sf, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
781a9f1d5b2SMichael Lange   ierr = PetscSFBcastEnd(sf, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
782a9f1d5b2SMichael Lange   /* Count recevied points in each stratum and compute the internal strata shift */
7837fab53ddSMatthew G. Knepley   ierr = PetscMalloc3(depth+1, &depthRecv, depth+1, &depthShift, depth+1, &depthIdx);CHKERRQ(ierr);
7847fab53ddSMatthew G. Knepley   for (d = 0; d < depth+1; ++d) depthRecv[d] = 0;
7857fab53ddSMatthew G. Knepley   for (p = 0; p < nleaves; ++p) depthRecv[remoteDepths[p]]++;
7867fab53ddSMatthew G. Knepley   depthShift[depth] = 0;
7877fab53ddSMatthew G. Knepley   for (d = 0; d < depth; ++d) depthShift[d] = depthRecv[depth];
7887fab53ddSMatthew G. Knepley   for (d = 1; d < depth; ++d) depthShift[d] += depthRecv[0];
7897fab53ddSMatthew G. Knepley   for (d = depth-2; d > 0; --d) depthShift[d] += depthRecv[d+1];
7907fab53ddSMatthew G. Knepley   for (d = 0; d < depth+1; ++d) {depthIdx[d] = 0;}
791a9f1d5b2SMichael Lange   /* Derive a new local permutation based on stratified indices */
792a9f1d5b2SMichael Lange   ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr);
7937fab53ddSMatthew G. Knepley   for (p = 0; p < nleaves; ++p) {
7947fab53ddSMatthew G. Knepley     const PetscInt dep = remoteDepths[p];
7957fab53ddSMatthew G. Knepley 
7967fab53ddSMatthew G. Knepley     ilocal[p] = depthShift[dep] + depthIdx[dep];
7977fab53ddSMatthew G. Knepley     depthIdx[dep]++;
798a9f1d5b2SMichael Lange   }
799a9f1d5b2SMichael Lange   ierr = PetscSFCreate(comm, migrationSF);CHKERRQ(ierr);
800a9f1d5b2SMichael Lange   ierr = PetscObjectSetName((PetscObject) *migrationSF, "Migration SF");CHKERRQ(ierr);
801a9f1d5b2SMichael Lange   ierr = PetscSFSetGraph(*migrationSF, nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_COPY_VALUES);CHKERRQ(ierr);
802a9f1d5b2SMichael Lange   ierr = PetscFree2(pointDepths,remoteDepths);CHKERRQ(ierr);
803a9f1d5b2SMichael Lange   ierr = PetscFree3(depthRecv, depthShift, depthIdx);CHKERRQ(ierr);
804a9f1d5b2SMichael Lange   PetscFunctionReturn(0);
805a9f1d5b2SMichael Lange }
806a9f1d5b2SMichael Lange 
807a9f1d5b2SMichael Lange #undef __FUNCT__
80870034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeField"
80970034214SMatthew G. Knepley /*@
81070034214SMatthew G. Knepley   DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
81170034214SMatthew G. Knepley 
81270034214SMatthew G. Knepley   Collective on DM
81370034214SMatthew G. Knepley 
81470034214SMatthew G. Knepley   Input Parameters:
81570034214SMatthew G. Knepley + dm - The DMPlex object
81670034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
81770034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
81870034214SMatthew G. Knepley - originalVec - The existing data
81970034214SMatthew G. Knepley 
82070034214SMatthew G. Knepley   Output Parameters:
82170034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout
82270034214SMatthew G. Knepley - newVec - The new data
82370034214SMatthew G. Knepley 
82470034214SMatthew G. Knepley   Level: developer
82570034214SMatthew G. Knepley 
82670034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeFieldIS(), DMPlexDistributeData()
82770034214SMatthew G. Knepley @*/
82870034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec)
82970034214SMatthew G. Knepley {
83070034214SMatthew G. Knepley   PetscSF        fieldSF;
83170034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
83270034214SMatthew G. Knepley   PetscScalar   *originalValues, *newValues;
83370034214SMatthew G. Knepley   PetscErrorCode ierr;
83470034214SMatthew G. Knepley 
83570034214SMatthew G. Knepley   PetscFunctionBegin;
83670034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
83770034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
83870034214SMatthew G. Knepley 
83970034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
84070034214SMatthew G. Knepley   ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr);
84170034214SMatthew G. Knepley   ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr);
84270034214SMatthew G. Knepley 
84370034214SMatthew G. Knepley   ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr);
84470034214SMatthew G. Knepley   ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr);
84570034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);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);
892aaf8c182SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr);
893aaf8c182SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr);
89470034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
89570034214SMatthew G. Knepley   ierr = ISRestoreIndices(originalIS, &originalValues);CHKERRQ(ierr);
89670034214SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) pointSF), fieldSize, newValues, PETSC_OWN_POINTER, newIS);CHKERRQ(ierr);
89770034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
89870034214SMatthew G. Knepley   PetscFunctionReturn(0);
89970034214SMatthew G. Knepley }
90070034214SMatthew G. Knepley 
90170034214SMatthew G. Knepley #undef __FUNCT__
90270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeData"
90370034214SMatthew G. Knepley /*@
90470034214SMatthew G. Knepley   DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
90570034214SMatthew G. Knepley 
90670034214SMatthew G. Knepley   Collective on DM
90770034214SMatthew G. Knepley 
90870034214SMatthew G. Knepley   Input Parameters:
90970034214SMatthew G. Knepley + dm - The DMPlex object
91070034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
91170034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
91270034214SMatthew G. Knepley . datatype - The type of data
91370034214SMatthew G. Knepley - originalData - The existing data
91470034214SMatthew G. Knepley 
91570034214SMatthew G. Knepley   Output Parameters:
91670034214SMatthew G. Knepley + newSection - The PetscSection describing the new data layout
91770034214SMatthew G. Knepley - newData - The new data
91870034214SMatthew G. Knepley 
91970034214SMatthew G. Knepley   Level: developer
92070034214SMatthew G. Knepley 
92170034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField()
92270034214SMatthew G. Knepley @*/
92370034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData)
92470034214SMatthew G. Knepley {
92570034214SMatthew G. Knepley   PetscSF        fieldSF;
92670034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
92770034214SMatthew G. Knepley   PetscMPIInt    dataSize;
92870034214SMatthew G. Knepley   PetscErrorCode ierr;
92970034214SMatthew G. Knepley 
93070034214SMatthew G. Knepley   PetscFunctionBegin;
93170034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
93270034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
93370034214SMatthew G. Knepley 
93470034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
93570034214SMatthew G. Knepley   ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr);
93670034214SMatthew G. Knepley   ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr);
93770034214SMatthew G. Knepley 
93870034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
93970034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
94070034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
94170034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
94270034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
94370034214SMatthew G. Knepley   PetscFunctionReturn(0);
94470034214SMatthew G. Knepley }
94570034214SMatthew G. Knepley 
94670034214SMatthew G. Knepley #undef __FUNCT__
947cc71bff1SMichael Lange #define __FUNCT__ "DMPlexDistributeCones"
948ac04eaf7SToby Isaac PetscErrorCode DMPlexDistributeCones(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping original, ISLocalToGlobalMapping renumbering, DM dmParallel)
949cc71bff1SMichael Lange {
950f5bf2dbfSMichael Lange   DM_Plex               *mesh  = (DM_Plex*) dm->data;
951cc71bff1SMichael Lange   DM_Plex               *pmesh = (DM_Plex*) (dmParallel)->data;
952cc71bff1SMichael Lange   MPI_Comm               comm;
953cc71bff1SMichael Lange   PetscSF                coneSF;
954cc71bff1SMichael Lange   PetscSection           originalConeSection, newConeSection;
955ac04eaf7SToby Isaac   PetscInt              *remoteOffsets, *cones, *globCones, *newCones, newConesSize;
956cc71bff1SMichael Lange   PetscBool              flg;
957cc71bff1SMichael Lange   PetscErrorCode         ierr;
958cc71bff1SMichael Lange 
959cc71bff1SMichael Lange   PetscFunctionBegin;
960cc71bff1SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
961cc71bff1SMichael Lange   PetscValidPointer(dmParallel,4);
962cc71bff1SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
963cc71bff1SMichael Lange 
964cc71bff1SMichael Lange   /* Distribute cone section */
965cc71bff1SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
966cc71bff1SMichael Lange   ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr);
967cc71bff1SMichael Lange   ierr = DMPlexGetConeSection(dmParallel, &newConeSection);CHKERRQ(ierr);
968cc71bff1SMichael Lange   ierr = PetscSFDistributeSection(migrationSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr);
969cc71bff1SMichael Lange   ierr = DMSetUp(dmParallel);CHKERRQ(ierr);
970cc71bff1SMichael Lange   {
971cc71bff1SMichael Lange     PetscInt pStart, pEnd, p;
972cc71bff1SMichael Lange 
973cc71bff1SMichael Lange     ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr);
974cc71bff1SMichael Lange     for (p = pStart; p < pEnd; ++p) {
975cc71bff1SMichael Lange       PetscInt coneSize;
976cc71bff1SMichael Lange       ierr               = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr);
977cc71bff1SMichael Lange       pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize);
978cc71bff1SMichael Lange     }
979cc71bff1SMichael Lange   }
980cc71bff1SMichael Lange   /* Communicate and renumber cones */
981cc71bff1SMichael Lange   ierr = PetscSFCreateSectionSF(migrationSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr);
982cc71bff1SMichael Lange   ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
983ac04eaf7SToby Isaac   if (original) {
984ac04eaf7SToby Isaac     PetscInt numCones;
985ac04eaf7SToby Isaac 
986ac04eaf7SToby Isaac     ierr = PetscSectionGetStorageSize(originalConeSection,&numCones);CHKERRQ(ierr); ierr = PetscMalloc1(numCones,&globCones);CHKERRQ(ierr);
987ac04eaf7SToby Isaac     ierr = ISLocalToGlobalMappingApplyBlock(original, numCones, cones, globCones);CHKERRQ(ierr);
988ac04eaf7SToby Isaac   }
989ac04eaf7SToby Isaac   else {
990ac04eaf7SToby Isaac     globCones = cones;
991ac04eaf7SToby Isaac   }
992cc71bff1SMichael Lange   ierr = DMPlexGetCones(dmParallel, &newCones);CHKERRQ(ierr);
993ac04eaf7SToby Isaac   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, globCones, newCones);CHKERRQ(ierr);
994ac04eaf7SToby Isaac   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, globCones, newCones);CHKERRQ(ierr);
995ac04eaf7SToby Isaac   if (original) {
996ac04eaf7SToby Isaac     ierr = PetscFree(globCones);CHKERRQ(ierr);
997ac04eaf7SToby Isaac   }
998cc71bff1SMichael Lange   ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr);
999cc71bff1SMichael Lange   ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr);
10003533c52bSMatthew G. Knepley #if PETSC_USE_DEBUG
10013533c52bSMatthew G. Knepley   {
10023533c52bSMatthew G. Knepley     PetscInt  p;
10033533c52bSMatthew G. Knepley     PetscBool valid = PETSC_TRUE;
10043533c52bSMatthew G. Knepley     for (p = 0; p < newConesSize; ++p) {
10053533c52bSMatthew G. Knepley       if (newCones[p] < 0) {valid = PETSC_FALSE; ierr = PetscPrintf(PETSC_COMM_SELF, "Point %d not in overlap SF\n", p);CHKERRQ(ierr);}
10063533c52bSMatthew G. Knepley     }
10073533c52bSMatthew G. Knepley     if (!valid) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map");
10083533c52bSMatthew G. Knepley   }
10093533c52bSMatthew G. Knepley #endif
1010cc71bff1SMichael Lange   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr);
1011cc71bff1SMichael Lange   if (flg) {
1012cc71bff1SMichael Lange     ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr);
1013cc71bff1SMichael Lange     ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1014cc71bff1SMichael Lange     ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr);
1015cc71bff1SMichael Lange     ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1016cc71bff1SMichael Lange     ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr);
1017cc71bff1SMichael Lange   }
1018cc71bff1SMichael Lange   ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr);
1019cc71bff1SMichael Lange   ierr = DMPlexGetConeOrientations(dmParallel, &newCones);CHKERRQ(ierr);
1020cc71bff1SMichael Lange   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
1021cc71bff1SMichael Lange   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
1022cc71bff1SMichael Lange   ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr);
1023cc71bff1SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
1024cc71bff1SMichael Lange   /* Create supports and stratify sieve */
1025cc71bff1SMichael Lange   {
1026cc71bff1SMichael Lange     PetscInt pStart, pEnd;
1027cc71bff1SMichael Lange 
1028cc71bff1SMichael Lange     ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr);
1029cc71bff1SMichael Lange     ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr);
1030cc71bff1SMichael Lange   }
1031cc71bff1SMichael Lange   ierr = DMPlexSymmetrize(dmParallel);CHKERRQ(ierr);
1032cc71bff1SMichael Lange   ierr = DMPlexStratify(dmParallel);CHKERRQ(ierr);
1033f5bf2dbfSMichael Lange   pmesh->useCone    = mesh->useCone;
1034f5bf2dbfSMichael Lange   pmesh->useClosure = mesh->useClosure;
1035cc71bff1SMichael Lange   PetscFunctionReturn(0);
1036cc71bff1SMichael Lange }
1037cc71bff1SMichael Lange 
1038cc71bff1SMichael Lange #undef __FUNCT__
10390df0e737SMichael Lange #define __FUNCT__ "DMPlexDistributeCoordinates"
10400df0e737SMichael Lange PetscErrorCode DMPlexDistributeCoordinates(DM dm, PetscSF migrationSF, DM dmParallel)
10410df0e737SMichael Lange {
10420df0e737SMichael Lange   MPI_Comm         comm;
10430df0e737SMichael Lange   PetscSection     originalCoordSection, newCoordSection;
10440df0e737SMichael Lange   Vec              originalCoordinates, newCoordinates;
10450df0e737SMichael Lange   PetscInt         bs;
10460df0e737SMichael Lange   const char      *name;
10470df0e737SMichael Lange   const PetscReal *maxCell, *L;
10485dc8c3f7SMatthew G. Knepley   const DMBoundaryType *bd;
10490df0e737SMichael Lange   PetscErrorCode   ierr;
10500df0e737SMichael Lange 
10510df0e737SMichael Lange   PetscFunctionBegin;
10520df0e737SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10530df0e737SMichael Lange   PetscValidPointer(dmParallel, 3);
10540df0e737SMichael Lange 
10550df0e737SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
10560df0e737SMichael Lange   ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr);
10570df0e737SMichael Lange   ierr = DMGetCoordinateSection(dmParallel, &newCoordSection);CHKERRQ(ierr);
10580df0e737SMichael Lange   ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr);
10590df0e737SMichael Lange   if (originalCoordinates) {
10600df0e737SMichael Lange     ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr);
10610df0e737SMichael Lange     ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr);
10620df0e737SMichael Lange     ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr);
10630df0e737SMichael Lange 
10640df0e737SMichael Lange     ierr = DMPlexDistributeField(dm, migrationSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr);
10650df0e737SMichael Lange     ierr = DMSetCoordinatesLocal(dmParallel, newCoordinates);CHKERRQ(ierr);
10660df0e737SMichael Lange     ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr);
10670df0e737SMichael Lange     ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr);
10680df0e737SMichael Lange     ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr);
10690df0e737SMichael Lange   }
10705dc8c3f7SMatthew G. Knepley   ierr = DMGetPeriodicity(dm, &maxCell, &L, &bd);CHKERRQ(ierr);
10715dc8c3f7SMatthew G. Knepley   if (L) {ierr = DMSetPeriodicity(dmParallel, maxCell, L, bd);CHKERRQ(ierr);}
10720df0e737SMichael Lange   PetscFunctionReturn(0);
10730df0e737SMichael Lange }
10740df0e737SMichael Lange 
10750df0e737SMichael Lange #undef __FUNCT__
10760df0e737SMichael Lange #define __FUNCT__ "DMPlexDistributeLabels"
1077d995df53SMatthew G. Knepley /* Here we are assuming that process 0 always has everything */
10782eb1fa14SMichael Lange PetscErrorCode DMPlexDistributeLabels(DM dm, PetscSF migrationSF, DM dmParallel)
10790df0e737SMichael Lange {
10800df0e737SMichael Lange   MPI_Comm       comm;
10810df0e737SMichael Lange   PetscMPIInt    rank;
1082d995df53SMatthew G. Knepley   PetscInt       numLabels, numLocalLabels, l;
1083d995df53SMatthew G. Knepley   PetscBool      hasLabels = PETSC_FALSE;
10840df0e737SMichael Lange   PetscErrorCode ierr;
10850df0e737SMichael Lange 
10860df0e737SMichael Lange   PetscFunctionBegin;
10870df0e737SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10882eb1fa14SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
10890df0e737SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
10900df0e737SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
10910df0e737SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
10920df0e737SMichael Lange 
1093d995df53SMatthew G. Knepley   /* Everyone must have either the same number of labels, or none */
1094d995df53SMatthew G. Knepley   ierr = DMPlexGetNumLabels(dm, &numLocalLabels);CHKERRQ(ierr);
1095d995df53SMatthew G. Knepley   numLabels = numLocalLabels;
10960df0e737SMichael Lange   ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1097627847f0SMatthew G. Knepley   if (numLabels == numLocalLabels) hasLabels = PETSC_TRUE;
1098bdd2d751SMichael Lange   for (l = numLabels-1; l >= 0; --l) {
1099bdd2d751SMichael Lange     DMLabel     label = NULL, labelNew = NULL;
11000df0e737SMichael Lange     PetscBool   isdepth;
11010df0e737SMichael Lange 
1102d995df53SMatthew G. Knepley     if (hasLabels) {
1103bdd2d751SMichael Lange       ierr = DMPlexGetLabelByNum(dm, l, &label);CHKERRQ(ierr);
11040df0e737SMichael Lange       /* Skip "depth" because it is recreated */
1105bdd2d751SMichael Lange       ierr = PetscStrcmp(label->name, "depth", &isdepth);CHKERRQ(ierr);
1106bdd2d751SMichael Lange     }
11070df0e737SMichael Lange     ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
1108bdd2d751SMichael Lange     if (isdepth) continue;
1109bdd2d751SMichael Lange     ierr = DMLabelDistribute(label, migrationSF, &labelNew);CHKERRQ(ierr);
1110bdd2d751SMichael Lange     ierr = DMPlexAddLabel(dmParallel, labelNew);CHKERRQ(ierr);
11110df0e737SMichael Lange   }
11120df0e737SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
11130df0e737SMichael Lange   PetscFunctionReturn(0);
11140df0e737SMichael Lange }
11150df0e737SMichael Lange 
11169734c634SMichael Lange #undef __FUNCT__
11179734c634SMichael Lange #define __FUNCT__ "DMPlexDistributeSetupHybrid"
11189734c634SMichael Lange PetscErrorCode DMPlexDistributeSetupHybrid(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
11199734c634SMichael Lange {
11209734c634SMichael Lange   DM_Plex        *mesh  = (DM_Plex*) dm->data;
11219734c634SMichael Lange   DM_Plex        *pmesh = (DM_Plex*) (dmParallel)->data;
11229734c634SMichael Lange   MPI_Comm        comm;
11239734c634SMichael Lange   const PetscInt *gpoints;
11249734c634SMichael Lange   PetscInt        dim, depth, n, d;
11259734c634SMichael Lange   PetscErrorCode  ierr;
11269734c634SMichael Lange 
11279734c634SMichael Lange   PetscFunctionBegin;
11289734c634SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11299734c634SMichael Lange   PetscValidPointer(dmParallel, 4);
11309734c634SMichael Lange 
11319734c634SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
11329734c634SMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
11339734c634SMichael Lange 
11349734c634SMichael Lange   /* Setup hybrid structure */
11359734c634SMichael Lange   for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];}
11369734c634SMichael Lange   ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr);
11379734c634SMichael Lange   ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr);
11389734c634SMichael Lange   ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr);
11399734c634SMichael Lange   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
11409734c634SMichael Lange   for (d = 0; d <= dim; ++d) {
11419734c634SMichael Lange     PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p;
11429734c634SMichael Lange 
11439734c634SMichael Lange     if (pmax < 0) continue;
11449734c634SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr);
11459734c634SMichael Lange     ierr = DMPlexGetDepthStratum(dmParallel, d, NULL, &pEnd);CHKERRQ(ierr);
11469734c634SMichael Lange     ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr);
11479734c634SMichael Lange     for (p = 0; p < n; ++p) {
11489734c634SMichael Lange       const PetscInt point = gpoints[p];
11499734c634SMichael Lange 
11509734c634SMichael Lange       if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax;
11519734c634SMichael Lange     }
11529734c634SMichael Lange     if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax;
11539734c634SMichael Lange     else            pmesh->hybridPointMax[d] = -1;
11549734c634SMichael Lange   }
11559734c634SMichael Lange   ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr);
11569734c634SMichael Lange   PetscFunctionReturn(0);
11579734c634SMichael Lange }
11580df0e737SMichael Lange 
1159a6f36705SMichael Lange #undef __FUNCT__
1160a6f36705SMichael Lange #define __FUNCT__ "DMPlexDistributeSetupTree"
1161*08633170SToby Isaac PetscErrorCode DMPlexDistributeSetupTree(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping original, ISLocalToGlobalMapping renumbering, DM dmParallel)
1162a6f36705SMichael Lange {
116315078cd4SMichael Lange   DM_Plex        *mesh  = (DM_Plex*) dm->data;
116415078cd4SMichael Lange   DM_Plex        *pmesh = (DM_Plex*) (dmParallel)->data;
1165a6f36705SMichael Lange   MPI_Comm        comm;
1166a6f36705SMichael Lange   DM              refTree;
1167a6f36705SMichael Lange   PetscSection    origParentSection, newParentSection;
1168a6f36705SMichael Lange   PetscInt        *origParents, *origChildIDs;
1169a6f36705SMichael Lange   PetscBool       flg;
1170a6f36705SMichael Lange   PetscErrorCode  ierr;
1171a6f36705SMichael Lange 
1172a6f36705SMichael Lange   PetscFunctionBegin;
1173a6f36705SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1174a6f36705SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 4);
1175a6f36705SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
1176a6f36705SMichael Lange 
1177a6f36705SMichael Lange   /* Set up tree */
1178a6f36705SMichael Lange   ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr);
1179a6f36705SMichael Lange   ierr = DMPlexSetReferenceTree(dmParallel,refTree);CHKERRQ(ierr);
1180a6f36705SMichael Lange   ierr = DMPlexGetTree(dm,&origParentSection,&origParents,&origChildIDs,NULL,NULL);CHKERRQ(ierr);
1181a6f36705SMichael Lange   if (origParentSection) {
1182a6f36705SMichael Lange     PetscInt        pStart, pEnd;
1183*08633170SToby Isaac     PetscInt        *newParents, *newChildIDs, *globParents;
1184a6f36705SMichael Lange     PetscInt        *remoteOffsetsParents, newParentSize;
1185a6f36705SMichael Lange     PetscSF         parentSF;
1186a6f36705SMichael Lange 
1187a6f36705SMichael Lange     ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
1188a6f36705SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dmParallel),&newParentSection);CHKERRQ(ierr);
1189a6f36705SMichael Lange     ierr = PetscSectionSetChart(newParentSection,pStart,pEnd);CHKERRQ(ierr);
1190a6f36705SMichael Lange     ierr = PetscSFDistributeSection(migrationSF, origParentSection, &remoteOffsetsParents, newParentSection);CHKERRQ(ierr);
1191a6f36705SMichael Lange     ierr = PetscSFCreateSectionSF(migrationSF, origParentSection, remoteOffsetsParents, newParentSection, &parentSF);CHKERRQ(ierr);
1192a6f36705SMichael Lange     ierr = PetscSectionGetStorageSize(newParentSection,&newParentSize);CHKERRQ(ierr);
1193a6f36705SMichael Lange     ierr = PetscMalloc2(newParentSize,&newParents,newParentSize,&newChildIDs);CHKERRQ(ierr);
1194*08633170SToby Isaac     if (original) {
1195*08633170SToby Isaac       PetscInt numParents;
1196*08633170SToby Isaac 
1197*08633170SToby Isaac       ierr = PetscSectionGetStorageSize(origParentSection,&numParents);CHKERRQ(ierr);
1198*08633170SToby Isaac       ierr = PetscMalloc1(numParents,&globParents);CHKERRQ(ierr);
1199*08633170SToby Isaac       ierr = ISLocalToGlobalMappingApplyBlock(original, numParents, origParents, globParents);CHKERRQ(ierr);
1200*08633170SToby Isaac     }
1201*08633170SToby Isaac     else {
1202*08633170SToby Isaac       globParents = origParents;
1203*08633170SToby Isaac     }
1204*08633170SToby Isaac     ierr = PetscSFBcastBegin(parentSF, MPIU_INT, globParents, newParents);CHKERRQ(ierr);
1205*08633170SToby Isaac     ierr = PetscSFBcastEnd(parentSF, MPIU_INT, globParents, newParents);CHKERRQ(ierr);
1206*08633170SToby Isaac     if (original) {
1207*08633170SToby Isaac       ierr = PetscFree(globParents);CHKERRQ(ierr);
1208*08633170SToby Isaac     }
1209a6f36705SMichael Lange     ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1210a6f36705SMichael Lange     ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1211a6f36705SMichael Lange     ierr = ISGlobalToLocalMappingApplyBlock(renumbering,IS_GTOLM_MASK, newParentSize, newParents, NULL, newParents);CHKERRQ(ierr);
12124a54e071SToby Isaac #if PETSC_USE_DEBUG
12134a54e071SToby Isaac     {
12144a54e071SToby Isaac       PetscInt  p;
12154a54e071SToby Isaac       PetscBool valid = PETSC_TRUE;
12164a54e071SToby Isaac       for (p = 0; p < newParentSize; ++p) {
12174a54e071SToby Isaac         if (newParents[p] < 0) {valid = PETSC_FALSE; ierr = PetscPrintf(PETSC_COMM_SELF, "Point %d not in overlap SF\n", p);CHKERRQ(ierr);}
12184a54e071SToby Isaac       }
12194a54e071SToby Isaac       if (!valid) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map");
12204a54e071SToby Isaac     }
12214a54e071SToby Isaac #endif
1222a6f36705SMichael Lange     ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-parents_view", &flg);CHKERRQ(ierr);
1223a6f36705SMichael Lange     if (flg) {
1224a6f36705SMichael Lange       ierr = PetscPrintf(comm, "Serial Parent Section: \n");CHKERRQ(ierr);
1225a6f36705SMichael Lange       ierr = PetscSectionView(origParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1226a6f36705SMichael Lange       ierr = PetscPrintf(comm, "Parallel Parent Section: \n");CHKERRQ(ierr);
1227a6f36705SMichael Lange       ierr = PetscSectionView(newParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1228a6f36705SMichael Lange       ierr = PetscSFView(parentSF, NULL);CHKERRQ(ierr);
1229a6f36705SMichael Lange     }
1230a6f36705SMichael Lange     ierr = DMPlexSetTree(dmParallel,newParentSection,newParents,newChildIDs);CHKERRQ(ierr);
1231a6f36705SMichael Lange     ierr = PetscSectionDestroy(&newParentSection);CHKERRQ(ierr);
1232a6f36705SMichael Lange     ierr = PetscFree2(newParents,newChildIDs);CHKERRQ(ierr);
1233a6f36705SMichael Lange     ierr = PetscSFDestroy(&parentSF);CHKERRQ(ierr);
1234a6f36705SMichael Lange   }
123515078cd4SMichael Lange   pmesh->useAnchors = mesh->useAnchors;
1236a6f36705SMichael Lange   PetscFunctionReturn(0);
1237a6f36705SMichael Lange }
12380df0e737SMichael Lange 
12390df0e737SMichael Lange #undef __FUNCT__
12404eca1733SMichael Lange #define __FUNCT__ "DMPlexDistributeSF"
1241f8987ae8SMichael Lange PetscErrorCode DMPlexDistributeSF(DM dm, PetscSF migrationSF, DM dmParallel)
12424eca1733SMichael Lange {
12434eca1733SMichael Lange   DM_Plex               *mesh  = (DM_Plex*) dm->data;
12444eca1733SMichael Lange   DM_Plex               *pmesh = (DM_Plex*) (dmParallel)->data;
12454eca1733SMichael Lange   PetscMPIInt            rank, numProcs;
12464eca1733SMichael Lange   MPI_Comm               comm;
12474eca1733SMichael Lange   PetscErrorCode         ierr;
12484eca1733SMichael Lange 
12494eca1733SMichael Lange   PetscFunctionBegin;
12504eca1733SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12514eca1733SMichael Lange   PetscValidPointer(dmParallel,7);
12524eca1733SMichael Lange 
12534eca1733SMichael Lange   /* Create point SF for parallel mesh */
12544eca1733SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
12554eca1733SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
12564eca1733SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
12574eca1733SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
12584eca1733SMichael Lange   {
12594eca1733SMichael Lange     const PetscInt *leaves;
12604eca1733SMichael Lange     PetscSFNode    *remotePoints, *rowners, *lowners;
12614eca1733SMichael Lange     PetscInt        numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints;
12624eca1733SMichael Lange     PetscInt        pStart, pEnd;
12634eca1733SMichael Lange 
12644eca1733SMichael Lange     ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
12654eca1733SMichael Lange     ierr = PetscSFGetGraph(migrationSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr);
12664eca1733SMichael Lange     ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr);
12674eca1733SMichael Lange     for (p=0; p<numRoots; p++) {
12684eca1733SMichael Lange       rowners[p].rank  = -1;
12694eca1733SMichael Lange       rowners[p].index = -1;
12704eca1733SMichael Lange     }
12714eca1733SMichael Lange     ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12724eca1733SMichael Lange     ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12734eca1733SMichael Lange     for (p = 0; p < numLeaves; ++p) {
12744eca1733SMichael Lange       if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */
12754eca1733SMichael Lange         lowners[p].rank  = rank;
12764eca1733SMichael Lange         lowners[p].index = leaves ? leaves[p] : p;
12774eca1733SMichael Lange       } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */
12784eca1733SMichael Lange         lowners[p].rank  = -2;
12794eca1733SMichael Lange         lowners[p].index = -2;
12804eca1733SMichael Lange       }
12814eca1733SMichael Lange     }
12824eca1733SMichael Lange     for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */
12834eca1733SMichael Lange       rowners[p].rank  = -3;
12844eca1733SMichael Lange       rowners[p].index = -3;
12854eca1733SMichael Lange     }
12864eca1733SMichael Lange     ierr = PetscSFReduceBegin(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
12874eca1733SMichael Lange     ierr = PetscSFReduceEnd(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
12884eca1733SMichael Lange     ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12894eca1733SMichael Lange     ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12904eca1733SMichael Lange     for (p = 0; p < numLeaves; ++p) {
12914eca1733SMichael Lange       if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed");
12924eca1733SMichael Lange       if (lowners[p].rank != rank) ++numGhostPoints;
12934eca1733SMichael Lange     }
12944eca1733SMichael Lange     ierr = PetscMalloc1(numGhostPoints, &ghostPoints);CHKERRQ(ierr);
12954eca1733SMichael Lange     ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr);
12964eca1733SMichael Lange     for (p = 0, gp = 0; p < numLeaves; ++p) {
12974eca1733SMichael Lange       if (lowners[p].rank != rank) {
12984eca1733SMichael Lange         ghostPoints[gp]        = leaves ? leaves[p] : p;
12994eca1733SMichael Lange         remotePoints[gp].rank  = lowners[p].rank;
13004eca1733SMichael Lange         remotePoints[gp].index = lowners[p].index;
13014eca1733SMichael Lange         ++gp;
13024eca1733SMichael Lange       }
13034eca1733SMichael Lange     }
13044eca1733SMichael Lange     ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr);
13054eca1733SMichael Lange     ierr = PetscSFSetGraph((dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
13064eca1733SMichael Lange     ierr = PetscSFSetFromOptions((dmParallel)->sf);CHKERRQ(ierr);
13074eca1733SMichael Lange   }
13084eca1733SMichael Lange   pmesh->useCone    = mesh->useCone;
13094eca1733SMichael Lange   pmesh->useClosure = mesh->useClosure;
13104eca1733SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
13114eca1733SMichael Lange   PetscFunctionReturn(0);
13124eca1733SMichael Lange }
13134eca1733SMichael Lange 
13144eca1733SMichael Lange #undef __FUNCT__
1315f5bf2dbfSMichael Lange #define __FUNCT__ "DMPlexCreatePointSF"
1316f5bf2dbfSMichael Lange /*@C
1317f5bf2dbfSMichael Lange   DMPlexDerivePointSF - Build a point SF from an SF describing a point migration
1318f5bf2dbfSMichael Lange 
1319f5bf2dbfSMichael Lange   Input Parameter:
1320f5bf2dbfSMichael Lange + dm          - The source DMPlex object
1321f5bf2dbfSMichael Lange . migrationSF - The star forest that describes the parallel point remapping
13221627f6ccSMichael Lange . ownership   - Flag causing a vote to determine point ownership
1323f5bf2dbfSMichael Lange 
1324f5bf2dbfSMichael Lange   Output Parameter:
1325f5bf2dbfSMichael Lange - pointSF     - The star forest describing the point overlap in the remapped DM
1326f5bf2dbfSMichael Lange 
1327f5bf2dbfSMichael Lange   Level: developer
1328f5bf2dbfSMichael Lange 
1329f5bf2dbfSMichael Lange .seealso: DMPlexDistribute(), DMPlexDistributeOverlap()
1330f5bf2dbfSMichael Lange @*/
1331f5bf2dbfSMichael Lange PetscErrorCode DMPlexCreatePointSF(DM dm, PetscSF migrationSF, PetscBool ownership, PetscSF *pointSF)
1332f5bf2dbfSMichael Lange {
1333f5bf2dbfSMichael Lange   PetscMPIInt        rank;
13341627f6ccSMichael Lange   PetscInt           p, nroots, nleaves, idx, npointLeaves;
1335f5bf2dbfSMichael Lange   PetscInt          *pointLocal;
1336f5bf2dbfSMichael Lange   const PetscInt    *leaves;
1337f5bf2dbfSMichael Lange   const PetscSFNode *roots;
1338f5bf2dbfSMichael Lange   PetscSFNode       *rootNodes, *leafNodes, *pointRemote;
1339f5bf2dbfSMichael Lange   PetscErrorCode     ierr;
1340f5bf2dbfSMichael Lange 
1341f5bf2dbfSMichael Lange   PetscFunctionBegin;
1342f5bf2dbfSMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1343f5bf2dbfSMichael Lange   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
1344f5bf2dbfSMichael Lange 
1345f5bf2dbfSMichael Lange   ierr = PetscSFGetGraph(migrationSF, &nroots, &nleaves, &leaves, &roots);CHKERRQ(ierr);
1346f5bf2dbfSMichael Lange   ierr = PetscMalloc2(nroots, &rootNodes, nleaves, &leafNodes);CHKERRQ(ierr);
1347f5bf2dbfSMichael Lange   if (ownership) {
13481627f6ccSMichael Lange     /* Point ownership vote: Process with highest rank ownes shared points */
1349f5bf2dbfSMichael Lange     for (p = 0; p < nleaves; ++p) {
1350f5bf2dbfSMichael Lange       /* Either put in a bid or we know we own it */
1351f5bf2dbfSMichael Lange       leafNodes[p].rank  = rank;
135243f7d02bSMichael Lange       leafNodes[p].index = p;
1353f5bf2dbfSMichael Lange     }
1354f5bf2dbfSMichael Lange     for (p = 0; p < nroots; p++) {
13551627f6ccSMichael Lange       /* Root must not participate in the reduction, flag so that MAXLOC does not use */
1356f5bf2dbfSMichael Lange       rootNodes[p].rank  = -3;
1357f5bf2dbfSMichael Lange       rootNodes[p].index = -3;
1358f5bf2dbfSMichael Lange     }
1359f5bf2dbfSMichael Lange     ierr = PetscSFReduceBegin(migrationSF, MPIU_2INT, leafNodes, rootNodes, MPI_MAXLOC);CHKERRQ(ierr);
1360f5bf2dbfSMichael Lange     ierr = PetscSFReduceEnd(migrationSF, MPIU_2INT, leafNodes, rootNodes, MPI_MAXLOC);CHKERRQ(ierr);
1361f5bf2dbfSMichael Lange   } else {
1362f5bf2dbfSMichael Lange     for (p = 0; p < nroots; p++) {
1363f5bf2dbfSMichael Lange       rootNodes[p].index = -1;
1364f5bf2dbfSMichael Lange       rootNodes[p].rank = rank;
1365f5bf2dbfSMichael Lange     };
1366f5bf2dbfSMichael Lange     for (p = 0; p < nleaves; p++) {
1367f5bf2dbfSMichael Lange       /* Write new local id into old location */
1368f5bf2dbfSMichael Lange       if (roots[p].rank == rank) {
1369f5bf2dbfSMichael Lange         rootNodes[roots[p].index].index = leaves[p];
1370f5bf2dbfSMichael Lange       }
1371f5bf2dbfSMichael Lange     }
1372f5bf2dbfSMichael Lange   }
1373f5bf2dbfSMichael Lange   ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rootNodes, leafNodes);CHKERRQ(ierr);
1374f5bf2dbfSMichael Lange   ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rootNodes, leafNodes);CHKERRQ(ierr);
1375f5bf2dbfSMichael Lange 
1376f5bf2dbfSMichael Lange   for (npointLeaves = 0, p = 0; p < nleaves; p++) {if (leafNodes[p].rank != rank) npointLeaves++;}
13771627f6ccSMichael Lange   ierr = PetscMalloc1(npointLeaves, &pointLocal);CHKERRQ(ierr);
13781627f6ccSMichael Lange   ierr = PetscMalloc1(npointLeaves, &pointRemote);CHKERRQ(ierr);
1379f5bf2dbfSMichael Lange   for (idx = 0, p = 0; p < nleaves; p++) {
1380f5bf2dbfSMichael Lange     if (leafNodes[p].rank != rank) {
1381f5bf2dbfSMichael Lange       pointLocal[idx] = p;
1382f5bf2dbfSMichael Lange       pointRemote[idx] = leafNodes[p];
1383f5bf2dbfSMichael Lange       idx++;
1384f5bf2dbfSMichael Lange     }
1385f5bf2dbfSMichael Lange   }
1386f5bf2dbfSMichael Lange   ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), pointSF);CHKERRQ(ierr);
13871627f6ccSMichael Lange   ierr = PetscSFSetFromOptions(*pointSF);CHKERRQ(ierr);
1388f5bf2dbfSMichael Lange   ierr = PetscSFSetGraph(*pointSF, nleaves, npointLeaves, pointLocal, PETSC_OWN_POINTER, pointRemote, PETSC_OWN_POINTER);CHKERRQ(ierr);
1389f5bf2dbfSMichael Lange   ierr = PetscFree2(rootNodes, leafNodes);CHKERRQ(ierr);
1390f5bf2dbfSMichael Lange   PetscFunctionReturn(0);
1391f5bf2dbfSMichael Lange }
1392f5bf2dbfSMichael Lange 
1393f5bf2dbfSMichael Lange #undef __FUNCT__
139415078cd4SMichael Lange #define __FUNCT__ "DMPlexMigrate"
139515078cd4SMichael Lange /*@C
139615078cd4SMichael Lange   DMPlexMigrate  - Migrates internal DM data over the supplied star forest
139715078cd4SMichael Lange 
139815078cd4SMichael Lange   Input Parameter:
139915078cd4SMichael Lange + dm       - The source DMPlex object
14001627f6ccSMichael Lange . sf       - The star forest communication context describing the migration pattern
140115078cd4SMichael Lange 
140215078cd4SMichael Lange   Output Parameter:
140315078cd4SMichael Lange - targetDM - The target DMPlex object
140415078cd4SMichael Lange 
1405b9f40539SMichael Lange   Level: intermediate
140615078cd4SMichael Lange 
140715078cd4SMichael Lange .seealso: DMPlexDistribute(), DMPlexDistributeOverlap()
140815078cd4SMichael Lange @*/
1409b9f40539SMichael Lange PetscErrorCode DMPlexMigrate(DM dm, PetscSF sf, DM targetDM)
141015078cd4SMichael Lange {
1411b9f40539SMichael Lange   MPI_Comm               comm;
1412b9f40539SMichael Lange   PetscInt               dim, nroots;
1413b9f40539SMichael Lange   PetscSF                sfPoint;
141415078cd4SMichael Lange   ISLocalToGlobalMapping ltogMigration;
1415ac04eaf7SToby Isaac   ISLocalToGlobalMapping ltogOriginal = NULL;
1416b9f40539SMichael Lange   PetscBool              flg;
141715078cd4SMichael Lange   PetscErrorCode         ierr;
141815078cd4SMichael Lange 
141915078cd4SMichael Lange   PetscFunctionBegin;
142015078cd4SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14211b858b30SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_Migrate, dm, 0, 0, 0);CHKERRQ(ierr);
1422b9f40539SMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
142315078cd4SMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
142415078cd4SMichael Lange   ierr = DMSetDimension(targetDM, dim);CHKERRQ(ierr);
142515078cd4SMichael Lange 
1426bfb0467fSMichael Lange   /* Check for a one-to-all distribution pattern */
1427b9f40539SMichael Lange   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
1428b9f40539SMichael Lange   ierr = PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1429bfb0467fSMichael Lange   if (nroots >= 0) {
1430b9f40539SMichael Lange     IS                     isOriginal;
1431ac04eaf7SToby Isaac     PetscInt               n, size, nleaves;
1432ac04eaf7SToby Isaac     PetscInt              *numbering_orig, *numbering_new;
1433b9f40539SMichael Lange     /* Get the original point numbering */
1434b9f40539SMichael Lange     ierr = DMPlexCreatePointNumbering(dm, &isOriginal);CHKERRQ(ierr);
1435b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingCreateIS(isOriginal, &ltogOriginal);CHKERRQ(ierr);
1436b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingGetSize(ltogOriginal, &size);CHKERRQ(ierr);
1437b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingGetBlockIndices(ltogOriginal, (const PetscInt**)&numbering_orig);CHKERRQ(ierr);
1438b9f40539SMichael Lange     /* Convert to positive global numbers */
1439b9f40539SMichael Lange     for (n=0; n<size; n++) {if (numbering_orig[n] < 0) numbering_orig[n] = -(numbering_orig[n]+1);}
1440b9f40539SMichael Lange     /* Derive the new local-to-global mapping from the old one */
1441b9f40539SMichael Lange     ierr = PetscSFGetGraph(sf, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr);
1442b9f40539SMichael Lange     ierr = PetscMalloc1(nleaves, &numbering_new);CHKERRQ(ierr);
1443b9f40539SMichael Lange     ierr = PetscSFBcastBegin(sf, MPIU_INT, (PetscInt *) numbering_orig, numbering_new);CHKERRQ(ierr);
1444b9f40539SMichael Lange     ierr = PetscSFBcastEnd(sf, MPIU_INT, (PetscInt *) numbering_orig, numbering_new);CHKERRQ(ierr);
1445b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingCreate(comm, 1, nleaves, (const PetscInt*) numbering_new, PETSC_OWN_POINTER, &ltogMigration);CHKERRQ(ierr);
1446b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingRestoreIndices(ltogOriginal, (const PetscInt**)&numbering_orig);CHKERRQ(ierr);
1447b9f40539SMichael Lange     ierr = ISDestroy(&isOriginal);CHKERRQ(ierr);
144815078cd4SMichael Lange   } else {
1449bfb0467fSMichael Lange     /* One-to-all distribution pattern: We can derive LToG from SF */
145015078cd4SMichael Lange     ierr = ISLocalToGlobalMappingCreateSF(sf, 0, &ltogMigration);CHKERRQ(ierr);
145115078cd4SMichael Lange   }
1452b9f40539SMichael Lange   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
1453b9f40539SMichael Lange   if (flg) {
1454b9f40539SMichael Lange     ierr = PetscPrintf(comm, "Point renumbering for DM migration:\n");CHKERRQ(ierr);
1455b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingView(ltogMigration, NULL);CHKERRQ(ierr);
1456b9f40539SMichael Lange   }
145715078cd4SMichael Lange   /* Migrate DM data to target DM */
1458ac04eaf7SToby Isaac   ierr = DMPlexDistributeCones(dm, sf, ltogOriginal, ltogMigration, targetDM);CHKERRQ(ierr);
145915078cd4SMichael Lange   ierr = DMPlexDistributeCoordinates(dm, sf, targetDM);CHKERRQ(ierr);
146015078cd4SMichael Lange   ierr = DMPlexDistributeLabels(dm, sf, targetDM);CHKERRQ(ierr);
146115078cd4SMichael Lange   ierr = DMPlexDistributeSetupHybrid(dm, sf, ltogMigration, targetDM);CHKERRQ(ierr);
1462*08633170SToby Isaac   ierr = DMPlexDistributeSetupTree(dm, sf, ltogOriginal, ltogMigration, targetDM);CHKERRQ(ierr);
1463ac04eaf7SToby Isaac   ierr = ISLocalToGlobalMappingDestroy(&ltogOriginal);CHKERRQ(ierr);
1464302440fdSBarry Smith   ierr = ISLocalToGlobalMappingDestroy(&ltogMigration);CHKERRQ(ierr);
14651b858b30SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_Migrate, dm, 0, 0, 0);CHKERRQ(ierr);
146615078cd4SMichael Lange   PetscFunctionReturn(0);
146715078cd4SMichael Lange }
146815078cd4SMichael Lange 
146915078cd4SMichael Lange #undef __FUNCT__
147070034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistribute"
14713b8f15a2SMatthew G. Knepley /*@C
147270034214SMatthew G. Knepley   DMPlexDistribute - Distributes the mesh and any associated sections.
147370034214SMatthew G. Knepley 
147470034214SMatthew G. Knepley   Not Collective
147570034214SMatthew G. Knepley 
147670034214SMatthew G. Knepley   Input Parameter:
147770034214SMatthew G. Knepley + dm  - The original DMPlex object
147870034214SMatthew G. Knepley - overlap - The overlap of partitions, 0 is the default
147970034214SMatthew G. Knepley 
148070034214SMatthew G. Knepley   Output Parameter:
148170034214SMatthew G. Knepley + sf - The PetscSF used for point distribution
148270034214SMatthew G. Knepley - parallelMesh - The distributed DMPlex object, or NULL
148370034214SMatthew G. Knepley 
148470034214SMatthew G. Knepley   Note: If the mesh was not distributed, the return value is NULL.
148570034214SMatthew G. Knepley 
148670034214SMatthew G. Knepley   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
148770034214SMatthew G. Knepley   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
148870034214SMatthew G. Knepley   representation on the mesh.
148970034214SMatthew G. Knepley 
149070034214SMatthew G. Knepley   Level: intermediate
149170034214SMatthew G. Knepley 
149270034214SMatthew G. Knepley .keywords: mesh, elements
149370034214SMatthew G. Knepley .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
149470034214SMatthew G. Knepley @*/
149580cf41d5SMatthew G. Knepley PetscErrorCode DMPlexDistribute(DM dm, PetscInt overlap, PetscSF *sf, DM *dmParallel)
149670034214SMatthew G. Knepley {
149770034214SMatthew G. Knepley   MPI_Comm               comm;
149815078cd4SMichael Lange   PetscPartitioner       partitioner;
1499f8987ae8SMichael Lange   IS                     cellPart;
1500f8987ae8SMichael Lange   PetscSection           cellPartSection;
1501cf86098cSMatthew G. Knepley   DM                     dmCoord;
1502f8987ae8SMichael Lange   DMLabel                lblPartition, lblMigration;
150343f7d02bSMichael Lange   PetscSF                sfProcess, sfMigration, sfStratified, sfPoint;
150470034214SMatthew G. Knepley   PetscBool              flg;
150570034214SMatthew G. Knepley   PetscMPIInt            rank, numProcs, p;
150670034214SMatthew G. Knepley   PetscErrorCode         ierr;
150770034214SMatthew G. Knepley 
150870034214SMatthew G. Knepley   PetscFunctionBegin;
150970034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
151070034214SMatthew G. Knepley   if (sf) PetscValidPointer(sf,4);
151170034214SMatthew G. Knepley   PetscValidPointer(dmParallel,5);
151270034214SMatthew G. Knepley 
151370034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
151470034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
151570034214SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
151670034214SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
151770034214SMatthew G. Knepley 
151870034214SMatthew G. Knepley   *dmParallel = NULL;
151970034214SMatthew G. Knepley   if (numProcs == 1) PetscFunctionReturn(0);
152070034214SMatthew G. Knepley 
152115078cd4SMichael Lange   /* Create cell partition */
152277623264SMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
152377623264SMatthew G. Knepley   ierr = PetscSectionCreate(comm, &cellPartSection);CHKERRQ(ierr);
152415078cd4SMichael Lange   ierr = DMPlexGetPartitioner(dm, &partitioner);CHKERRQ(ierr);
152515078cd4SMichael Lange   ierr = PetscPartitionerPartition(partitioner, dm, cellPartSection, &cellPart);CHKERRQ(ierr);
1526f8987ae8SMichael Lange   {
1527f8987ae8SMichael Lange     /* Convert partition to DMLabel */
1528f8987ae8SMichael Lange     PetscInt proc, pStart, pEnd, npoints, poffset;
1529f8987ae8SMichael Lange     const PetscInt *points;
1530f8987ae8SMichael Lange     ierr = DMLabelCreate("Point Partition", &lblPartition);CHKERRQ(ierr);
1531f8987ae8SMichael Lange     ierr = ISGetIndices(cellPart, &points);CHKERRQ(ierr);
1532f8987ae8SMichael Lange     ierr = PetscSectionGetChart(cellPartSection, &pStart, &pEnd);CHKERRQ(ierr);
1533f8987ae8SMichael Lange     for (proc = pStart; proc < pEnd; proc++) {
1534f8987ae8SMichael Lange       ierr = PetscSectionGetDof(cellPartSection, proc, &npoints);CHKERRQ(ierr);
1535f8987ae8SMichael Lange       ierr = PetscSectionGetOffset(cellPartSection, proc, &poffset);CHKERRQ(ierr);
1536f8987ae8SMichael Lange       for (p = poffset; p < poffset+npoints; p++) {
1537f8987ae8SMichael Lange         ierr = DMLabelSetValue(lblPartition, points[p], proc);CHKERRQ(ierr);
153870034214SMatthew G. Knepley       }
1539f8987ae8SMichael Lange     }
1540f8987ae8SMichael Lange     ierr = ISRestoreIndices(cellPart, &points);CHKERRQ(ierr);
1541f8987ae8SMichael Lange   }
1542f8987ae8SMichael Lange   ierr = DMPlexPartitionLabelClosure(dm, lblPartition);CHKERRQ(ierr);
1543f8987ae8SMichael Lange   {
1544f8987ae8SMichael Lange     /* Build a global process SF */
1545f8987ae8SMichael Lange     PetscSFNode *remoteProc;
1546f8987ae8SMichael Lange     ierr = PetscMalloc1(numProcs, &remoteProc);CHKERRQ(ierr);
1547f8987ae8SMichael Lange     for (p = 0; p < numProcs; ++p) {
1548f8987ae8SMichael Lange       remoteProc[p].rank  = p;
1549f8987ae8SMichael Lange       remoteProc[p].index = rank;
1550f8987ae8SMichael Lange     }
1551f8987ae8SMichael Lange     ierr = PetscSFCreate(comm, &sfProcess);CHKERRQ(ierr);
1552f8987ae8SMichael Lange     ierr = PetscObjectSetName((PetscObject) sfProcess, "Process SF");CHKERRQ(ierr);
1553f8987ae8SMichael Lange     ierr = PetscSFSetGraph(sfProcess, numProcs, numProcs, NULL, PETSC_OWN_POINTER, remoteProc, PETSC_OWN_POINTER);CHKERRQ(ierr);
1554f8987ae8SMichael Lange   }
1555f8987ae8SMichael Lange   ierr = DMLabelCreate("Point migration", &lblMigration);CHKERRQ(ierr);
1556f8987ae8SMichael Lange   ierr = DMPlexPartitionLabelInvert(dm, lblPartition, sfProcess, lblMigration);CHKERRQ(ierr);
1557302440fdSBarry Smith   ierr = DMPlexPartitionLabelCreateSF(dm, lblMigration, &sfMigration);CHKERRQ(ierr);
155843f7d02bSMichael Lange   /* Stratify the SF in case we are migrating an already parallel plex */
155943f7d02bSMichael Lange   ierr = DMPlexStratifyMigrationSF(dm, sfMigration, &sfStratified);CHKERRQ(ierr);
156043f7d02bSMichael Lange   ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr);
156143f7d02bSMichael Lange   sfMigration = sfStratified;
1562f8987ae8SMichael Lange   ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
156370034214SMatthew G. Knepley   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
156470034214SMatthew G. Knepley   if (flg) {
1565f8987ae8SMichael Lange     ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr);
1566f8987ae8SMichael Lange     ierr = DMLabelView(lblPartition, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1567f8987ae8SMichael Lange     ierr = PetscSFView(sfMigration, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
156870034214SMatthew G. Knepley   }
1569f8987ae8SMichael Lange 
157015078cd4SMichael Lange   /* Create non-overlapping parallel DM and migrate internal data */
157170034214SMatthew G. Knepley   ierr = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr);
157270034214SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr);
1573b9f40539SMichael Lange   ierr = DMPlexMigrate(dm, sfMigration, *dmParallel);CHKERRQ(ierr);
157470034214SMatthew G. Knepley 
1575a157612eSMichael Lange   /* Build the point SF without overlap */
15761627f6ccSMichael Lange   ierr = DMPlexCreatePointSF(*dmParallel, sfMigration, PETSC_TRUE, &sfPoint);CHKERRQ(ierr);
1577f5bf2dbfSMichael Lange   ierr = DMSetPointSF(*dmParallel, sfPoint);CHKERRQ(ierr);
1578cf86098cSMatthew G. Knepley   ierr = DMGetCoordinateDM(*dmParallel, &dmCoord);CHKERRQ(ierr);
1579cf86098cSMatthew G. Knepley   if (dmCoord) {ierr = DMSetPointSF(dmCoord, sfPoint);CHKERRQ(ierr);}
1580f5bf2dbfSMichael Lange   if (flg) {ierr = PetscSFView(sfPoint, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);}
158170034214SMatthew G. Knepley 
1582a157612eSMichael Lange   if (overlap > 0) {
158315078cd4SMichael Lange     DM                 dmOverlap;
158415078cd4SMichael Lange     PetscInt           nroots, nleaves;
158515078cd4SMichael Lange     PetscSFNode       *newRemote;
158615078cd4SMichael Lange     const PetscSFNode *oldRemote;
158715078cd4SMichael Lange     PetscSF            sfOverlap, sfOverlapPoint;
1588a157612eSMichael Lange     /* Add the partition overlap to the distributed DM */
1589b9f40539SMichael Lange     ierr = DMPlexDistributeOverlap(*dmParallel, overlap, &sfOverlap, &dmOverlap);CHKERRQ(ierr);
1590a157612eSMichael Lange     ierr = DMDestroy(dmParallel);CHKERRQ(ierr);
1591a157612eSMichael Lange     *dmParallel = dmOverlap;
1592c389ea9fSToby Isaac     if (flg) {
15933d822a50SMichael Lange       ierr = PetscPrintf(comm, "Overlap Migration SF:\n");CHKERRQ(ierr);
159415078cd4SMichael Lange       ierr = PetscSFView(sfOverlap, NULL);CHKERRQ(ierr);
1595c389ea9fSToby Isaac     }
159643331d4aSMichael Lange 
1597f8987ae8SMichael Lange     /* Re-map the migration SF to establish the full migration pattern */
1598f8987ae8SMichael Lange     ierr = PetscSFGetGraph(sfMigration, &nroots, NULL, NULL, &oldRemote);CHKERRQ(ierr);
159915078cd4SMichael Lange     ierr = PetscSFGetGraph(sfOverlap, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr);
160043331d4aSMichael Lange     ierr = PetscMalloc1(nleaves, &newRemote);CHKERRQ(ierr);
160115078cd4SMichael Lange     ierr = PetscSFBcastBegin(sfOverlap, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr);
160215078cd4SMichael Lange     ierr = PetscSFBcastEnd(sfOverlap, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr);
160315078cd4SMichael Lange     ierr = PetscSFCreate(comm, &sfOverlapPoint);CHKERRQ(ierr);
160415078cd4SMichael Lange     ierr = PetscSFSetGraph(sfOverlapPoint, nroots, nleaves, NULL, PETSC_OWN_POINTER, newRemote, PETSC_OWN_POINTER);CHKERRQ(ierr);
160515078cd4SMichael Lange     ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr);
1606f8987ae8SMichael Lange     ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr);
160715078cd4SMichael Lange     sfMigration = sfOverlapPoint;
1608c389ea9fSToby Isaac   }
1609bf5244c7SMatthew G. Knepley   /* Cleanup Partition */
1610f8987ae8SMichael Lange   ierr = PetscSFDestroy(&sfProcess);CHKERRQ(ierr);
1611f8987ae8SMichael Lange   ierr = DMLabelDestroy(&lblPartition);CHKERRQ(ierr);
1612f8987ae8SMichael Lange   ierr = DMLabelDestroy(&lblMigration);CHKERRQ(ierr);
16134eca1733SMichael Lange   ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr);
16144eca1733SMichael Lange   ierr = ISDestroy(&cellPart);CHKERRQ(ierr);
1615721cbd76SMatthew G. Knepley   /* Copy BC */
1616721cbd76SMatthew G. Knepley   ierr = DMPlexCopyBoundary(dm, *dmParallel);CHKERRQ(ierr);
1617721cbd76SMatthew G. Knepley   /* Cleanup */
1618f8987ae8SMichael Lange   if (sf) {*sf = sfMigration;}
1619f8987ae8SMichael Lange   else    {ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr);}
1620f5bf2dbfSMichael Lange   ierr = PetscSFDestroy(&sfPoint);CHKERRQ(ierr);
162170034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
162270034214SMatthew G. Knepley   PetscFunctionReturn(0);
162370034214SMatthew G. Knepley }
1624a157612eSMichael Lange 
1625a157612eSMichael Lange #undef __FUNCT__
1626a157612eSMichael Lange #define __FUNCT__ "DMPlexDistributeOverlap"
1627a157612eSMichael Lange /*@C
1628a157612eSMichael Lange   DMPlexDistribute - Add partition overlap to a distributed non-overlapping DM.
1629a157612eSMichael Lange 
1630a157612eSMichael Lange   Not Collective
1631a157612eSMichael Lange 
1632a157612eSMichael Lange   Input Parameter:
1633a157612eSMichael Lange + dm  - The non-overlapping distrbuted DMPlex object
1634a157612eSMichael Lange - overlap - The overlap of partitions, 0 is the default
1635a157612eSMichael Lange 
1636a157612eSMichael Lange   Output Parameter:
1637a157612eSMichael Lange + sf - The PetscSF used for point distribution
1638a157612eSMichael Lange - dmOverlap - The overlapping distributed DMPlex object, or NULL
1639a157612eSMichael Lange 
1640a157612eSMichael Lange   Note: If the mesh was not distributed, the return value is NULL.
1641a157612eSMichael Lange 
1642a157612eSMichael Lange   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
1643a157612eSMichael Lange   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
1644a157612eSMichael Lange   representation on the mesh.
1645a157612eSMichael Lange 
1646a157612eSMichael Lange   Level: intermediate
1647a157612eSMichael Lange 
1648a157612eSMichael Lange .keywords: mesh, elements
1649a157612eSMichael Lange .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
1650a157612eSMichael Lange @*/
1651b9f40539SMichael Lange PetscErrorCode DMPlexDistributeOverlap(DM dm, PetscInt overlap, PetscSF *sf, DM *dmOverlap)
1652a157612eSMichael Lange {
1653a157612eSMichael Lange   MPI_Comm               comm;
1654a157612eSMichael Lange   PetscMPIInt            rank;
1655a157612eSMichael Lange   PetscSection           rootSection, leafSection;
1656a157612eSMichael Lange   IS                     rootrank, leafrank;
1657cf86098cSMatthew G. Knepley   DM                     dmCoord;
1658a9f1d5b2SMichael Lange   DMLabel                lblOverlap;
1659f5bf2dbfSMichael Lange   PetscSF                sfOverlap, sfStratified, sfPoint;
1660a157612eSMichael Lange   PetscErrorCode         ierr;
1661a157612eSMichael Lange 
1662a157612eSMichael Lange   PetscFunctionBegin;
1663a157612eSMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1664a157612eSMichael Lange   if (sf) PetscValidPointer(sf, 3);
1665a157612eSMichael Lange   PetscValidPointer(dmOverlap, 4);
1666a157612eSMichael Lange 
16671b858b30SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeOverlap, dm, 0, 0, 0);CHKERRQ(ierr);
1668a157612eSMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1669a157612eSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1670a157612eSMichael Lange 
1671a157612eSMichael Lange   /* Compute point overlap with neighbouring processes on the distributed DM */
1672a157612eSMichael Lange   ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1673a157612eSMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
1674a157612eSMichael Lange   ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr);
1675a157612eSMichael Lange   ierr = DMPlexDistributeOwnership(dm, rootSection, &rootrank, leafSection, &leafrank);CHKERRQ(ierr);
1676a9f1d5b2SMichael Lange   ierr = DMPlexCreateOverlap(dm, overlap, rootSection, rootrank, leafSection, leafrank, &lblOverlap);CHKERRQ(ierr);
1677a9f1d5b2SMichael Lange   /* Convert overlap label to stratified migration SF */
1678a9f1d5b2SMichael Lange   ierr = DMPlexPartitionLabelCreateSF(dm, lblOverlap, &sfOverlap);CHKERRQ(ierr);
1679a9f1d5b2SMichael Lange   ierr = DMPlexStratifyMigrationSF(dm, sfOverlap, &sfStratified);CHKERRQ(ierr);
1680a9f1d5b2SMichael Lange   ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr);
1681a9f1d5b2SMichael Lange   sfOverlap = sfStratified;
1682a9f1d5b2SMichael Lange   ierr = PetscObjectSetName((PetscObject) sfOverlap, "Overlap SF");CHKERRQ(ierr);
1683a9f1d5b2SMichael Lange   ierr = PetscSFSetFromOptions(sfOverlap);CHKERRQ(ierr);
1684a9f1d5b2SMichael Lange 
168515fff7beSMatthew G. Knepley   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
168615fff7beSMatthew G. Knepley   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
168715fff7beSMatthew G. Knepley   ierr = ISDestroy(&rootrank);CHKERRQ(ierr);
168815fff7beSMatthew G. Knepley   ierr = ISDestroy(&leafrank);CHKERRQ(ierr);
1689a157612eSMichael Lange   ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1690a157612eSMichael Lange 
1691a157612eSMichael Lange   /* Build the overlapping DM */
1692a157612eSMichael Lange   ierr = DMPlexCreate(comm, dmOverlap);CHKERRQ(ierr);
1693a157612eSMichael Lange   ierr = PetscObjectSetName((PetscObject) *dmOverlap, "Parallel Mesh");CHKERRQ(ierr);
1694a9f1d5b2SMichael Lange   ierr = DMPlexMigrate(dm, sfOverlap, *dmOverlap);CHKERRQ(ierr);
1695f5bf2dbfSMichael Lange   /* Build the new point SF */
16961627f6ccSMichael Lange   ierr = DMPlexCreatePointSF(*dmOverlap, sfOverlap, PETSC_FALSE, &sfPoint);CHKERRQ(ierr);
1697f5bf2dbfSMichael Lange   ierr = DMSetPointSF(*dmOverlap, sfPoint);CHKERRQ(ierr);
1698cf86098cSMatthew G. Knepley   ierr = DMGetCoordinateDM(*dmOverlap, &dmCoord);CHKERRQ(ierr);
1699cf86098cSMatthew G. Knepley   if (dmCoord) {ierr = DMSetPointSF(dmCoord, sfPoint);CHKERRQ(ierr);}
1700f5bf2dbfSMichael Lange   ierr = PetscSFDestroy(&sfPoint);CHKERRQ(ierr);
1701a157612eSMichael Lange   /* Cleanup overlap partition */
1702a9f1d5b2SMichael Lange   ierr = DMLabelDestroy(&lblOverlap);CHKERRQ(ierr);
1703a9f1d5b2SMichael Lange   if (sf) *sf = sfOverlap;
1704a9f1d5b2SMichael Lange   else    {ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr);}
17051b858b30SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeOverlap, dm, 0, 0, 0);CHKERRQ(ierr);
1706a157612eSMichael Lange   PetscFunctionReturn(0);
1707a157612eSMichael Lange }
1708