xref: /petsc/src/dm/impls/plex/plexdistribute.c (revision a17985de68acc6928b22292ec164853d2c38efbf)
170034214SMatthew G. Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
270034214SMatthew G. Knepley #include <petscsf.h>
370034214SMatthew G. Knepley 
470034214SMatthew G. Knepley #undef __FUNCT__
570034214SMatthew G. Knepley #define __FUNCT__ "DMPlexSetAdjacencyUseCone"
670034214SMatthew G. Knepley /*@
770034214SMatthew G. Knepley   DMPlexSetAdjacencyUseCone - Define adjacency in the mesh using either the cone or the support first
870034214SMatthew G. Knepley 
970034214SMatthew G. Knepley   Input Parameters:
1070034214SMatthew G. Knepley + dm      - The DM object
1170034214SMatthew G. Knepley - useCone - Flag to use the cone first
1270034214SMatthew G. Knepley 
1370034214SMatthew G. Knepley   Level: intermediate
1470034214SMatthew G. Knepley 
1570034214SMatthew G. Knepley   Notes:
1670034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE
1770034214SMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in star(cone(p)),    useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
1870034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
1970034214SMatthew G. Knepley 
2070034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
2170034214SMatthew G. Knepley @*/
2270034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseCone(DM dm, PetscBool useCone)
2370034214SMatthew G. Knepley {
2470034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
2570034214SMatthew G. Knepley 
2670034214SMatthew G. Knepley   PetscFunctionBegin;
2770034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2870034214SMatthew G. Knepley   mesh->useCone = useCone;
2970034214SMatthew G. Knepley   PetscFunctionReturn(0);
3070034214SMatthew G. Knepley }
3170034214SMatthew G. Knepley 
3270034214SMatthew G. Knepley #undef __FUNCT__
3370034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacencyUseCone"
3470034214SMatthew G. Knepley /*@
3570034214SMatthew G. Knepley   DMPlexGetAdjacencyUseCone - Query whether adjacency in the mesh uses the cone or the support first
3670034214SMatthew G. Knepley 
3770034214SMatthew G. Knepley   Input Parameter:
3870034214SMatthew G. Knepley . dm      - The DM object
3970034214SMatthew G. Knepley 
4070034214SMatthew G. Knepley   Output Parameter:
4170034214SMatthew G. Knepley . useCone - Flag to use the cone first
4270034214SMatthew G. Knepley 
4370034214SMatthew G. Knepley   Level: intermediate
4470034214SMatthew G. Knepley 
4570034214SMatthew G. Knepley   Notes:
4670034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE
4770034214SMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in star(cone(p)),    useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
4870034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
4970034214SMatthew G. Knepley 
5070034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
5170034214SMatthew G. Knepley @*/
5270034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseCone(DM dm, PetscBool *useCone)
5370034214SMatthew G. Knepley {
5470034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
5570034214SMatthew G. Knepley 
5670034214SMatthew G. Knepley   PetscFunctionBegin;
5770034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5870034214SMatthew G. Knepley   PetscValidIntPointer(useCone, 2);
5970034214SMatthew G. Knepley   *useCone = mesh->useCone;
6070034214SMatthew G. Knepley   PetscFunctionReturn(0);
6170034214SMatthew G. Knepley }
6270034214SMatthew G. Knepley 
6370034214SMatthew G. Knepley #undef __FUNCT__
6470034214SMatthew G. Knepley #define __FUNCT__ "DMPlexSetAdjacencyUseClosure"
6570034214SMatthew G. Knepley /*@
6670034214SMatthew G. Knepley   DMPlexSetAdjacencyUseClosure - Define adjacency in the mesh using the transitive closure
6770034214SMatthew G. Knepley 
6870034214SMatthew G. Knepley   Input Parameters:
6970034214SMatthew G. Knepley + dm      - The DM object
7070034214SMatthew G. Knepley - useClosure - Flag to use the closure
7170034214SMatthew G. Knepley 
7270034214SMatthew G. Knepley   Level: intermediate
7370034214SMatthew G. Knepley 
7470034214SMatthew G. Knepley   Notes:
7570034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE
7670034214SMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in star(cone(p)),    useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
7770034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
7870034214SMatthew G. Knepley 
7970034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
8070034214SMatthew G. Knepley @*/
8170034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseClosure(DM dm, PetscBool useClosure)
8270034214SMatthew G. Knepley {
8370034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
8470034214SMatthew G. Knepley 
8570034214SMatthew G. Knepley   PetscFunctionBegin;
8670034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8770034214SMatthew G. Knepley   mesh->useClosure = useClosure;
8870034214SMatthew G. Knepley   PetscFunctionReturn(0);
8970034214SMatthew G. Knepley }
9070034214SMatthew G. Knepley 
9170034214SMatthew G. Knepley #undef __FUNCT__
9270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacencyUseClosure"
9370034214SMatthew G. Knepley /*@
9470034214SMatthew G. Knepley   DMPlexGetAdjacencyUseClosure - Query whether adjacency in the mesh uses the transitive closure
9570034214SMatthew G. Knepley 
9670034214SMatthew G. Knepley   Input Parameter:
9770034214SMatthew G. Knepley . dm      - The DM object
9870034214SMatthew G. Knepley 
9970034214SMatthew G. Knepley   Output Parameter:
10070034214SMatthew G. Knepley . useClosure - Flag to use the closure
10170034214SMatthew G. Knepley 
10270034214SMatthew G. Knepley   Level: intermediate
10370034214SMatthew G. Knepley 
10470034214SMatthew G. Knepley   Notes:
10570034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE
10670034214SMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in star(cone(p)),    useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
10770034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
10870034214SMatthew G. Knepley 
10970034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
11070034214SMatthew G. Knepley @*/
11170034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseClosure(DM dm, PetscBool *useClosure)
11270034214SMatthew G. Knepley {
11370034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
11470034214SMatthew G. Knepley 
11570034214SMatthew G. Knepley   PetscFunctionBegin;
11670034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11770034214SMatthew G. Knepley   PetscValidIntPointer(useClosure, 2);
11870034214SMatthew G. Knepley   *useClosure = mesh->useClosure;
11970034214SMatthew G. Knepley   PetscFunctionReturn(0);
12070034214SMatthew G. Knepley }
12170034214SMatthew G. Knepley 
12270034214SMatthew G. Knepley #undef __FUNCT__
123*a17985deSToby Isaac #define __FUNCT__ "DMPlexSetAdjacencyUseAnchors"
1248b0b4c70SToby Isaac /*@
125*a17985deSToby Isaac   DMPlexSetAdjacencyUseAnchors - Define adjacency in the mesh using the point-to-point constraints.
1268b0b4c70SToby Isaac 
1278b0b4c70SToby Isaac   Input Parameters:
1288b0b4c70SToby Isaac + dm      - The DM object
1298b0b4c70SToby Isaac - useConstraints - Flag to use the constraints.  If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place.
1308b0b4c70SToby Isaac 
1318b0b4c70SToby Isaac   Level: intermediate
1328b0b4c70SToby Isaac 
133*a17985deSToby Isaac .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
1348b0b4c70SToby Isaac @*/
135*a17985deSToby Isaac PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM dm, PetscBool useConstraints)
1368b0b4c70SToby Isaac {
1378b0b4c70SToby Isaac   DM_Plex *mesh = (DM_Plex *) dm->data;
1388b0b4c70SToby Isaac 
1398b0b4c70SToby Isaac   PetscFunctionBegin;
1408b0b4c70SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1418b0b4c70SToby Isaac   mesh->useConstraints = useConstraints;
1428b0b4c70SToby Isaac   PetscFunctionReturn(0);
1438b0b4c70SToby Isaac }
1448b0b4c70SToby Isaac 
1458b0b4c70SToby Isaac #undef __FUNCT__
146*a17985deSToby Isaac #define __FUNCT__ "DMPlexGetAdjacencyUseAnchors"
1478b0b4c70SToby Isaac /*@
148*a17985deSToby Isaac   DMPlexGetAdjacencyUseAnchors - Query whether adjacency in the mesh uses the point-to-point constraints.
1498b0b4c70SToby Isaac 
1508b0b4c70SToby Isaac   Input Parameter:
1518b0b4c70SToby Isaac . dm      - The DM object
1528b0b4c70SToby Isaac 
1538b0b4c70SToby Isaac   Output Parameter:
1548b0b4c70SToby Isaac . useConstraints - Flag to use the closure.  If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place.
1558b0b4c70SToby Isaac 
1568b0b4c70SToby Isaac   Level: intermediate
1578b0b4c70SToby Isaac 
158*a17985deSToby Isaac .seealso: DMPlexSetAdjacencyUseAnchors(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
1598b0b4c70SToby Isaac @*/
160*a17985deSToby Isaac PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM dm, PetscBool *useConstraints)
1618b0b4c70SToby Isaac {
1628b0b4c70SToby Isaac   DM_Plex *mesh = (DM_Plex *) dm->data;
1638b0b4c70SToby Isaac 
1648b0b4c70SToby Isaac   PetscFunctionBegin;
1658b0b4c70SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1668b0b4c70SToby Isaac   PetscValidIntPointer(useConstraints, 2);
1678b0b4c70SToby Isaac   *useConstraints = mesh->useConstraints;
1688b0b4c70SToby Isaac   PetscFunctionReturn(0);
1698b0b4c70SToby Isaac }
1708b0b4c70SToby Isaac 
1718b0b4c70SToby Isaac #undef __FUNCT__
17270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Cone_Internal"
17370034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Cone_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
17470034214SMatthew G. Knepley {
17570034214SMatthew G. Knepley   const PetscInt *cone = NULL;
17670034214SMatthew G. Knepley   PetscInt        numAdj = 0, maxAdjSize = *adjSize, coneSize, c;
17770034214SMatthew G. Knepley   PetscErrorCode  ierr;
17870034214SMatthew G. Knepley 
17970034214SMatthew G. Knepley   PetscFunctionBeginHot;
18070034214SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
18170034214SMatthew G. Knepley   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
18270034214SMatthew G. Knepley   for (c = 0; c < coneSize; ++c) {
18370034214SMatthew G. Knepley     const PetscInt *support = NULL;
18470034214SMatthew G. Knepley     PetscInt        supportSize, s, q;
18570034214SMatthew G. Knepley 
18670034214SMatthew G. Knepley     ierr = DMPlexGetSupportSize(dm, cone[c], &supportSize);CHKERRQ(ierr);
18770034214SMatthew G. Knepley     ierr = DMPlexGetSupport(dm, cone[c], &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);
21070034214SMatthew G. Knepley   for (s = 0; s < supportSize; ++s) {
21170034214SMatthew G. Knepley     const PetscInt *cone = NULL;
21270034214SMatthew G. Knepley     PetscInt        coneSize, c, q;
21370034214SMatthew G. Knepley 
21470034214SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, support[s], &coneSize);CHKERRQ(ierr);
21570034214SMatthew G. Knepley     ierr = DMPlexGetCone(dm, support[s], &cone);CHKERRQ(ierr);
21670034214SMatthew G. Knepley     for (c = 0; c < coneSize; ++c) {
21770034214SMatthew G. Knepley       for (q = 0; q < numAdj || (adj[numAdj++] = cone[c],0); ++q) {
21870034214SMatthew G. Knepley         if (cone[c] == adj[q]) break;
21970034214SMatthew G. Knepley       }
22070034214SMatthew G. Knepley       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
22170034214SMatthew G. Knepley     }
22270034214SMatthew G. Knepley   }
22370034214SMatthew G. Knepley   *adjSize = numAdj;
22470034214SMatthew G. Knepley   PetscFunctionReturn(0);
22570034214SMatthew G. Knepley }
22670034214SMatthew G. Knepley 
22770034214SMatthew G. Knepley #undef __FUNCT__
22870034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Transitive_Internal"
22970034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Transitive_Internal(DM dm, PetscInt p, PetscBool useClosure, PetscInt *adjSize, PetscInt adj[])
23070034214SMatthew G. Knepley {
23170034214SMatthew G. Knepley   PetscInt      *star = NULL;
23270034214SMatthew G. Knepley   PetscInt       numAdj = 0, maxAdjSize = *adjSize, starSize, s;
23370034214SMatthew G. Knepley   PetscErrorCode ierr;
23470034214SMatthew G. Knepley 
23570034214SMatthew G. Knepley   PetscFunctionBeginHot;
23670034214SMatthew G. Knepley   ierr = DMPlexGetTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr);
23770034214SMatthew G. Knepley   for (s = 0; s < starSize*2; s += 2) {
23870034214SMatthew G. Knepley     const PetscInt *closure = NULL;
23970034214SMatthew G. Knepley     PetscInt        closureSize, c, q;
24070034214SMatthew G. Knepley 
24170034214SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
24270034214SMatthew G. Knepley     for (c = 0; c < closureSize*2; c += 2) {
24370034214SMatthew G. Knepley       for (q = 0; q < numAdj || (adj[numAdj++] = closure[c],0); ++q) {
24470034214SMatthew G. Knepley         if (closure[c] == adj[q]) break;
24570034214SMatthew G. Knepley       }
24670034214SMatthew G. Knepley       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
24770034214SMatthew G. Knepley     }
24870034214SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
24970034214SMatthew G. Knepley   }
25070034214SMatthew G. Knepley   ierr = DMPlexRestoreTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr);
25170034214SMatthew G. Knepley   *adjSize = numAdj;
25270034214SMatthew G. Knepley   PetscFunctionReturn(0);
25370034214SMatthew G. Knepley }
25470034214SMatthew G. Knepley 
25570034214SMatthew G. Knepley #undef __FUNCT__
25670034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Internal"
2578b0b4c70SToby Isaac PetscErrorCode DMPlexGetAdjacency_Internal(DM dm, PetscInt p, PetscBool useCone, PetscBool useTransitiveClosure, PetscBool useConstraints, PetscInt *adjSize, PetscInt *adj[])
25870034214SMatthew G. Knepley {
25979bad331SMatthew G. Knepley   static PetscInt asiz = 0;
2608b0b4c70SToby Isaac   PetscInt maxAnchors = 1;
2618b0b4c70SToby Isaac   PetscInt aStart = -1, aEnd = -1;
2628b0b4c70SToby Isaac   PetscInt maxAdjSize;
2638b0b4c70SToby Isaac   PetscSection aSec = NULL;
2648b0b4c70SToby Isaac   IS aIS = NULL;
2658b0b4c70SToby Isaac   const PetscInt *anchors;
26670034214SMatthew G. Knepley   PetscErrorCode  ierr;
26770034214SMatthew G. Knepley 
26870034214SMatthew G. Knepley   PetscFunctionBeginHot;
2698b0b4c70SToby Isaac   if (useConstraints) {
270*a17985deSToby Isaac     ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr);
2718b0b4c70SToby Isaac     if (aSec) {
2728b0b4c70SToby Isaac       ierr = PetscSectionGetMaxDof(aSec,&maxAnchors);CHKERRQ(ierr);
27324c766afSToby Isaac       maxAnchors = PetscMax(1,maxAnchors);
2748b0b4c70SToby Isaac       ierr = PetscSectionGetChart(aSec,&aStart,&aEnd);CHKERRQ(ierr);
2758b0b4c70SToby Isaac       ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr);
2768b0b4c70SToby Isaac     }
2778b0b4c70SToby Isaac   }
27879bad331SMatthew G. Knepley   if (!*adj) {
27924c766afSToby Isaac     PetscInt depth, coneSeries, supportSeries, maxC, maxS, pStart, pEnd;
28079bad331SMatthew G. Knepley 
28124c766afSToby Isaac     ierr  = DMPlexGetChart(dm, &pStart,&pEnd);CHKERRQ(ierr);
28279bad331SMatthew G. Knepley     ierr  = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
28324c766afSToby Isaac     ierr  = DMPlexGetMaxSizes(dm, &maxC, &maxS);CHKERRQ(ierr);
28424c766afSToby Isaac     coneSeries    = (maxC > 1) ? ((PetscPowInt(maxC,depth+1)-1)/(maxC-1)) : depth+1;
28524c766afSToby Isaac     supportSeries = (maxS > 1) ? ((PetscPowInt(maxS,depth+1)-1)/(maxS-1)) : depth+1;
28624c766afSToby Isaac     asiz  = PetscMax(PetscPowInt(maxS,depth)*coneSeries,PetscPowInt(maxC,depth)*supportSeries);
2878b0b4c70SToby Isaac     asiz *= maxAnchors;
28824c766afSToby Isaac     asiz  = PetscMin(asiz,pEnd-pStart);
28979bad331SMatthew G. Knepley     ierr  = PetscMalloc1(asiz,adj);CHKERRQ(ierr);
29079bad331SMatthew G. Knepley   }
29179bad331SMatthew G. Knepley   if (*adjSize < 0) *adjSize = asiz;
2928b0b4c70SToby Isaac   maxAdjSize = *adjSize;
29370034214SMatthew G. Knepley   if (useTransitiveClosure) {
29479bad331SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Transitive_Internal(dm, p, useCone, adjSize, *adj);CHKERRQ(ierr);
29570034214SMatthew G. Knepley   } else if (useCone) {
29679bad331SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Cone_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
29770034214SMatthew G. Knepley   } else {
29879bad331SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Support_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
29970034214SMatthew G. Knepley   }
3008b0b4c70SToby Isaac   if (useConstraints && aSec) {
3018b0b4c70SToby Isaac     PetscInt origSize = *adjSize;
3028b0b4c70SToby Isaac     PetscInt numAdj = origSize;
3038b0b4c70SToby Isaac     PetscInt i = 0, j;
3048b0b4c70SToby Isaac     PetscInt *orig = *adj;
3058b0b4c70SToby Isaac 
3068b0b4c70SToby Isaac     while (i < origSize) {
3078b0b4c70SToby Isaac       PetscInt p = orig[i];
3088b0b4c70SToby Isaac       PetscInt aDof = 0;
3098b0b4c70SToby Isaac 
3108b0b4c70SToby Isaac       if (p >= aStart && p < aEnd) {
3118b0b4c70SToby Isaac         ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr);
3128b0b4c70SToby Isaac       }
3138b0b4c70SToby Isaac       if (aDof) {
3148b0b4c70SToby Isaac         PetscInt aOff;
3158b0b4c70SToby Isaac         PetscInt s, q;
3168b0b4c70SToby Isaac 
3178b0b4c70SToby Isaac         for (j = i + 1; j < numAdj; j++) {
3188b0b4c70SToby Isaac           orig[j - 1] = orig[j];
3198b0b4c70SToby Isaac         }
3208b0b4c70SToby Isaac         origSize--;
3218b0b4c70SToby Isaac         numAdj--;
3228b0b4c70SToby Isaac         ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr);
3238b0b4c70SToby Isaac         for (s = 0; s < aDof; ++s) {
3248b0b4c70SToby Isaac           for (q = 0; q < numAdj || (orig[numAdj++] = anchors[aOff+s],0); ++q) {
3258b0b4c70SToby Isaac             if (anchors[aOff+s] == orig[q]) break;
3268b0b4c70SToby Isaac           }
3278b0b4c70SToby Isaac           if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
3288b0b4c70SToby Isaac         }
3298b0b4c70SToby Isaac       }
3308b0b4c70SToby Isaac       else {
3318b0b4c70SToby Isaac         i++;
3328b0b4c70SToby Isaac       }
3338b0b4c70SToby Isaac     }
3348b0b4c70SToby Isaac     *adjSize = numAdj;
3358b0b4c70SToby Isaac     ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr);
3368b0b4c70SToby Isaac   }
33770034214SMatthew G. Knepley   PetscFunctionReturn(0);
33870034214SMatthew G. Knepley }
33970034214SMatthew G. Knepley 
34070034214SMatthew G. Knepley #undef __FUNCT__
34170034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency"
34270034214SMatthew G. Knepley /*@
34370034214SMatthew G. Knepley   DMPlexGetAdjacency - Return all points adjacent to the given point
34470034214SMatthew G. Knepley 
34570034214SMatthew G. Knepley   Input Parameters:
34670034214SMatthew G. Knepley + dm - The DM object
34770034214SMatthew G. Knepley . p  - The point
34870034214SMatthew G. Knepley . adjSize - The maximum size of adj if it is non-NULL, or PETSC_DETERMINE
34970034214SMatthew G. Knepley - adj - Either NULL so that the array is allocated, or an existing array with size adjSize
35070034214SMatthew G. Knepley 
35170034214SMatthew G. Knepley   Output Parameters:
35270034214SMatthew G. Knepley + adjSize - The number of adjacent points
35370034214SMatthew G. Knepley - adj - The adjacent points
35470034214SMatthew G. Knepley 
35570034214SMatthew G. Knepley   Level: advanced
35670034214SMatthew G. Knepley 
35770034214SMatthew G. Knepley   Notes: The user must PetscFree the adj array if it was not passed in.
35870034214SMatthew G. Knepley 
35970034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexDistribute(), DMCreateMatrix(), DMPlexPreallocateOperator()
36070034214SMatthew G. Knepley @*/
36170034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacency(DM dm, PetscInt p, PetscInt *adjSize, PetscInt *adj[])
36270034214SMatthew G. Knepley {
36370034214SMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
36470034214SMatthew G. Knepley   PetscErrorCode ierr;
36570034214SMatthew G. Knepley 
36670034214SMatthew G. Knepley   PetscFunctionBeginHot;
36770034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
36870034214SMatthew G. Knepley   PetscValidPointer(adjSize,3);
36970034214SMatthew G. Knepley   PetscValidPointer(adj,4);
3708b0b4c70SToby Isaac   ierr = DMPlexGetAdjacency_Internal(dm, p, mesh->useCone, mesh->useClosure, mesh->useConstraints, adjSize, adj);CHKERRQ(ierr);
37170034214SMatthew G. Knepley   PetscFunctionReturn(0);
37270034214SMatthew G. Knepley }
37370034214SMatthew G. Knepley 
37470034214SMatthew G. Knepley #undef __FUNCT__
37570034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeField"
37670034214SMatthew G. Knepley /*@
37770034214SMatthew G. Knepley   DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
37870034214SMatthew G. Knepley 
37970034214SMatthew G. Knepley   Collective on DM
38070034214SMatthew G. Knepley 
38170034214SMatthew G. Knepley   Input Parameters:
38270034214SMatthew G. Knepley + dm - The DMPlex object
38370034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
38470034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
38570034214SMatthew G. Knepley - originalVec - The existing data
38670034214SMatthew G. Knepley 
38770034214SMatthew G. Knepley   Output Parameters:
38870034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout
38970034214SMatthew G. Knepley - newVec - The new data
39070034214SMatthew G. Knepley 
39170034214SMatthew G. Knepley   Level: developer
39270034214SMatthew G. Knepley 
39370034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeData()
39470034214SMatthew G. Knepley @*/
39570034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec)
39670034214SMatthew G. Knepley {
39770034214SMatthew G. Knepley   PetscSF        fieldSF;
39870034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
39970034214SMatthew G. Knepley   PetscScalar   *originalValues, *newValues;
40070034214SMatthew G. Knepley   PetscErrorCode ierr;
40170034214SMatthew G. Knepley 
40270034214SMatthew G. Knepley   PetscFunctionBegin;
40370034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
40470034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
40570034214SMatthew G. Knepley 
40670034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
40770034214SMatthew G. Knepley   ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr);
40870034214SMatthew G. Knepley   ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr);
40970034214SMatthew G. Knepley 
41070034214SMatthew G. Knepley   ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr);
41170034214SMatthew G. Knepley   ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr);
41270034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
41370034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
41470034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
41570034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
41670034214SMatthew G. Knepley   ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr);
41770034214SMatthew G. Knepley   ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr);
41870034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
41970034214SMatthew G. Knepley   PetscFunctionReturn(0);
42070034214SMatthew G. Knepley }
42170034214SMatthew G. Knepley 
42270034214SMatthew G. Knepley #undef __FUNCT__
42370034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeData"
42470034214SMatthew G. Knepley /*@
42570034214SMatthew G. Knepley   DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
42670034214SMatthew G. Knepley 
42770034214SMatthew G. Knepley   Collective on DM
42870034214SMatthew G. Knepley 
42970034214SMatthew G. Knepley   Input Parameters:
43070034214SMatthew G. Knepley + dm - The DMPlex object
43170034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
43270034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
43370034214SMatthew G. Knepley . datatype - The type of data
43470034214SMatthew G. Knepley - originalData - The existing data
43570034214SMatthew G. Knepley 
43670034214SMatthew G. Knepley   Output Parameters:
43770034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout
43870034214SMatthew G. Knepley - newData - The new data
43970034214SMatthew G. Knepley 
44070034214SMatthew G. Knepley   Level: developer
44170034214SMatthew G. Knepley 
44270034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField()
44370034214SMatthew G. Knepley @*/
44470034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData)
44570034214SMatthew G. Knepley {
44670034214SMatthew G. Knepley   PetscSF        fieldSF;
44770034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
44870034214SMatthew G. Knepley   PetscMPIInt    dataSize;
44970034214SMatthew G. Knepley   PetscErrorCode ierr;
45070034214SMatthew G. Knepley 
45170034214SMatthew G. Knepley   PetscFunctionBegin;
45270034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
45370034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
45470034214SMatthew G. Knepley 
45570034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
45670034214SMatthew G. Knepley   ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr);
45770034214SMatthew G. Knepley   ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr);
45870034214SMatthew G. Knepley 
45970034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
46070034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
46170034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
46270034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
46370034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
46470034214SMatthew G. Knepley   PetscFunctionReturn(0);
46570034214SMatthew G. Knepley }
46670034214SMatthew G. Knepley 
46770034214SMatthew G. Knepley #undef __FUNCT__
46870034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistribute"
46970034214SMatthew G. Knepley /*@C
47070034214SMatthew G. Knepley   DMPlexDistribute - Distributes the mesh and any associated sections.
47170034214SMatthew G. Knepley 
47270034214SMatthew G. Knepley   Not Collective
47370034214SMatthew G. Knepley 
47470034214SMatthew G. Knepley   Input Parameter:
47570034214SMatthew G. Knepley + dm  - The original DMPlex object
47670034214SMatthew G. Knepley . partitioner - The partitioning package, or NULL for the default
47770034214SMatthew G. Knepley - overlap - The overlap of partitions, 0 is the default
47870034214SMatthew G. Knepley 
47970034214SMatthew G. Knepley   Output Parameter:
48070034214SMatthew G. Knepley + sf - The PetscSF used for point distribution
48170034214SMatthew G. Knepley - parallelMesh - The distributed DMPlex object, or NULL
48270034214SMatthew G. Knepley 
483a9c22940SMatthew G. Knepley   Note: If the mesh was not distributed, the return value is NULL.
484a9c22940SMatthew G. Knepley 
485a9c22940SMatthew G. Knepley   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
486a9c22940SMatthew G. Knepley   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
487a9c22940SMatthew G. Knepley   representation on the mesh.
48870034214SMatthew G. Knepley 
48970034214SMatthew G. Knepley   Level: intermediate
49070034214SMatthew G. Knepley 
49170034214SMatthew G. Knepley .keywords: mesh, elements
492a9c22940SMatthew G. Knepley .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
49370034214SMatthew G. Knepley @*/
49470034214SMatthew G. Knepley PetscErrorCode DMPlexDistribute(DM dm, const char partitioner[], PetscInt overlap, PetscSF *sf, DM *dmParallel)
49570034214SMatthew G. Knepley {
49670034214SMatthew G. Knepley   DM_Plex               *mesh   = (DM_Plex*) dm->data, *pmesh;
49770034214SMatthew G. Knepley   MPI_Comm               comm;
49870034214SMatthew G. Knepley   const PetscInt         height = 0;
49970034214SMatthew G. Knepley   PetscInt               dim, numRemoteRanks;
50070034214SMatthew G. Knepley   IS                     origCellPart,        origPart,        cellPart,        part;
50170034214SMatthew G. Knepley   PetscSection           origCellPartSection, origPartSection, cellPartSection, partSection;
50270034214SMatthew G. Knepley   PetscSFNode           *remoteRanks;
50370034214SMatthew G. Knepley   PetscSF                partSF, pointSF, coneSF;
50470034214SMatthew G. Knepley   ISLocalToGlobalMapping renumbering;
50570034214SMatthew G. Knepley   PetscSection           originalConeSection, newConeSection;
50670034214SMatthew G. Knepley   PetscInt              *remoteOffsets;
50770034214SMatthew G. Knepley   PetscInt              *cones, *newCones, newConesSize;
50870034214SMatthew G. Knepley   PetscBool              flg;
50970034214SMatthew G. Knepley   PetscMPIInt            rank, numProcs, p;
51070034214SMatthew G. Knepley   PetscErrorCode         ierr;
51170034214SMatthew G. Knepley 
51270034214SMatthew G. Knepley   PetscFunctionBegin;
51370034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
51470034214SMatthew G. Knepley   if (sf) PetscValidPointer(sf,4);
51570034214SMatthew G. Knepley   PetscValidPointer(dmParallel,5);
51670034214SMatthew G. Knepley 
51770034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
51870034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
51970034214SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
52070034214SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
52170034214SMatthew G. Knepley 
52270034214SMatthew G. Knepley   *dmParallel = NULL;
52370034214SMatthew G. Knepley   if (numProcs == 1) PetscFunctionReturn(0);
52470034214SMatthew G. Knepley 
525c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
52670034214SMatthew G. Knepley   /* Create cell partition - We need to rewrite to use IS, use the MatPartition stuff */
52770034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Partition,dm,0,0,0);CHKERRQ(ierr);
52870034214SMatthew G. Knepley   if (overlap > 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Overlap > 1 not yet implemented");
52970034214SMatthew G. Knepley   ierr = DMPlexCreatePartition(dm, partitioner, height, overlap > 0 ? PETSC_TRUE : PETSC_FALSE, &cellPartSection, &cellPart, &origCellPartSection, &origCellPart);CHKERRQ(ierr);
53070034214SMatthew G. Knepley   /* Create SF assuming a serial partition for all processes: Could check for IS length here */
53170034214SMatthew G. Knepley   if (!rank) numRemoteRanks = numProcs;
53270034214SMatthew G. Knepley   else       numRemoteRanks = 0;
53370034214SMatthew G. Knepley   ierr = PetscMalloc1(numRemoteRanks, &remoteRanks);CHKERRQ(ierr);
53470034214SMatthew G. Knepley   for (p = 0; p < numRemoteRanks; ++p) {
53570034214SMatthew G. Knepley     remoteRanks[p].rank  = p;
53670034214SMatthew G. Knepley     remoteRanks[p].index = 0;
53770034214SMatthew G. Knepley   }
53870034214SMatthew G. Knepley   ierr = PetscSFCreate(comm, &partSF);CHKERRQ(ierr);
53970034214SMatthew G. Knepley   ierr = PetscSFSetGraph(partSF, 1, numRemoteRanks, NULL, PETSC_OWN_POINTER, remoteRanks, PETSC_OWN_POINTER);CHKERRQ(ierr);
54070034214SMatthew G. Knepley   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
54170034214SMatthew G. Knepley   if (flg) {
54270034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Cell Partition:\n");CHKERRQ(ierr);
54370034214SMatthew G. Knepley     ierr = PetscSectionView(cellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
54470034214SMatthew G. Knepley     ierr = ISView(cellPart, NULL);CHKERRQ(ierr);
54570034214SMatthew G. Knepley     if (origCellPart) {
54670034214SMatthew G. Knepley       ierr = PetscPrintf(comm, "Original Cell Partition:\n");CHKERRQ(ierr);
54770034214SMatthew G. Knepley       ierr = PetscSectionView(origCellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
54870034214SMatthew G. Knepley       ierr = ISView(origCellPart, NULL);CHKERRQ(ierr);
54970034214SMatthew G. Knepley     }
55070034214SMatthew G. Knepley     ierr = PetscSFView(partSF, NULL);CHKERRQ(ierr);
55170034214SMatthew G. Knepley   }
55270034214SMatthew G. Knepley   /* Close the partition over the mesh */
55370034214SMatthew G. Knepley   ierr = DMPlexCreatePartitionClosure(dm, cellPartSection, cellPart, &partSection, &part);CHKERRQ(ierr);
55470034214SMatthew G. Knepley   ierr = ISDestroy(&cellPart);CHKERRQ(ierr);
55570034214SMatthew G. Knepley   ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr);
55670034214SMatthew G. Knepley   /* Create new mesh */
55770034214SMatthew G. Knepley   ierr  = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr);
558c73cfb54SMatthew G. Knepley   ierr  = DMSetDimension(*dmParallel, dim);CHKERRQ(ierr);
55970034214SMatthew G. Knepley   ierr  = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr);
56070034214SMatthew G. Knepley   pmesh = (DM_Plex*) (*dmParallel)->data;
56170034214SMatthew G. Knepley   /* Distribute sieve points and the global point numbering (replaces creating remote bases) */
56270034214SMatthew G. Knepley   ierr = PetscSFConvertPartition(partSF, partSection, part, &renumbering, &pointSF);CHKERRQ(ierr);
56370034214SMatthew G. Knepley   if (flg) {
56470034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Point Partition:\n");CHKERRQ(ierr);
56570034214SMatthew G. Knepley     ierr = PetscSectionView(partSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
56670034214SMatthew G. Knepley     ierr = ISView(part, NULL);CHKERRQ(ierr);
56770034214SMatthew G. Knepley     ierr = PetscSFView(pointSF, NULL);CHKERRQ(ierr);
56870034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Point Renumbering after partition:\n");CHKERRQ(ierr);
56970034214SMatthew G. Knepley     ierr = ISLocalToGlobalMappingView(renumbering, NULL);CHKERRQ(ierr);
57070034214SMatthew G. Knepley   }
57170034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Partition,dm,0,0,0);CHKERRQ(ierr);
57270034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
57370034214SMatthew G. Knepley   /* Distribute cone section */
57470034214SMatthew G. Knepley   ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr);
57570034214SMatthew G. Knepley   ierr = DMPlexGetConeSection(*dmParallel, &newConeSection);CHKERRQ(ierr);
57670034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr);
57770034214SMatthew G. Knepley   ierr = DMSetUp(*dmParallel);CHKERRQ(ierr);
57870034214SMatthew G. Knepley   {
57970034214SMatthew G. Knepley     PetscInt pStart, pEnd, p;
58070034214SMatthew G. Knepley 
58170034214SMatthew G. Knepley     ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr);
58270034214SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
58370034214SMatthew G. Knepley       PetscInt coneSize;
58470034214SMatthew G. Knepley       ierr               = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr);
58570034214SMatthew G. Knepley       pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize);
58670034214SMatthew G. Knepley     }
58770034214SMatthew G. Knepley   }
58870034214SMatthew G. Knepley   /* Communicate and renumber cones */
58970034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr);
59070034214SMatthew G. Knepley   ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
59170034214SMatthew G. Knepley   ierr = DMPlexGetCones(*dmParallel, &newCones);CHKERRQ(ierr);
59270034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
59370034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
59470034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr);
5959d90f715SBarry Smith   ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr);
59670034214SMatthew G. Knepley   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr);
59770034214SMatthew G. Knepley   if (flg) {
59870034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr);
59970034214SMatthew G. Knepley     ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
60070034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr);
60170034214SMatthew G. Knepley     ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
60270034214SMatthew G. Knepley     ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr);
60370034214SMatthew G. Knepley   }
60470034214SMatthew G. Knepley   ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr);
60570034214SMatthew G. Knepley   ierr = DMPlexGetConeOrientations(*dmParallel, &newCones);CHKERRQ(ierr);
60670034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
60770034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
60870034214SMatthew G. Knepley   ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr);
60970034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
61070034214SMatthew G. Knepley   /* Create supports and stratify sieve */
61170034214SMatthew G. Knepley   {
61270034214SMatthew G. Knepley     PetscInt pStart, pEnd;
61370034214SMatthew G. Knepley 
61470034214SMatthew G. Knepley     ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr);
61570034214SMatthew G. Knepley     ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr);
61670034214SMatthew G. Knepley   }
61770034214SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dmParallel);CHKERRQ(ierr);
61870034214SMatthew G. Knepley   ierr = DMPlexStratify(*dmParallel);CHKERRQ(ierr);
61970034214SMatthew G. Knepley   /* Create point SF for parallel mesh */
62070034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
62170034214SMatthew G. Knepley   {
62270034214SMatthew G. Knepley     const PetscInt *leaves;
62370034214SMatthew G. Knepley     PetscSFNode    *remotePoints, *rowners, *lowners;
62470034214SMatthew G. Knepley     PetscInt        numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints;
62570034214SMatthew G. Knepley     PetscInt        pStart, pEnd;
62670034214SMatthew G. Knepley 
62770034214SMatthew G. Knepley     ierr = DMPlexGetChart(*dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
62870034214SMatthew G. Knepley     ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr);
62970034214SMatthew G. Knepley     ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr);
63070034214SMatthew G. Knepley     for (p=0; p<numRoots; p++) {
63170034214SMatthew G. Knepley       rowners[p].rank  = -1;
63270034214SMatthew G. Knepley       rowners[p].index = -1;
63370034214SMatthew G. Knepley     }
63470034214SMatthew G. Knepley     if (origCellPart) {
63570034214SMatthew G. Knepley       /* Make sure points in the original partition are not assigned to other procs */
63670034214SMatthew G. Knepley       const PetscInt *origPoints;
63770034214SMatthew G. Knepley 
63870034214SMatthew G. Knepley       ierr = DMPlexCreatePartitionClosure(dm, origCellPartSection, origCellPart, &origPartSection, &origPart);CHKERRQ(ierr);
63970034214SMatthew G. Knepley       ierr = ISGetIndices(origPart, &origPoints);CHKERRQ(ierr);
64070034214SMatthew G. Knepley       for (p = 0; p < numProcs; ++p) {
64170034214SMatthew G. Knepley         PetscInt dof, off, d;
64270034214SMatthew G. Knepley 
64370034214SMatthew G. Knepley         ierr = PetscSectionGetDof(origPartSection, p, &dof);CHKERRQ(ierr);
64470034214SMatthew G. Knepley         ierr = PetscSectionGetOffset(origPartSection, p, &off);CHKERRQ(ierr);
64570034214SMatthew G. Knepley         for (d = off; d < off+dof; ++d) {
64670034214SMatthew G. Knepley           rowners[origPoints[d]].rank = p;
64770034214SMatthew G. Knepley         }
64870034214SMatthew G. Knepley       }
64970034214SMatthew G. Knepley       ierr = ISRestoreIndices(origPart, &origPoints);CHKERRQ(ierr);
65070034214SMatthew G. Knepley       ierr = ISDestroy(&origPart);CHKERRQ(ierr);
65170034214SMatthew G. Knepley       ierr = PetscSectionDestroy(&origPartSection);CHKERRQ(ierr);
65270034214SMatthew G. Knepley     }
65370034214SMatthew G. Knepley     ierr = ISDestroy(&origCellPart);CHKERRQ(ierr);
65470034214SMatthew G. Knepley     ierr = PetscSectionDestroy(&origCellPartSection);CHKERRQ(ierr);
65570034214SMatthew G. Knepley 
65670034214SMatthew G. Knepley     ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
65770034214SMatthew G. Knepley     ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
65870034214SMatthew G. Knepley     for (p = 0; p < numLeaves; ++p) {
65970034214SMatthew G. Knepley       if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */
66070034214SMatthew G. Knepley         lowners[p].rank  = rank;
66170034214SMatthew G. Knepley         lowners[p].index = leaves ? leaves[p] : p;
66270034214SMatthew G. Knepley       } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */
66370034214SMatthew G. Knepley         lowners[p].rank  = -2;
66470034214SMatthew G. Knepley         lowners[p].index = -2;
66570034214SMatthew G. Knepley       }
66670034214SMatthew G. Knepley     }
66770034214SMatthew G. Knepley     for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */
66870034214SMatthew G. Knepley       rowners[p].rank  = -3;
66970034214SMatthew G. Knepley       rowners[p].index = -3;
67070034214SMatthew G. Knepley     }
67170034214SMatthew G. Knepley     ierr = PetscSFReduceBegin(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
67270034214SMatthew G. Knepley     ierr = PetscSFReduceEnd(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
67370034214SMatthew G. Knepley     ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
67470034214SMatthew G. Knepley     ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
67570034214SMatthew G. Knepley     for (p = 0; p < numLeaves; ++p) {
67670034214SMatthew G. Knepley       if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed");
67770034214SMatthew G. Knepley       if (lowners[p].rank != rank) ++numGhostPoints;
67870034214SMatthew G. Knepley     }
67970034214SMatthew G. Knepley     ierr = PetscMalloc1(numGhostPoints,    &ghostPoints);CHKERRQ(ierr);
68070034214SMatthew G. Knepley     ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr);
68170034214SMatthew G. Knepley     for (p = 0, gp = 0; p < numLeaves; ++p) {
68270034214SMatthew G. Knepley       if (lowners[p].rank != rank) {
68370034214SMatthew G. Knepley         ghostPoints[gp]        = leaves ? leaves[p] : p;
68470034214SMatthew G. Knepley         remotePoints[gp].rank  = lowners[p].rank;
68570034214SMatthew G. Knepley         remotePoints[gp].index = lowners[p].index;
68670034214SMatthew G. Knepley         ++gp;
68770034214SMatthew G. Knepley       }
68870034214SMatthew G. Knepley     }
68970034214SMatthew G. Knepley     ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr);
69070034214SMatthew G. Knepley     ierr = PetscSFSetGraph((*dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
69170034214SMatthew G. Knepley     ierr = PetscSFSetFromOptions((*dmParallel)->sf);CHKERRQ(ierr);
69270034214SMatthew G. Knepley   }
693a42b08eeSMatthew G. Knepley   pmesh->useCone        = mesh->useCone;
694a42b08eeSMatthew G. Knepley   pmesh->useClosure     = mesh->useClosure;
695c389ea9fSToby Isaac   pmesh->useConstraints = mesh->useConstraints;
69670034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
697bf5244c7SMatthew G. Knepley   /* Distribute Coordinates */
698bf5244c7SMatthew G. Knepley   {
699bf5244c7SMatthew G. Knepley     PetscSection     originalCoordSection, newCoordSection;
700bf5244c7SMatthew G. Knepley     Vec              originalCoordinates, newCoordinates;
701bf5244c7SMatthew G. Knepley     PetscInt         bs;
702bf5244c7SMatthew G. Knepley     const char      *name;
703bf5244c7SMatthew G. Knepley     const PetscReal *maxCell, *L;
704bf5244c7SMatthew G. Knepley 
705bf5244c7SMatthew G. Knepley     ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr);
706bf5244c7SMatthew G. Knepley     ierr = DMGetCoordinateSection(*dmParallel, &newCoordSection);CHKERRQ(ierr);
707bf5244c7SMatthew G. Knepley     ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr);
708bf5244c7SMatthew G. Knepley     ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr);
709bf5244c7SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr);
710bf5244c7SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr);
711bf5244c7SMatthew G. Knepley 
712bf5244c7SMatthew G. Knepley     ierr = DMPlexDistributeField(dm, pointSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr);
713bf5244c7SMatthew G. Knepley     ierr = DMSetCoordinatesLocal(*dmParallel, newCoordinates);CHKERRQ(ierr);
714bf5244c7SMatthew G. Knepley     ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr);
715bf5244c7SMatthew G. Knepley     ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr);
716bf5244c7SMatthew G. Knepley     ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr);
717bf5244c7SMatthew G. Knepley     ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr);
718bf5244c7SMatthew G. Knepley     if (L) {ierr = DMSetPeriodicity(*dmParallel, maxCell, L);CHKERRQ(ierr);}
719bf5244c7SMatthew G. Knepley   }
720bf5244c7SMatthew G. Knepley   /* Distribute labels */
721bf5244c7SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
722bf5244c7SMatthew G. Knepley   {
723bf5244c7SMatthew G. Knepley     DMLabel  next      = mesh->labels, newNext = pmesh->labels;
724bf5244c7SMatthew G. Knepley     PetscInt numLabels = 0, l;
725bf5244c7SMatthew G. Knepley 
726bf5244c7SMatthew G. Knepley     /* Bcast number of labels */
727bf5244c7SMatthew G. Knepley     while (next) {++numLabels; next = next->next;}
728bf5244c7SMatthew G. Knepley     ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
729bf5244c7SMatthew G. Knepley     next = mesh->labels;
730bf5244c7SMatthew G. Knepley     for (l = 0; l < numLabels; ++l) {
731bf5244c7SMatthew G. Knepley       DMLabel   labelNew;
732bf5244c7SMatthew G. Knepley       PetscBool isdepth;
733bf5244c7SMatthew G. Knepley 
734bf5244c7SMatthew G. Knepley       /* Skip "depth" because it is recreated */
735bf5244c7SMatthew G. Knepley       if (!rank) {ierr = PetscStrcmp(next->name, "depth", &isdepth);CHKERRQ(ierr);}
736bf5244c7SMatthew G. Knepley       ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
737bf5244c7SMatthew G. Knepley       if (isdepth) {if (!rank) next = next->next; continue;}
738bf5244c7SMatthew G. Knepley       ierr = DMLabelDistribute(next, partSection, part, renumbering, &labelNew);CHKERRQ(ierr);
739bf5244c7SMatthew G. Knepley       /* Insert into list */
740bf5244c7SMatthew G. Knepley       if (newNext) newNext->next = labelNew;
741bf5244c7SMatthew G. Knepley       else         pmesh->labels = labelNew;
742bf5244c7SMatthew G. Knepley       newNext = labelNew;
743bf5244c7SMatthew G. Knepley       if (!rank) next = next->next;
744bf5244c7SMatthew G. Knepley     }
745bf5244c7SMatthew G. Knepley   }
746bf5244c7SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
747bf5244c7SMatthew G. Knepley   /* Setup hybrid structure */
748bf5244c7SMatthew G. Knepley   {
749bf5244c7SMatthew G. Knepley     const PetscInt *gpoints;
750bf5244c7SMatthew G. Knepley     PetscInt        depth, n, d;
751bf5244c7SMatthew G. Knepley 
752bf5244c7SMatthew G. Knepley     for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];}
753bf5244c7SMatthew G. Knepley     ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr);
754bf5244c7SMatthew G. Knepley     ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr);
755bf5244c7SMatthew G. Knepley     ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr);
756bf5244c7SMatthew G. Knepley     ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
757bf5244c7SMatthew G. Knepley     for (d = 0; d <= dim; ++d) {
758bf5244c7SMatthew G. Knepley       PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p;
759bf5244c7SMatthew G. Knepley 
760bf5244c7SMatthew G. Knepley       if (pmax < 0) continue;
761bf5244c7SMatthew G. Knepley       ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr);
762bf5244c7SMatthew G. Knepley       ierr = DMPlexGetDepthStratum(*dmParallel, d, NULL, &pEnd);CHKERRQ(ierr);
763bf5244c7SMatthew G. Knepley       ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr);
764bf5244c7SMatthew G. Knepley       for (p = 0; p < n; ++p) {
765bf5244c7SMatthew G. Knepley         const PetscInt point = gpoints[p];
766bf5244c7SMatthew G. Knepley 
767bf5244c7SMatthew G. Knepley         if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax;
768bf5244c7SMatthew G. Knepley       }
769bf5244c7SMatthew G. Knepley       if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax;
770bf5244c7SMatthew G. Knepley       else            pmesh->hybridPointMax[d] = -1;
771bf5244c7SMatthew G. Knepley     }
772bf5244c7SMatthew G. Knepley     ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr);
773bf5244c7SMatthew G. Knepley   }
774c389ea9fSToby Isaac   /* Set up tree */
775c389ea9fSToby Isaac   {
776c389ea9fSToby Isaac     DM              refTree;
777c389ea9fSToby Isaac     PetscSection    origParentSection, newParentSection;
778c389ea9fSToby Isaac     PetscInt        *origParents, *origChildIDs;
779c389ea9fSToby Isaac 
780c389ea9fSToby Isaac     ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr);
781c389ea9fSToby Isaac     ierr = DMPlexSetReferenceTree(*dmParallel,refTree);CHKERRQ(ierr);
782c389ea9fSToby Isaac     ierr = DMPlexGetTree(dm,&origParentSection,&origParents,&origChildIDs,NULL,NULL);CHKERRQ(ierr);
783c389ea9fSToby Isaac     if (origParentSection) {
784c389ea9fSToby Isaac       PetscInt        pStart, pEnd;
785c389ea9fSToby Isaac       PetscInt        *newParents, *newChildIDs;
786c389ea9fSToby Isaac       PetscInt        *remoteOffsetsParents, newParentSize;
787c389ea9fSToby Isaac       PetscSF         parentSF;
788c389ea9fSToby Isaac 
789c389ea9fSToby Isaac       ierr = DMPlexGetChart(*dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
790c389ea9fSToby Isaac       ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dmParallel),&newParentSection);CHKERRQ(ierr);
791c389ea9fSToby Isaac       ierr = PetscSectionSetChart(newParentSection,pStart,pEnd);CHKERRQ(ierr);
792c389ea9fSToby Isaac       ierr = PetscSFDistributeSection(pointSF, origParentSection, &remoteOffsetsParents, newParentSection);CHKERRQ(ierr);
793c389ea9fSToby Isaac       ierr = PetscSFCreateSectionSF(pointSF, origParentSection, remoteOffsetsParents, newParentSection, &parentSF);CHKERRQ(ierr);
794c389ea9fSToby Isaac       ierr = PetscSectionGetStorageSize(newParentSection,&newParentSize);CHKERRQ(ierr);
795c389ea9fSToby Isaac       ierr = PetscMalloc2(newParentSize,&newParents,newParentSize,&newChildIDs);CHKERRQ(ierr);
796c389ea9fSToby Isaac       ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
797c389ea9fSToby Isaac       ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
798c389ea9fSToby Isaac       ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
799c389ea9fSToby Isaac       ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
800c389ea9fSToby Isaac       ierr = ISGlobalToLocalMappingApplyBlock(renumbering,IS_GTOLM_MASK, newParentSize, newParents, NULL, newParents);CHKERRQ(ierr);
801c389ea9fSToby Isaac       ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-parents_view", &flg);CHKERRQ(ierr);
802c389ea9fSToby Isaac       if (flg) {
803c389ea9fSToby Isaac         ierr = PetscPrintf(comm, "Serial Parent Section: \n");CHKERRQ(ierr);
804c389ea9fSToby Isaac         ierr = PetscSectionView(origParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
805c389ea9fSToby Isaac         ierr = PetscPrintf(comm, "Parallel Parent Section: \n");CHKERRQ(ierr);
806c389ea9fSToby Isaac         ierr = PetscSectionView(newParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
807c389ea9fSToby Isaac         ierr = PetscSFView(parentSF, NULL);CHKERRQ(ierr);
808c389ea9fSToby Isaac       }
809c389ea9fSToby Isaac       ierr = DMPlexSetTree(*dmParallel,newParentSection,newParents,newChildIDs);CHKERRQ(ierr);
810c389ea9fSToby Isaac       ierr = PetscSectionDestroy(&newParentSection);CHKERRQ(ierr);
811c389ea9fSToby Isaac       ierr = PetscFree2(newParents,newChildIDs);CHKERRQ(ierr);
812c389ea9fSToby Isaac       ierr = PetscSFDestroy(&parentSF);CHKERRQ(ierr);
813c389ea9fSToby Isaac     }
814c389ea9fSToby Isaac   }
815bf5244c7SMatthew G. Knepley   /* Cleanup Partition */
816bf5244c7SMatthew G. Knepley   ierr = ISLocalToGlobalMappingDestroy(&renumbering);CHKERRQ(ierr);
817bf5244c7SMatthew G. Knepley   ierr = PetscSFDestroy(&partSF);CHKERRQ(ierr);
818bf5244c7SMatthew G. Knepley   ierr = PetscSectionDestroy(&partSection);CHKERRQ(ierr);
819bf5244c7SMatthew G. Knepley   ierr = ISDestroy(&part);CHKERRQ(ierr);
82086719b60SMatthew G. Knepley   /* Copy BC */
82186719b60SMatthew G. Knepley   ierr = DMPlexCopyBoundary(dm, *dmParallel);CHKERRQ(ierr);
82270034214SMatthew G. Knepley   /* Cleanup */
82370034214SMatthew G. Knepley   if (sf) {*sf = pointSF;}
82470034214SMatthew G. Knepley   else    {ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr);}
82570034214SMatthew G. Knepley   ierr = DMSetFromOptions(*dmParallel);CHKERRQ(ierr);
82670034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
82770034214SMatthew G. Knepley   PetscFunctionReturn(0);
82870034214SMatthew G. Knepley }
829