xref: /petsc/src/dm/impls/plex/plexdistribute.c (revision bdd2d7519b0c721ec0e026c4dbaee90bb0134da7)
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__
123a17985deSToby Isaac #define __FUNCT__ "DMPlexSetAdjacencyUseAnchors"
1248b0b4c70SToby Isaac /*@
125a17985deSToby Isaac   DMPlexSetAdjacencyUseAnchors - Define adjacency in the mesh using the point-to-point constraints.
1268b0b4c70SToby Isaac 
1278b0b4c70SToby Isaac   Input Parameters:
1288b0b4c70SToby Isaac + dm      - The DM object
1295b317d89SToby Isaac - useAnchors - Flag to use the constraints.  If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place.
1308b0b4c70SToby Isaac 
1318b0b4c70SToby Isaac   Level: intermediate
1328b0b4c70SToby Isaac 
133a17985deSToby Isaac .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
1348b0b4c70SToby Isaac @*/
1355b317d89SToby Isaac PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM dm, PetscBool useAnchors)
1368b0b4c70SToby Isaac {
1378b0b4c70SToby Isaac   DM_Plex *mesh = (DM_Plex *) dm->data;
1388b0b4c70SToby Isaac 
1398b0b4c70SToby Isaac   PetscFunctionBegin;
1408b0b4c70SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1415b317d89SToby Isaac   mesh->useAnchors = useAnchors;
1428b0b4c70SToby Isaac   PetscFunctionReturn(0);
1438b0b4c70SToby Isaac }
1448b0b4c70SToby Isaac 
1458b0b4c70SToby Isaac #undef __FUNCT__
146a17985deSToby Isaac #define __FUNCT__ "DMPlexGetAdjacencyUseAnchors"
1478b0b4c70SToby Isaac /*@
148a17985deSToby Isaac   DMPlexGetAdjacencyUseAnchors - Query whether adjacency in the mesh uses the point-to-point constraints.
1498b0b4c70SToby Isaac 
1508b0b4c70SToby Isaac   Input Parameter:
1518b0b4c70SToby Isaac . dm      - The DM object
1528b0b4c70SToby Isaac 
1538b0b4c70SToby Isaac   Output Parameter:
1545b317d89SToby Isaac . useAnchors - Flag to use the closure.  If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place.
1558b0b4c70SToby Isaac 
1568b0b4c70SToby Isaac   Level: intermediate
1578b0b4c70SToby Isaac 
158a17985deSToby Isaac .seealso: DMPlexSetAdjacencyUseAnchors(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
1598b0b4c70SToby Isaac @*/
1605b317d89SToby Isaac PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM dm, PetscBool *useAnchors)
1618b0b4c70SToby Isaac {
1628b0b4c70SToby Isaac   DM_Plex *mesh = (DM_Plex *) dm->data;
1638b0b4c70SToby Isaac 
1648b0b4c70SToby Isaac   PetscFunctionBegin;
1658b0b4c70SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1665b317d89SToby Isaac   PetscValidIntPointer(useAnchors, 2);
1675b317d89SToby Isaac   *useAnchors = mesh->useAnchors;
1688b0b4c70SToby Isaac   PetscFunctionReturn(0);
1698b0b4c70SToby Isaac }
1708b0b4c70SToby Isaac 
1718b0b4c70SToby Isaac #undef __FUNCT__
17270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Cone_Internal"
17370034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Cone_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
17470034214SMatthew G. Knepley {
17570034214SMatthew G. Knepley   const PetscInt *cone = NULL;
17670034214SMatthew G. Knepley   PetscInt        numAdj = 0, maxAdjSize = *adjSize, coneSize, c;
17770034214SMatthew G. Knepley   PetscErrorCode  ierr;
17870034214SMatthew G. Knepley 
17970034214SMatthew G. Knepley   PetscFunctionBeginHot;
18070034214SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
18170034214SMatthew G. Knepley   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
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"
2575b317d89SToby Isaac PetscErrorCode DMPlexGetAdjacency_Internal(DM dm, PetscInt p, PetscBool useCone, PetscBool useTransitiveClosure, PetscBool useAnchors, 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;
2695b317d89SToby Isaac   if (useAnchors) {
270a17985deSToby 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   }
3005b317d89SToby Isaac   if (useAnchors && 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);
3705b317d89SToby Isaac   ierr = DMPlexGetAdjacency_Internal(dm, p, mesh->useCone, mesh->useClosure, mesh->useAnchors, adjSize, adj);CHKERRQ(ierr);
37170034214SMatthew G. Knepley   PetscFunctionReturn(0);
37270034214SMatthew G. Knepley }
373b0a623aaSMatthew G. Knepley #undef __FUNCT__
374b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexCreateTwoSidedProcessSF"
375b0a623aaSMatthew G. Knepley /*@
376b0a623aaSMatthew G. Knepley   DMPlexCreateTwoSidedProcessSF - Create an SF which just has process connectivity
377b0a623aaSMatthew G. Knepley 
378b0a623aaSMatthew G. Knepley   Collective on DM
379b0a623aaSMatthew G. Knepley 
380b0a623aaSMatthew G. Knepley   Input Parameters:
381b0a623aaSMatthew G. Knepley + dm      - The DM
382b0a623aaSMatthew G. Knepley - sfPoint - The PetscSF which encodes point connectivity
383b0a623aaSMatthew G. Knepley 
384b0a623aaSMatthew G. Knepley   Output Parameters:
385b0a623aaSMatthew G. Knepley + processRanks - A list of process neighbors, or NULL
386b0a623aaSMatthew G. Knepley - sfProcess    - An SF encoding the two-sided process connectivity, or NULL
387b0a623aaSMatthew G. Knepley 
388b0a623aaSMatthew G. Knepley   Level: developer
389b0a623aaSMatthew G. Knepley 
390b0a623aaSMatthew G. Knepley .seealso: PetscSFCreate(), DMPlexCreateProcessSF()
391b0a623aaSMatthew G. Knepley @*/
392b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM dm, PetscSF sfPoint, PetscSection rootRankSection, IS rootRanks, PetscSection leafRankSection, IS leafRanks, IS *processRanks, PetscSF *sfProcess)
393b0a623aaSMatthew G. Knepley {
394b0a623aaSMatthew G. Knepley   const PetscSFNode *remotePoints;
395b0a623aaSMatthew G. Knepley   PetscInt          *localPointsNew;
396b0a623aaSMatthew G. Knepley   PetscSFNode       *remotePointsNew;
397b0a623aaSMatthew G. Knepley   const PetscInt    *nranks;
398b0a623aaSMatthew G. Knepley   PetscInt          *ranksNew;
399b0a623aaSMatthew G. Knepley   PetscBT            neighbors;
400b0a623aaSMatthew G. Knepley   PetscInt           pStart, pEnd, p, numLeaves, l, numNeighbors, n;
401b0a623aaSMatthew G. Knepley   PetscMPIInt        numProcs, proc, rank;
402b0a623aaSMatthew G. Knepley   PetscErrorCode     ierr;
403b0a623aaSMatthew G. Knepley 
404b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
405b0a623aaSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
406b0a623aaSMatthew G. Knepley   PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2);
407b0a623aaSMatthew G. Knepley   if (processRanks) {PetscValidPointer(processRanks, 3);}
408b0a623aaSMatthew G. Knepley   if (sfProcess)    {PetscValidPointer(sfProcess, 4);}
409b0a623aaSMatthew G. Knepley   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr);
410b0a623aaSMatthew G. Knepley   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
411b0a623aaSMatthew G. Knepley   ierr = PetscSFGetGraph(sfPoint, NULL, &numLeaves, NULL, &remotePoints);CHKERRQ(ierr);
412b0a623aaSMatthew G. Knepley   ierr = PetscBTCreate(numProcs, &neighbors);CHKERRQ(ierr);
413b0a623aaSMatthew G. Knepley   ierr = PetscBTMemzero(numProcs, neighbors);CHKERRQ(ierr);
414b0a623aaSMatthew G. Knepley   /* Compute root-to-leaf process connectivity */
415b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetChart(rootRankSection, &pStart, &pEnd);CHKERRQ(ierr);
416b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(rootRanks, &nranks);CHKERRQ(ierr);
417b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
418b0a623aaSMatthew G. Knepley     PetscInt ndof, noff, n;
419b0a623aaSMatthew G. Knepley 
420b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetDof(rootRankSection, p, &ndof);CHKERRQ(ierr);
421b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(rootRankSection, p, &noff);CHKERRQ(ierr);
422b0a623aaSMatthew G. Knepley     for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);}
423b0a623aaSMatthew G. Knepley   }
424b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(rootRanks, &nranks);CHKERRQ(ierr);
425b0a623aaSMatthew G. Knepley   /* Compute leaf-to-neighbor process connectivity */
426b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetChart(leafRankSection, &pStart, &pEnd);CHKERRQ(ierr);
427b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(leafRanks, &nranks);CHKERRQ(ierr);
428b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
429b0a623aaSMatthew G. Knepley     PetscInt ndof, noff, n;
430b0a623aaSMatthew G. Knepley 
431b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetDof(leafRankSection, p, &ndof);CHKERRQ(ierr);
432b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(leafRankSection, p, &noff);CHKERRQ(ierr);
433b0a623aaSMatthew G. Knepley     for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);}
434b0a623aaSMatthew G. Knepley   }
435b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(leafRanks, &nranks);CHKERRQ(ierr);
436b0a623aaSMatthew G. Knepley   /* Compute leaf-to-root process connectivity */
437b0a623aaSMatthew G. Knepley   for (l = 0; l < numLeaves; ++l) {PetscBTSet(neighbors, remotePoints[l].rank);}
438b0a623aaSMatthew G. Knepley   /* Calculate edges */
439b0a623aaSMatthew G. Knepley   PetscBTClear(neighbors, rank);
440b0a623aaSMatthew G. Knepley   for(proc = 0, numNeighbors = 0; proc < numProcs; ++proc) {if (PetscBTLookup(neighbors, proc)) ++numNeighbors;}
441b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(numNeighbors, &ranksNew);CHKERRQ(ierr);
442b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(numNeighbors, &localPointsNew);CHKERRQ(ierr);
443b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(numNeighbors, &remotePointsNew);CHKERRQ(ierr);
444b0a623aaSMatthew G. Knepley   for(proc = 0, n = 0; proc < numProcs; ++proc) {
445b0a623aaSMatthew G. Knepley     if (PetscBTLookup(neighbors, proc)) {
446b0a623aaSMatthew G. Knepley       ranksNew[n]              = proc;
447b0a623aaSMatthew G. Knepley       localPointsNew[n]        = proc;
448b0a623aaSMatthew G. Knepley       remotePointsNew[n].index = rank;
449b0a623aaSMatthew G. Knepley       remotePointsNew[n].rank  = proc;
450b0a623aaSMatthew G. Knepley       ++n;
451b0a623aaSMatthew G. Knepley     }
452b0a623aaSMatthew G. Knepley   }
453b0a623aaSMatthew G. Knepley   ierr = PetscBTDestroy(&neighbors);CHKERRQ(ierr);
454b0a623aaSMatthew G. Knepley   if (processRanks) {ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), numNeighbors, ranksNew, PETSC_OWN_POINTER, processRanks);CHKERRQ(ierr);}
455b0a623aaSMatthew G. Knepley   else              {ierr = PetscFree(ranksNew);CHKERRQ(ierr);}
456b0a623aaSMatthew G. Knepley   if (sfProcess) {
457b0a623aaSMatthew G. Knepley     ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess);CHKERRQ(ierr);
458b0a623aaSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) *sfProcess, "Two-Sided Process SF");CHKERRQ(ierr);
459b0a623aaSMatthew G. Knepley     ierr = PetscSFSetFromOptions(*sfProcess);CHKERRQ(ierr);
460b0a623aaSMatthew G. Knepley     ierr = PetscSFSetGraph(*sfProcess, numProcs, numNeighbors, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr);
461b0a623aaSMatthew G. Knepley   }
462b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
463b0a623aaSMatthew G. Knepley }
464b0a623aaSMatthew G. Knepley 
465b0a623aaSMatthew G. Knepley #undef __FUNCT__
466b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeOwnership"
467b0a623aaSMatthew G. Knepley /*@
468b0a623aaSMatthew G. Knepley   DMPlexDistributeOwnership - Compute owner information for shared points. This basically gets two-sided for an SF.
469b0a623aaSMatthew G. Knepley 
470b0a623aaSMatthew G. Knepley   Collective on DM
471b0a623aaSMatthew G. Knepley 
472b0a623aaSMatthew G. Knepley   Input Parameter:
473b0a623aaSMatthew G. Knepley . dm - The DM
474b0a623aaSMatthew G. Knepley 
475b0a623aaSMatthew G. Knepley   Output Parameters:
476b0a623aaSMatthew G. Knepley + rootSection - The number of leaves for a given root point
477b0a623aaSMatthew G. Knepley . rootrank    - The rank of each edge into the root point
478b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point
479b0a623aaSMatthew G. Knepley - leafrank    - The rank of each process sharing a leaf point
480b0a623aaSMatthew G. Knepley 
481b0a623aaSMatthew G. Knepley   Level: developer
482b0a623aaSMatthew G. Knepley 
483b0a623aaSMatthew G. Knepley .seealso: DMPlexCreateOverlap()
484b0a623aaSMatthew G. Knepley @*/
485b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexDistributeOwnership(DM dm, PetscSection rootSection, IS *rootrank, PetscSection leafSection, IS *leafrank)
486b0a623aaSMatthew G. Knepley {
487b0a623aaSMatthew G. Knepley   MPI_Comm        comm;
488b0a623aaSMatthew G. Knepley   PetscSF         sfPoint;
489b0a623aaSMatthew G. Knepley   const PetscInt *rootdegree;
490b0a623aaSMatthew G. Knepley   PetscInt       *myrank, *remoterank;
491b0a623aaSMatthew G. Knepley   PetscInt        pStart, pEnd, p, nedges;
492b0a623aaSMatthew G. Knepley   PetscMPIInt     rank;
493b0a623aaSMatthew G. Knepley   PetscErrorCode  ierr;
494b0a623aaSMatthew G. Knepley 
495b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
496b0a623aaSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
497b0a623aaSMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
498b0a623aaSMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
499b0a623aaSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
500b0a623aaSMatthew G. Knepley   /* Compute number of leaves for each root */
501b0a623aaSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) rootSection, "Root Section");CHKERRQ(ierr);
502b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetChart(rootSection, pStart, pEnd);CHKERRQ(ierr);
503b0a623aaSMatthew G. Knepley   ierr = PetscSFComputeDegreeBegin(sfPoint, &rootdegree);CHKERRQ(ierr);
504b0a623aaSMatthew G. Knepley   ierr = PetscSFComputeDegreeEnd(sfPoint, &rootdegree);CHKERRQ(ierr);
505b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {ierr = PetscSectionSetDof(rootSection, p, rootdegree[p-pStart]);CHKERRQ(ierr);}
506b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
507b0a623aaSMatthew G. Knepley   /* Gather rank of each leaf to root */
508b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(rootSection, &nedges);CHKERRQ(ierr);
509b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(pEnd-pStart, &myrank);CHKERRQ(ierr);
510b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(nedges,  &remoterank);CHKERRQ(ierr);
511b0a623aaSMatthew G. Knepley   for (p = 0; p < pEnd-pStart; ++p) myrank[p] = rank;
512b0a623aaSMatthew G. Knepley   ierr = PetscSFGatherBegin(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr);
513b0a623aaSMatthew G. Knepley   ierr = PetscSFGatherEnd(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr);
514b0a623aaSMatthew G. Knepley   ierr = PetscFree(myrank);CHKERRQ(ierr);
515b0a623aaSMatthew G. Knepley   ierr = ISCreateGeneral(comm, nedges, remoterank, PETSC_OWN_POINTER, rootrank);CHKERRQ(ierr);
516b0a623aaSMatthew G. Knepley   /* Distribute remote ranks to leaves */
517b0a623aaSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) leafSection, "Leaf Section");CHKERRQ(ierr);
518b0a623aaSMatthew G. Knepley   ierr = DMPlexDistributeFieldIS(dm, sfPoint, rootSection, *rootrank, leafSection, leafrank);CHKERRQ(ierr);
519b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
520b0a623aaSMatthew G. Knepley }
521b0a623aaSMatthew G. Knepley 
522b0a623aaSMatthew G. Knepley #undef __FUNCT__
523b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexCreateOverlap"
524278397a0SMatthew G. Knepley /*@C
525b0a623aaSMatthew G. Knepley   DMPlexCreateOverlap - Compute owner information for shared points. This basically gets two-sided for an SF.
526b0a623aaSMatthew G. Knepley 
527b0a623aaSMatthew G. Knepley   Collective on DM
528b0a623aaSMatthew G. Knepley 
529b0a623aaSMatthew G. Knepley   Input Parameters:
530b0a623aaSMatthew G. Knepley + dm          - The DM
531b0a623aaSMatthew G. Knepley . rootSection - The number of leaves for a given root point
532b0a623aaSMatthew G. Knepley . rootrank    - The rank of each edge into the root point
533b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point
534b0a623aaSMatthew G. Knepley - leafrank    - The rank of each process sharing a leaf point
535b0a623aaSMatthew G. Knepley 
536b0a623aaSMatthew G. Knepley   Output Parameters:
537b0a623aaSMatthew G. Knepley + ovRootSection - The number of new overlap points for each neighboring process
538b0a623aaSMatthew G. Knepley . ovRootPoints  - The new overlap points for each neighboring process
539b0a623aaSMatthew G. Knepley . ovLeafSection - The number of new overlap points from each neighboring process
540b0a623aaSMatthew G. Knepley - ovLeafPoints  - The new overlap points from each neighboring process
541b0a623aaSMatthew G. Knepley 
542b0a623aaSMatthew G. Knepley   Level: developer
543b0a623aaSMatthew G. Knepley 
544b0a623aaSMatthew G. Knepley .seealso: DMPlexDistributeOwnership()
545b0a623aaSMatthew G. Knepley @*/
546e540f424SMichael Lange PetscErrorCode DMPlexCreateOverlap(DM dm, PetscSection rootSection, IS rootrank, PetscSection leafSection, IS leafrank, PetscSF *overlapSF)
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   IS                 valueIS;
552e540f424SMichael Lange   DMLabel            ovLeafLabel;
553b0a623aaSMatthew G. Knepley   const PetscSFNode *remote;
554b0a623aaSMatthew G. Knepley   const PetscInt    *local;
555b0a623aaSMatthew G. Knepley   const PetscInt    *nrank, *rrank, *neighbors;
556e540f424SMichael Lange   PetscSFNode       *ovRootPoints, *ovLeafPoints, *remotePoints;
557e540f424SMichael Lange   PetscSection       ovRootSection, ovLeafSection;
558b0a623aaSMatthew G. Knepley   PetscInt          *adj = NULL;
559b0a623aaSMatthew G. Knepley   PetscInt           pStart, pEnd, p, sStart, sEnd, nleaves, l, numNeighbors, n, ovSize;
560e540f424SMichael Lange   PetscInt           idx, numRemote;
561b0a623aaSMatthew G. Knepley   PetscMPIInt        rank, numProcs;
562e540f424SMichael Lange   PetscBool          flg;
563b0a623aaSMatthew G. Knepley   PetscErrorCode     ierr;
564b0a623aaSMatthew G. Knepley 
565b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
566e540f424SMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
567e540f424SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
568e540f424SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
569b0a623aaSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
570b0a623aaSMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
571b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetChart(leafSection, &sStart, &sEnd);CHKERRQ(ierr);
572b0a623aaSMatthew G. Knepley   ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr);
573b0a623aaSMatthew G. Knepley   ierr = DMLabelCreate("Overlap adjacency", &ovAdjByRank);CHKERRQ(ierr);
574b0a623aaSMatthew G. Knepley   /* Handle leaves: shared with the root point */
575b0a623aaSMatthew G. Knepley   for (l = 0; l < nleaves; ++l) {
576b0a623aaSMatthew G. Knepley     PetscInt adjSize = PETSC_DETERMINE, a;
577b0a623aaSMatthew G. Knepley 
578b0a623aaSMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, local[l], &adjSize, &adj);CHKERRQ(ierr);
579b0a623aaSMatthew G. Knepley     for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remote[l].rank);CHKERRQ(ierr);}
580b0a623aaSMatthew G. Knepley   }
581b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(rootrank, &rrank);CHKERRQ(ierr);
582b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(leafrank, &nrank);CHKERRQ(ierr);
583b0a623aaSMatthew G. Knepley   /* Handle roots */
584b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
585b0a623aaSMatthew G. Knepley     PetscInt adjSize = PETSC_DETERMINE, neighbors = 0, noff, n, a;
586b0a623aaSMatthew G. Knepley 
587b0a623aaSMatthew G. Knepley     if ((p >= sStart) && (p < sEnd)) {
588b0a623aaSMatthew G. Knepley       /* Some leaves share a root with other leaves on different processes */
589b0a623aaSMatthew G. Knepley       ierr = PetscSectionGetDof(leafSection, p, &neighbors);CHKERRQ(ierr);
590b0a623aaSMatthew G. Knepley       if (neighbors) {
591b0a623aaSMatthew G. Knepley         ierr = PetscSectionGetOffset(leafSection, p, &noff);CHKERRQ(ierr);
592b0a623aaSMatthew G. Knepley         ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
593b0a623aaSMatthew G. Knepley         for (n = 0; n < neighbors; ++n) {
594b0a623aaSMatthew G. Knepley           const PetscInt remoteRank = nrank[noff+n];
595b0a623aaSMatthew G. Knepley 
596b0a623aaSMatthew G. Knepley           if (remoteRank == rank) continue;
597b0a623aaSMatthew G. Knepley           for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);}
598b0a623aaSMatthew G. Knepley         }
599b0a623aaSMatthew G. Knepley       }
600b0a623aaSMatthew G. Knepley     }
601b0a623aaSMatthew G. Knepley     /* Roots are shared with leaves */
602b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetDof(rootSection, p, &neighbors);CHKERRQ(ierr);
603b0a623aaSMatthew G. Knepley     if (!neighbors) continue;
604b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(rootSection, p, &noff);CHKERRQ(ierr);
605b0a623aaSMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
606b0a623aaSMatthew G. Knepley     for (n = 0; n < neighbors; ++n) {
607b0a623aaSMatthew G. Knepley       const PetscInt remoteRank = rrank[noff+n];
608b0a623aaSMatthew G. Knepley 
609b0a623aaSMatthew G. Knepley       if (remoteRank == rank) continue;
610b0a623aaSMatthew G. Knepley       for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);}
611b0a623aaSMatthew G. Knepley     }
612b0a623aaSMatthew G. Knepley   }
613b0a623aaSMatthew G. Knepley   ierr = PetscFree(adj);CHKERRQ(ierr);
614b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(rootrank, &rrank);CHKERRQ(ierr);
615b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(leafrank, &nrank);CHKERRQ(ierr);
616e540f424SMichael Lange   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-overlap_view", &flg);CHKERRQ(ierr);
617e540f424SMichael Lange   if (flg) {
618b0a623aaSMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr);
619b0a623aaSMatthew G. Knepley     ierr = DMLabelView(ovAdjByRank, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
620b0a623aaSMatthew G. Knepley   }
621b0a623aaSMatthew G. Knepley   /* Convert to (point, rank) and use actual owners */
622e540f424SMichael Lange   ierr = PetscSectionCreate(comm, &ovRootSection);CHKERRQ(ierr);
623e540f424SMichael Lange   ierr = PetscObjectSetName((PetscObject) ovRootSection, "Overlap Root Section");CHKERRQ(ierr);
624b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetChart(ovRootSection, 0, numProcs);CHKERRQ(ierr);
625b0a623aaSMatthew G. Knepley   ierr = DMLabelGetValueIS(ovAdjByRank, &valueIS);CHKERRQ(ierr);
626b0a623aaSMatthew G. Knepley   ierr = ISGetLocalSize(valueIS, &numNeighbors);CHKERRQ(ierr);
627b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(valueIS, &neighbors);CHKERRQ(ierr);
628b0a623aaSMatthew G. Knepley   for (n = 0; n < numNeighbors; ++n) {
629b0a623aaSMatthew G. Knepley     PetscInt numPoints;
630b0a623aaSMatthew G. Knepley 
631b0a623aaSMatthew G. Knepley     ierr = DMLabelGetStratumSize(ovAdjByRank, neighbors[n], &numPoints);CHKERRQ(ierr);
632b0a623aaSMatthew G. Knepley     ierr = PetscSectionAddDof(ovRootSection, neighbors[n], numPoints);CHKERRQ(ierr);
633b0a623aaSMatthew G. Knepley   }
634b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetUp(ovRootSection);CHKERRQ(ierr);
635b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(ovRootSection, &ovSize);CHKERRQ(ierr);
636e540f424SMichael Lange   ierr = PetscMalloc1(ovSize, &ovRootPoints);CHKERRQ(ierr);
637b0a623aaSMatthew G. Knepley   ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr);
638b0a623aaSMatthew G. Knepley   for (n = 0; n < numNeighbors; ++n) {
639b0a623aaSMatthew G. Knepley     IS              pointIS;
640b0a623aaSMatthew G. Knepley     const PetscInt *points;
641b0a623aaSMatthew G. Knepley     PetscInt        off, numPoints, p;
642b0a623aaSMatthew G. Knepley 
643b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(ovRootSection, neighbors[n], &off);CHKERRQ(ierr);
644b0a623aaSMatthew G. Knepley     ierr = DMLabelGetStratumIS(ovAdjByRank, neighbors[n], &pointIS);CHKERRQ(ierr);
645b0a623aaSMatthew G. Knepley     ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr);
646b0a623aaSMatthew G. Knepley     ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr);
647b0a623aaSMatthew G. Knepley     for (p = 0; p < numPoints; ++p) {
648b0a623aaSMatthew G. Knepley       ierr = PetscFindInt(points[p], nleaves, local, &l);CHKERRQ(ierr);
649e540f424SMichael Lange       if (l >= 0) {ovRootPoints[off+p] = remote[l];}
650e540f424SMichael Lange       else        {ovRootPoints[off+p].index = points[p]; ovRootPoints[off+p].rank = rank;}
651b0a623aaSMatthew G. Knepley     }
652b0a623aaSMatthew G. Knepley     ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr);
653b0a623aaSMatthew G. Knepley     ierr = ISDestroy(&pointIS);CHKERRQ(ierr);
654b0a623aaSMatthew G. Knepley   }
655b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(valueIS, &neighbors);CHKERRQ(ierr);
656b0a623aaSMatthew G. Knepley   ierr = ISDestroy(&valueIS);CHKERRQ(ierr);
657b0a623aaSMatthew G. Knepley   ierr = DMLabelDestroy(&ovAdjByRank);CHKERRQ(ierr);
658b0a623aaSMatthew G. Knepley   /* Make process SF */
659b0a623aaSMatthew G. Knepley   ierr = DMPlexCreateTwoSidedProcessSF(dm, sfPoint, rootSection, rootrank, leafSection, leafrank, NULL, &sfProc);CHKERRQ(ierr);
660e540f424SMichael Lange   if (flg) {
661b0a623aaSMatthew G. Knepley     ierr = PetscSFView(sfProc, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
662b0a623aaSMatthew G. Knepley   }
663b0a623aaSMatthew G. Knepley   /* Communicate overlap */
664e540f424SMichael Lange   ierr = PetscSectionCreate(comm, &ovLeafSection);CHKERRQ(ierr);
665e540f424SMichael Lange   ierr = PetscObjectSetName((PetscObject) ovLeafSection, "Overlap Leaf Section");CHKERRQ(ierr);
666e540f424SMichael Lange   ierr = DMPlexDistributeData(dm, sfProc, ovRootSection, MPIU_2INT, ovRootPoints, ovLeafSection, (void**) &ovLeafPoints);CHKERRQ(ierr);
667e540f424SMichael Lange   /* Filter remote contributions (ovLeafPoints) into the overlapSF */
668e540f424SMichael Lange   ierr = DMLabelCreate("ovLeafs", &ovLeafLabel);CHKERRQ(ierr);
669e540f424SMichael Lange   ierr = PetscSectionGetStorageSize(ovLeafSection, &ovSize);CHKERRQ(ierr);
670e540f424SMichael Lange   for (p = 0; p < ovSize; p++) {
671e540f424SMichael Lange     /* Don't import points from yourself */
672e540f424SMichael Lange     if (ovLeafPoints[p].rank == rank) continue;
673e540f424SMichael Lange     ierr = DMLabelSetValue(ovLeafLabel, ovLeafPoints[p].index, ovLeafPoints[p].rank);CHKERRQ(ierr);
674e540f424SMichael Lange   }
675e540f424SMichael Lange   /* Don't import points already in the pointSF */
676e540f424SMichael Lange   for (l = 0; l < nleaves; ++l) {
677e540f424SMichael Lange     ierr = DMLabelClearValue(ovLeafLabel, remote[l].index, remote[l].rank);CHKERRQ(ierr);
678e540f424SMichael Lange   }
679e540f424SMichael Lange   for (numRemote = 0, n = 0; n < numProcs; ++n) {
680e540f424SMichael Lange     PetscInt numPoints;
681e540f424SMichael Lange     ierr = DMLabelGetStratumSize(ovLeafLabel, n, &numPoints);CHKERRQ(ierr);
682e540f424SMichael Lange     numRemote += numPoints;
683e540f424SMichael Lange   }
684e540f424SMichael Lange   ierr = PetscMalloc1(numRemote, &remotePoints);CHKERRQ(ierr);
685e540f424SMichael Lange   for (idx = 0, n = 0; n < numProcs; ++n) {
686e540f424SMichael Lange     IS remoteRootIS;
687e540f424SMichael Lange     PetscInt numPoints;
688e540f424SMichael Lange     const PetscInt *remoteRoots;
689e540f424SMichael Lange     ierr = DMLabelGetStratumSize(ovLeafLabel, n, &numPoints);CHKERRQ(ierr);
690e540f424SMichael Lange     if (numPoints <= 0) continue;
691e540f424SMichael Lange     ierr = DMLabelGetStratumIS(ovLeafLabel, n, &remoteRootIS);CHKERRQ(ierr);
692e540f424SMichael Lange     ierr = ISGetIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
693e540f424SMichael Lange     for (p = 0; p < numPoints; p++) {
694e540f424SMichael Lange       remotePoints[idx].index = remoteRoots[p];
695e540f424SMichael Lange       remotePoints[idx].rank = n;
696e540f424SMichael Lange       idx++;
697e540f424SMichael Lange     }
698e540f424SMichael Lange     ierr = ISRestoreIndices(remoteRootIS, &remoteRoots);CHKERRQ(ierr);
699e540f424SMichael Lange   }
700e540f424SMichael Lange   ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), overlapSF);CHKERRQ(ierr);
701e540f424SMichael Lange   ierr = PetscObjectSetName((PetscObject) *overlapSF, "Overlap SF");CHKERRQ(ierr);
702e540f424SMichael Lange   ierr = PetscSFSetFromOptions(*overlapSF);CHKERRQ(ierr);
703e540f424SMichael Lange   ierr = PetscSFSetGraph(*overlapSF, pEnd-pStart, numRemote, NULL, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
704e540f424SMichael Lange   if (flg) {
705e540f424SMichael Lange     ierr = PetscPrintf(comm, "Overlap SF\n");CHKERRQ(ierr);
706e540f424SMichael Lange     ierr = PetscSFView(*overlapSF, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
707e540f424SMichael Lange   }
708e540f424SMichael Lange   /* Clean up */
709b0a623aaSMatthew G. Knepley   ierr = PetscSFDestroy(&sfProc);CHKERRQ(ierr);
710e540f424SMichael Lange   ierr = PetscFree(ovRootPoints);CHKERRQ(ierr);
711e540f424SMichael Lange   ierr = PetscSectionDestroy(&ovRootSection);CHKERRQ(ierr);
712e540f424SMichael Lange   ierr = PetscFree(ovLeafPoints);CHKERRQ(ierr);
713e540f424SMichael Lange   ierr = PetscSectionDestroy(&ovLeafSection);CHKERRQ(ierr);
714b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
715b0a623aaSMatthew G. Knepley }
71670034214SMatthew G. Knepley 
71770034214SMatthew G. Knepley #undef __FUNCT__
71846f9b1c3SMichael Lange #define __FUNCT__ "DMPlexCreateOverlapMigrationSF"
71946f9b1c3SMichael Lange PetscErrorCode DMPlexCreateOverlapMigrationSF(DM dm, PetscSF overlapSF, PetscSF *migrationSF)
72046f9b1c3SMichael Lange {
72146f9b1c3SMichael Lange   MPI_Comm           comm;
72246f9b1c3SMichael Lange   PetscMPIInt        rank, numProcs;
72346f9b1c3SMichael Lange   PetscInt           d, dim, p, pStart, pEnd, nroots, nleaves, newLeaves, point, numSharedPoints;
72446f9b1c3SMichael Lange   PetscInt          *pointDepths, *remoteDepths, *ilocal;
72546f9b1c3SMichael Lange   PetscInt          *depthRecv, *depthShift, *depthIdx;
72646f9b1c3SMichael Lange   PetscSFNode       *iremote;
72746f9b1c3SMichael Lange   PetscSF            pointSF;
72846f9b1c3SMichael Lange   const PetscInt    *sharedLocal;
72946f9b1c3SMichael Lange   const PetscSFNode *overlapRemote, *sharedRemote;
73046f9b1c3SMichael Lange   PetscErrorCode     ierr;
73146f9b1c3SMichael Lange 
73246f9b1c3SMichael Lange   PetscFunctionBegin;
73346f9b1c3SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
73446f9b1c3SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
73546f9b1c3SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
73646f9b1c3SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
73746f9b1c3SMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
73846f9b1c3SMichael Lange 
73946f9b1c3SMichael Lange   /* Before building the migration SF we need to know the new stratum offsets */
74046f9b1c3SMichael Lange   ierr = PetscSFGetGraph(overlapSF, &nroots, &nleaves, NULL, &overlapRemote);CHKERRQ(ierr);
74146f9b1c3SMichael Lange   ierr = PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);CHKERRQ(ierr);
74246f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
74346f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
74446f9b1c3SMichael Lange     for (p=pStart; p<pEnd; p++) pointDepths[p] = d;
74546f9b1c3SMichael Lange   }
74646f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) remoteDepths[p] = -1;
74746f9b1c3SMichael Lange   ierr = PetscSFBcastBegin(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
74846f9b1c3SMichael Lange   ierr = PetscSFBcastEnd(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
74946f9b1c3SMichael Lange 
75046f9b1c3SMichael Lange   /* Count recevied points in each stratum and compute the internal strata shift */
75146f9b1c3SMichael Lange   ierr = PetscMalloc3(dim+1, &depthRecv, dim+1, &depthShift, dim+1, &depthIdx);CHKERRQ(ierr);
75246f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) depthRecv[d]=0;
75346f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) depthRecv[remoteDepths[p]]++;
75446f9b1c3SMichael Lange   depthShift[dim] = 0;
75546f9b1c3SMichael Lange   for (d=0; d<dim; d++) depthShift[d] = depthRecv[dim];
75646f9b1c3SMichael Lange   for (d=1; d<dim; d++) depthShift[d] += depthRecv[0];
75746f9b1c3SMichael Lange   for (d=dim-2; d>0; d--) depthShift[d] += depthRecv[d+1];
75846f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
75946f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
76046f9b1c3SMichael Lange     depthIdx[d] = pStart + depthShift[d];
76146f9b1c3SMichael Lange   }
76246f9b1c3SMichael Lange 
76346f9b1c3SMichael Lange   /* Form the overlap SF build an SF that describes the full overlap migration SF */
76446f9b1c3SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
76546f9b1c3SMichael Lange   newLeaves = pEnd - pStart + nleaves;
76646f9b1c3SMichael Lange   ierr = PetscMalloc2(newLeaves, &ilocal, newLeaves, &iremote);CHKERRQ(ierr);
76746f9b1c3SMichael Lange   /* First map local points to themselves */
76846f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
76946f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
77046f9b1c3SMichael Lange     for (p=pStart; p<pEnd; p++) {
77146f9b1c3SMichael Lange       point = p + depthShift[d];
77246f9b1c3SMichael Lange       ilocal[point] = point;
77346f9b1c3SMichael Lange       iremote[point].index = p;
77446f9b1c3SMichael Lange       iremote[point].rank = rank;
77546f9b1c3SMichael Lange       depthIdx[d]++;
77646f9b1c3SMichael Lange     }
77746f9b1c3SMichael Lange   }
77846f9b1c3SMichael Lange 
77946f9b1c3SMichael Lange   /* Add in the remote roots for currently shared points */
78046f9b1c3SMichael Lange   ierr = DMGetPointSF(dm, &pointSF);CHKERRQ(ierr);
78146f9b1c3SMichael Lange   ierr = PetscSFGetGraph(pointSF, NULL, &numSharedPoints, &sharedLocal, &sharedRemote);CHKERRQ(ierr);
78246f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
78346f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
78446f9b1c3SMichael Lange     for (p=0; p<numSharedPoints; p++) {
78546f9b1c3SMichael Lange       if (pStart <= sharedLocal[p] && sharedLocal[p] < pEnd) {
78646f9b1c3SMichael Lange         point = sharedLocal[p] + depthShift[d];
78746f9b1c3SMichael Lange         iremote[point].index = sharedRemote[p].index;
78846f9b1c3SMichael Lange         iremote[point].rank = sharedRemote[p].rank;
78946f9b1c3SMichael Lange       }
79046f9b1c3SMichael Lange     }
79146f9b1c3SMichael Lange   }
79246f9b1c3SMichael Lange 
79346f9b1c3SMichael Lange   /* Now add the incoming overlap points */
79446f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) {
79546f9b1c3SMichael Lange     point = depthIdx[remoteDepths[p]];
79646f9b1c3SMichael Lange     ilocal[point] = point;
79746f9b1c3SMichael Lange     iremote[point].index = overlapRemote[p].index;
79846f9b1c3SMichael Lange     iremote[point].rank = overlapRemote[p].rank;
79946f9b1c3SMichael Lange     depthIdx[remoteDepths[p]]++;
80046f9b1c3SMichael Lange   }
80146f9b1c3SMichael Lange 
80246f9b1c3SMichael Lange   ierr = PetscSFCreate(comm, migrationSF);CHKERRQ(ierr);
80346f9b1c3SMichael Lange   ierr = PetscObjectSetName((PetscObject) *migrationSF, "Overlap Migration SF");CHKERRQ(ierr);
80446f9b1c3SMichael Lange   ierr = PetscSFSetFromOptions(*migrationSF);CHKERRQ(ierr);
80546f9b1c3SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
80646f9b1c3SMichael Lange   ierr = PetscSFSetGraph(*migrationSF, pEnd-pStart, newLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
80746f9b1c3SMichael Lange 
80846f9b1c3SMichael Lange   ierr = PetscFree3(depthRecv, depthShift, depthIdx);CHKERRQ(ierr);
80946f9b1c3SMichael Lange   PetscFunctionReturn(0);
81046f9b1c3SMichael Lange }
81146f9b1c3SMichael Lange 
81246f9b1c3SMichael Lange #undef __FUNCT__
81370034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeField"
81470034214SMatthew G. Knepley /*@
81570034214SMatthew G. Knepley   DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
81670034214SMatthew G. Knepley 
81770034214SMatthew G. Knepley   Collective on DM
81870034214SMatthew G. Knepley 
81970034214SMatthew G. Knepley   Input Parameters:
82070034214SMatthew G. Knepley + dm - The DMPlex object
82170034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
82270034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
82370034214SMatthew G. Knepley - originalVec - The existing data
82470034214SMatthew G. Knepley 
82570034214SMatthew G. Knepley   Output Parameters:
82670034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout
82770034214SMatthew G. Knepley - newVec - The new data
82870034214SMatthew G. Knepley 
82970034214SMatthew G. Knepley   Level: developer
83070034214SMatthew G. Knepley 
8311e8fc25dSMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeFieldIS(), DMPlexDistributeData()
83270034214SMatthew G. Knepley @*/
83370034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec)
83470034214SMatthew G. Knepley {
83570034214SMatthew G. Knepley   PetscSF        fieldSF;
83670034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
83770034214SMatthew G. Knepley   PetscScalar   *originalValues, *newValues;
83870034214SMatthew G. Knepley   PetscErrorCode ierr;
83970034214SMatthew G. Knepley 
84070034214SMatthew G. Knepley   PetscFunctionBegin;
84170034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
84270034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
84370034214SMatthew G. Knepley 
84470034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
84570034214SMatthew G. Knepley   ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr);
84670034214SMatthew G. Knepley   ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr);
84770034214SMatthew G. Knepley 
84870034214SMatthew G. Knepley   ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr);
84970034214SMatthew G. Knepley   ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr);
85070034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
85170034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
85270034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
85370034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
85470034214SMatthew G. Knepley   ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr);
85570034214SMatthew G. Knepley   ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr);
85670034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
85770034214SMatthew G. Knepley   PetscFunctionReturn(0);
85870034214SMatthew G. Knepley }
85970034214SMatthew G. Knepley 
86070034214SMatthew G. Knepley #undef __FUNCT__
8611e8fc25dSMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeFieldIS"
8621e8fc25dSMatthew G. Knepley /*@
8631e8fc25dSMatthew G. Knepley   DMPlexDistributeFieldIS - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
8641e8fc25dSMatthew G. Knepley 
8651e8fc25dSMatthew G. Knepley   Collective on DM
8661e8fc25dSMatthew G. Knepley 
8671e8fc25dSMatthew G. Knepley   Input Parameters:
8681e8fc25dSMatthew G. Knepley + dm - The DMPlex object
8691e8fc25dSMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
8701e8fc25dSMatthew G. Knepley . originalSection - The PetscSection for existing data layout
8711e8fc25dSMatthew G. Knepley - originalIS - The existing data
8721e8fc25dSMatthew G. Knepley 
8731e8fc25dSMatthew G. Knepley   Output Parameters:
8741e8fc25dSMatthew G. Knepley + newSection - The PetscSF describing the new data layout
8751e8fc25dSMatthew G. Knepley - newIS - The new data
8761e8fc25dSMatthew G. Knepley 
8771e8fc25dSMatthew G. Knepley   Level: developer
8781e8fc25dSMatthew G. Knepley 
8791e8fc25dSMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexDistributeData()
8801e8fc25dSMatthew G. Knepley @*/
8811e8fc25dSMatthew G. Knepley PetscErrorCode DMPlexDistributeFieldIS(DM dm, PetscSF pointSF, PetscSection originalSection, IS originalIS, PetscSection newSection, IS *newIS)
8821e8fc25dSMatthew G. Knepley {
8831e8fc25dSMatthew G. Knepley   PetscSF         fieldSF;
8841e8fc25dSMatthew G. Knepley   PetscInt       *newValues, *remoteOffsets, fieldSize;
8851e8fc25dSMatthew G. Knepley   const PetscInt *originalValues;
8861e8fc25dSMatthew G. Knepley   PetscErrorCode  ierr;
8871e8fc25dSMatthew G. Knepley 
8881e8fc25dSMatthew G. Knepley   PetscFunctionBegin;
8891e8fc25dSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
8901e8fc25dSMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
8911e8fc25dSMatthew G. Knepley 
8921e8fc25dSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
8931e8fc25dSMatthew G. Knepley   ierr = PetscMalloc(fieldSize * sizeof(PetscInt), &newValues);CHKERRQ(ierr);
8941e8fc25dSMatthew G. Knepley 
8951e8fc25dSMatthew G. Knepley   ierr = ISGetIndices(originalIS, &originalValues);CHKERRQ(ierr);
8961e8fc25dSMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
8971e8fc25dSMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, MPIU_INT, originalValues, newValues);CHKERRQ(ierr);
8981e8fc25dSMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, MPIU_INT, originalValues, newValues);CHKERRQ(ierr);
8991e8fc25dSMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
9001e8fc25dSMatthew G. Knepley   ierr = ISRestoreIndices(originalIS, &originalValues);CHKERRQ(ierr);
9011e8fc25dSMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) pointSF), fieldSize, newValues, PETSC_OWN_POINTER, newIS);CHKERRQ(ierr);
9021e8fc25dSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
9031e8fc25dSMatthew G. Knepley   PetscFunctionReturn(0);
9041e8fc25dSMatthew G. Knepley }
9051e8fc25dSMatthew G. Knepley 
9061e8fc25dSMatthew G. Knepley #undef __FUNCT__
90770034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeData"
90870034214SMatthew G. Knepley /*@
90970034214SMatthew G. Knepley   DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
91070034214SMatthew G. Knepley 
91170034214SMatthew G. Knepley   Collective on DM
91270034214SMatthew G. Knepley 
91370034214SMatthew G. Knepley   Input Parameters:
91470034214SMatthew G. Knepley + dm - The DMPlex object
91570034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
91670034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
91770034214SMatthew G. Knepley . datatype - The type of data
91870034214SMatthew G. Knepley - originalData - The existing data
91970034214SMatthew G. Knepley 
92070034214SMatthew G. Knepley   Output Parameters:
9211e8fc25dSMatthew G. Knepley + newSection - The PetscSection describing the new data layout
92270034214SMatthew G. Knepley - newData - The new data
92370034214SMatthew G. Knepley 
92470034214SMatthew G. Knepley   Level: developer
92570034214SMatthew G. Knepley 
92670034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField()
92770034214SMatthew G. Knepley @*/
92870034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData)
92970034214SMatthew G. Knepley {
93070034214SMatthew G. Knepley   PetscSF        fieldSF;
93170034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
93270034214SMatthew G. Knepley   PetscMPIInt    dataSize;
93370034214SMatthew G. Knepley   PetscErrorCode ierr;
93470034214SMatthew G. Knepley 
93570034214SMatthew G. Knepley   PetscFunctionBegin;
93670034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
93770034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
93870034214SMatthew G. Knepley 
93970034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
94070034214SMatthew G. Knepley   ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr);
94170034214SMatthew G. Knepley   ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr);
94270034214SMatthew G. Knepley 
94370034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
94470034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
94570034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
94670034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
94770034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
94870034214SMatthew G. Knepley   PetscFunctionReturn(0);
94970034214SMatthew G. Knepley }
95070034214SMatthew G. Knepley 
95170034214SMatthew G. Knepley #undef __FUNCT__
952cc71bff1SMichael Lange #define __FUNCT__ "DMPlexDistributeCones"
953cc71bff1SMichael Lange PetscErrorCode DMPlexDistributeCones(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
954cc71bff1SMichael Lange {
955cc71bff1SMichael Lange   DM_Plex               *pmesh   = (DM_Plex*) (dmParallel)->data;
956cc71bff1SMichael Lange   MPI_Comm               comm;
957cc71bff1SMichael Lange   PetscSF                coneSF;
958cc71bff1SMichael Lange   PetscSection           originalConeSection, newConeSection;
959cc71bff1SMichael Lange   PetscInt              *remoteOffsets, *cones, *newCones, newConesSize;
960cc71bff1SMichael Lange   PetscBool              flg;
961cc71bff1SMichael Lange   PetscErrorCode         ierr;
962cc71bff1SMichael Lange 
963cc71bff1SMichael Lange   PetscFunctionBegin;
964cc71bff1SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
965cc71bff1SMichael Lange   PetscValidPointer(dmParallel,4);
966cc71bff1SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
967cc71bff1SMichael Lange 
968cc71bff1SMichael Lange   /* Distribute cone section */
969cc71bff1SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
970cc71bff1SMichael Lange   ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr);
971cc71bff1SMichael Lange   ierr = DMPlexGetConeSection(dmParallel, &newConeSection);CHKERRQ(ierr);
972cc71bff1SMichael Lange   ierr = PetscSFDistributeSection(migrationSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr);
973cc71bff1SMichael Lange   ierr = DMSetUp(dmParallel);CHKERRQ(ierr);
974cc71bff1SMichael Lange   {
975cc71bff1SMichael Lange     PetscInt pStart, pEnd, p;
976cc71bff1SMichael Lange 
977cc71bff1SMichael Lange     ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr);
978cc71bff1SMichael Lange     for (p = pStart; p < pEnd; ++p) {
979cc71bff1SMichael Lange       PetscInt coneSize;
980cc71bff1SMichael Lange       ierr               = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr);
981cc71bff1SMichael Lange       pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize);
982cc71bff1SMichael Lange     }
983cc71bff1SMichael Lange   }
984cc71bff1SMichael Lange   /* Communicate and renumber cones */
985cc71bff1SMichael Lange   ierr = PetscSFCreateSectionSF(migrationSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr);
986cc71bff1SMichael Lange   ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
987cc71bff1SMichael Lange   ierr = DMPlexGetCones(dmParallel, &newCones);CHKERRQ(ierr);
988cc71bff1SMichael Lange   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
989cc71bff1SMichael Lange   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
990cc71bff1SMichael Lange   ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr);
991cc71bff1SMichael Lange   ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr);
992cc71bff1SMichael Lange   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr);
993cc71bff1SMichael Lange   if (flg) {
994cc71bff1SMichael Lange     ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr);
995cc71bff1SMichael Lange     ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
996cc71bff1SMichael Lange     ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr);
997cc71bff1SMichael Lange     ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
998cc71bff1SMichael Lange     ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr);
999cc71bff1SMichael Lange   }
1000cc71bff1SMichael Lange   ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr);
1001cc71bff1SMichael Lange   ierr = DMPlexGetConeOrientations(dmParallel, &newCones);CHKERRQ(ierr);
1002cc71bff1SMichael Lange   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
1003cc71bff1SMichael Lange   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
1004cc71bff1SMichael Lange   ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr);
1005cc71bff1SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
1006cc71bff1SMichael Lange   /* Create supports and stratify sieve */
1007cc71bff1SMichael Lange   {
1008cc71bff1SMichael Lange     PetscInt pStart, pEnd;
1009cc71bff1SMichael Lange 
1010cc71bff1SMichael Lange     ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr);
1011cc71bff1SMichael Lange     ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr);
1012cc71bff1SMichael Lange   }
1013cc71bff1SMichael Lange   ierr = DMPlexSymmetrize(dmParallel);CHKERRQ(ierr);
1014cc71bff1SMichael Lange   ierr = DMPlexStratify(dmParallel);CHKERRQ(ierr);
1015cc71bff1SMichael Lange   PetscFunctionReturn(0);
1016cc71bff1SMichael Lange }
1017cc71bff1SMichael Lange 
1018cc71bff1SMichael Lange #undef __FUNCT__
10190df0e737SMichael Lange #define __FUNCT__ "DMPlexDistributeCoordinates"
10200df0e737SMichael Lange PetscErrorCode DMPlexDistributeCoordinates(DM dm, PetscSF migrationSF, DM dmParallel)
10210df0e737SMichael Lange {
10220df0e737SMichael Lange   MPI_Comm         comm;
10230df0e737SMichael Lange   PetscSection     originalCoordSection, newCoordSection;
10240df0e737SMichael Lange   Vec              originalCoordinates, newCoordinates;
10250df0e737SMichael Lange   PetscInt         bs;
10260df0e737SMichael Lange   const char      *name;
10270df0e737SMichael Lange   const PetscReal *maxCell, *L;
10280df0e737SMichael Lange   PetscErrorCode   ierr;
10290df0e737SMichael Lange 
10300df0e737SMichael Lange   PetscFunctionBegin;
10310df0e737SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10320df0e737SMichael Lange   PetscValidPointer(dmParallel, 3);
10330df0e737SMichael Lange 
10340df0e737SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
10350df0e737SMichael Lange   ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr);
10360df0e737SMichael Lange   ierr = DMGetCoordinateSection(dmParallel, &newCoordSection);CHKERRQ(ierr);
10370df0e737SMichael Lange   ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr);
10380df0e737SMichael Lange   if (originalCoordinates) {
10390df0e737SMichael Lange     ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr);
10400df0e737SMichael Lange     ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr);
10410df0e737SMichael Lange     ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr);
10420df0e737SMichael Lange 
10430df0e737SMichael Lange     ierr = DMPlexDistributeField(dm, migrationSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr);
10440df0e737SMichael Lange     ierr = DMSetCoordinatesLocal(dmParallel, newCoordinates);CHKERRQ(ierr);
10450df0e737SMichael Lange     ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr);
10460df0e737SMichael Lange     ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr);
10470df0e737SMichael Lange     ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr);
10480df0e737SMichael Lange   }
10490df0e737SMichael Lange   ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr);
10500df0e737SMichael Lange   if (L) {ierr = DMSetPeriodicity(dmParallel, maxCell, L);CHKERRQ(ierr);}
10510df0e737SMichael Lange   PetscFunctionReturn(0);
10520df0e737SMichael Lange }
10530df0e737SMichael Lange 
10540df0e737SMichael Lange #undef __FUNCT__
10550df0e737SMichael Lange #define __FUNCT__ "DMPlexDistributeLabels"
10562eb1fa14SMichael Lange PetscErrorCode DMPlexDistributeLabels(DM dm, PetscSF migrationSF, DM dmParallel)
10570df0e737SMichael Lange {
10580df0e737SMichael Lange   MPI_Comm       comm;
10590df0e737SMichael Lange   PetscMPIInt    rank;
1060*bdd2d751SMichael Lange   PetscInt       numLabels, l;
10610df0e737SMichael Lange   PetscErrorCode ierr;
10620df0e737SMichael Lange 
10630df0e737SMichael Lange   PetscFunctionBegin;
10640df0e737SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10652eb1fa14SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
10660df0e737SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
10670df0e737SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
10680df0e737SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
10690df0e737SMichael Lange 
10700df0e737SMichael Lange   /* Bcast number of labels */
1071*bdd2d751SMichael Lange   ierr = DMPlexGetNumLabels(dm, &numLabels);CHKERRQ(ierr);
10720df0e737SMichael Lange   ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1073*bdd2d751SMichael Lange   for (l = numLabels-1; l >= 0; --l) {
1074*bdd2d751SMichael Lange     DMLabel     label = NULL, labelNew = NULL;
10750df0e737SMichael Lange     PetscBool   isdepth;
10760df0e737SMichael Lange 
1077*bdd2d751SMichael Lange     if (!rank) {
1078*bdd2d751SMichael Lange       ierr = DMPlexGetLabelByNum(dm, l, &label);CHKERRQ(ierr);
10790df0e737SMichael Lange       /* Skip "depth" because it is recreated */
1080*bdd2d751SMichael Lange       ierr = PetscStrcmp(label->name, "depth", &isdepth);CHKERRQ(ierr);
1081*bdd2d751SMichael Lange     }
10820df0e737SMichael Lange     ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
1083*bdd2d751SMichael Lange     if (isdepth) continue;
1084*bdd2d751SMichael Lange     ierr = DMLabelDistribute(label, migrationSF, &labelNew);CHKERRQ(ierr);
1085*bdd2d751SMichael Lange     ierr = DMPlexAddLabel(dmParallel, labelNew);CHKERRQ(ierr);
10860df0e737SMichael Lange   }
10870df0e737SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
10880df0e737SMichael Lange   PetscFunctionReturn(0);
10890df0e737SMichael Lange }
10900df0e737SMichael Lange 
10919734c634SMichael Lange #undef __FUNCT__
10929734c634SMichael Lange #define __FUNCT__ "DMPlexDistributeSetupHybrid"
10939734c634SMichael Lange PetscErrorCode DMPlexDistributeSetupHybrid(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
10949734c634SMichael Lange {
10959734c634SMichael Lange   DM_Plex        *mesh  = (DM_Plex*) dm->data;
10969734c634SMichael Lange   DM_Plex        *pmesh = (DM_Plex*) (dmParallel)->data;
10979734c634SMichael Lange   MPI_Comm        comm;
10989734c634SMichael Lange   const PetscInt *gpoints;
10999734c634SMichael Lange   PetscInt        dim, depth, n, d;
11009734c634SMichael Lange   PetscErrorCode  ierr;
11019734c634SMichael Lange 
11029734c634SMichael Lange   PetscFunctionBegin;
11039734c634SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11049734c634SMichael Lange   PetscValidPointer(dmParallel, 4);
11059734c634SMichael Lange 
11069734c634SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
11079734c634SMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
11089734c634SMichael Lange 
11099734c634SMichael Lange   /* Setup hybrid structure */
11109734c634SMichael Lange   for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];}
11119734c634SMichael Lange   ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr);
11129734c634SMichael Lange   ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr);
11139734c634SMichael Lange   ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr);
11149734c634SMichael Lange   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
11159734c634SMichael Lange   for (d = 0; d <= dim; ++d) {
11169734c634SMichael Lange     PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p;
11179734c634SMichael Lange 
11189734c634SMichael Lange     if (pmax < 0) continue;
11199734c634SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr);
11209734c634SMichael Lange     ierr = DMPlexGetDepthStratum(dmParallel, d, NULL, &pEnd);CHKERRQ(ierr);
11219734c634SMichael Lange     ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr);
11229734c634SMichael Lange     for (p = 0; p < n; ++p) {
11239734c634SMichael Lange       const PetscInt point = gpoints[p];
11249734c634SMichael Lange 
11259734c634SMichael Lange       if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax;
11269734c634SMichael Lange     }
11279734c634SMichael Lange     if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax;
11289734c634SMichael Lange     else            pmesh->hybridPointMax[d] = -1;
11299734c634SMichael Lange   }
11309734c634SMichael Lange   ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr);
11319734c634SMichael Lange   PetscFunctionReturn(0);
11329734c634SMichael Lange }
11330df0e737SMichael Lange 
1134a6f36705SMichael Lange #undef __FUNCT__
1135a6f36705SMichael Lange #define __FUNCT__ "DMPlexDistributeSetupTree"
1136a6f36705SMichael Lange PetscErrorCode DMPlexDistributeSetupTree(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
1137a6f36705SMichael Lange {
1138a6f36705SMichael Lange   MPI_Comm        comm;
1139a6f36705SMichael Lange   DM              refTree;
1140a6f36705SMichael Lange   PetscSection    origParentSection, newParentSection;
1141a6f36705SMichael Lange   PetscInt        *origParents, *origChildIDs;
1142a6f36705SMichael Lange   PetscBool       flg;
1143a6f36705SMichael Lange   PetscErrorCode  ierr;
1144a6f36705SMichael Lange 
1145a6f36705SMichael Lange   PetscFunctionBegin;
1146a6f36705SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1147a6f36705SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 4);
1148a6f36705SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
1149a6f36705SMichael Lange 
1150a6f36705SMichael Lange   /* Set up tree */
1151a6f36705SMichael Lange   ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr);
1152a6f36705SMichael Lange   ierr = DMPlexSetReferenceTree(dmParallel,refTree);CHKERRQ(ierr);
1153a6f36705SMichael Lange   ierr = DMPlexGetTree(dm,&origParentSection,&origParents,&origChildIDs,NULL,NULL);CHKERRQ(ierr);
1154a6f36705SMichael Lange   if (origParentSection) {
1155a6f36705SMichael Lange     PetscInt        pStart, pEnd;
1156a6f36705SMichael Lange     PetscInt        *newParents, *newChildIDs;
1157a6f36705SMichael Lange     PetscInt        *remoteOffsetsParents, newParentSize;
1158a6f36705SMichael Lange     PetscSF         parentSF;
1159a6f36705SMichael Lange 
1160a6f36705SMichael Lange     ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
1161a6f36705SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dmParallel),&newParentSection);CHKERRQ(ierr);
1162a6f36705SMichael Lange     ierr = PetscSectionSetChart(newParentSection,pStart,pEnd);CHKERRQ(ierr);
1163a6f36705SMichael Lange     ierr = PetscSFDistributeSection(migrationSF, origParentSection, &remoteOffsetsParents, newParentSection);CHKERRQ(ierr);
1164a6f36705SMichael Lange     ierr = PetscSFCreateSectionSF(migrationSF, origParentSection, remoteOffsetsParents, newParentSection, &parentSF);CHKERRQ(ierr);
1165a6f36705SMichael Lange     ierr = PetscSectionGetStorageSize(newParentSection,&newParentSize);CHKERRQ(ierr);
1166a6f36705SMichael Lange     ierr = PetscMalloc2(newParentSize,&newParents,newParentSize,&newChildIDs);CHKERRQ(ierr);
1167a6f36705SMichael Lange     ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
1168a6f36705SMichael Lange     ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
1169a6f36705SMichael Lange     ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1170a6f36705SMichael Lange     ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1171a6f36705SMichael Lange     ierr = ISGlobalToLocalMappingApplyBlock(renumbering,IS_GTOLM_MASK, newParentSize, newParents, NULL, newParents);CHKERRQ(ierr);
1172a6f36705SMichael Lange     ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-parents_view", &flg);CHKERRQ(ierr);
1173a6f36705SMichael Lange     if (flg) {
1174a6f36705SMichael Lange       ierr = PetscPrintf(comm, "Serial Parent Section: \n");CHKERRQ(ierr);
1175a6f36705SMichael Lange       ierr = PetscSectionView(origParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1176a6f36705SMichael Lange       ierr = PetscPrintf(comm, "Parallel Parent Section: \n");CHKERRQ(ierr);
1177a6f36705SMichael Lange       ierr = PetscSectionView(newParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1178a6f36705SMichael Lange       ierr = PetscSFView(parentSF, NULL);CHKERRQ(ierr);
1179a6f36705SMichael Lange     }
1180a6f36705SMichael Lange     ierr = DMPlexSetTree(dmParallel,newParentSection,newParents,newChildIDs);CHKERRQ(ierr);
1181a6f36705SMichael Lange     ierr = PetscSectionDestroy(&newParentSection);CHKERRQ(ierr);
1182a6f36705SMichael Lange     ierr = PetscFree2(newParents,newChildIDs);CHKERRQ(ierr);
1183a6f36705SMichael Lange     ierr = PetscSFDestroy(&parentSF);CHKERRQ(ierr);
1184a6f36705SMichael Lange   }
1185a6f36705SMichael Lange   PetscFunctionReturn(0);
1186a6f36705SMichael Lange }
11870df0e737SMichael Lange 
11880df0e737SMichael Lange #undef __FUNCT__
11894eca1733SMichael Lange #define __FUNCT__ "DMPlexDistributeSF"
11904eca1733SMichael Lange PetscErrorCode DMPlexDistributeSF(DM dm, PetscSF migrationSF, PetscSection partSection, IS part, PetscSection origPartSection, IS origPart, DM dmParallel)
11914eca1733SMichael Lange {
11924eca1733SMichael Lange   DM_Plex               *mesh  = (DM_Plex*) dm->data;
11934eca1733SMichael Lange   DM_Plex               *pmesh = (DM_Plex*) (dmParallel)->data;
11944eca1733SMichael Lange   PetscMPIInt            rank, numProcs;
11954eca1733SMichael Lange   MPI_Comm               comm;
11964eca1733SMichael Lange   PetscErrorCode         ierr;
11974eca1733SMichael Lange 
11984eca1733SMichael Lange   PetscFunctionBegin;
11994eca1733SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
12004eca1733SMichael Lange   PetscValidPointer(dmParallel,7);
12014eca1733SMichael Lange 
12024eca1733SMichael Lange   /* Create point SF for parallel mesh */
12034eca1733SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
12044eca1733SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
12054eca1733SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
12064eca1733SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
12074eca1733SMichael Lange   {
12084eca1733SMichael Lange     const PetscInt *leaves;
12094eca1733SMichael Lange     PetscSFNode    *remotePoints, *rowners, *lowners;
12104eca1733SMichael Lange     PetscInt        numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints;
12114eca1733SMichael Lange     PetscInt        pStart, pEnd;
12124eca1733SMichael Lange 
12134eca1733SMichael Lange     ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
12144eca1733SMichael Lange     ierr = PetscSFGetGraph(migrationSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr);
12154eca1733SMichael Lange     ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr);
12164eca1733SMichael Lange     for (p=0; p<numRoots; p++) {
12174eca1733SMichael Lange       rowners[p].rank  = -1;
12184eca1733SMichael Lange       rowners[p].index = -1;
12194eca1733SMichael Lange     }
12204eca1733SMichael Lange     if (origPart) {
12214eca1733SMichael Lange       /* Make sure points in the original partition are not assigned to other procs */
12224eca1733SMichael Lange       const PetscInt *origPoints;
12234eca1733SMichael Lange 
12244eca1733SMichael Lange       ierr = ISGetIndices(origPart, &origPoints);CHKERRQ(ierr);
12254eca1733SMichael Lange       for (p = 0; p < numProcs; ++p) {
12264eca1733SMichael Lange         PetscInt dof, off, d;
12274eca1733SMichael Lange 
12284eca1733SMichael Lange         ierr = PetscSectionGetDof(origPartSection, p, &dof);CHKERRQ(ierr);
12294eca1733SMichael Lange         ierr = PetscSectionGetOffset(origPartSection, p, &off);CHKERRQ(ierr);
12304eca1733SMichael Lange         for (d = off; d < off+dof; ++d) {
12314eca1733SMichael Lange           rowners[origPoints[d]].rank = p;
12324eca1733SMichael Lange         }
12334eca1733SMichael Lange       }
12344eca1733SMichael Lange       ierr = ISRestoreIndices(origPart, &origPoints);CHKERRQ(ierr);
12354eca1733SMichael Lange     }
12364eca1733SMichael Lange     ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12374eca1733SMichael Lange     ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12384eca1733SMichael Lange     for (p = 0; p < numLeaves; ++p) {
12394eca1733SMichael Lange       if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */
12404eca1733SMichael Lange         lowners[p].rank  = rank;
12414eca1733SMichael Lange         lowners[p].index = leaves ? leaves[p] : p;
12424eca1733SMichael Lange       } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */
12434eca1733SMichael Lange         lowners[p].rank  = -2;
12444eca1733SMichael Lange         lowners[p].index = -2;
12454eca1733SMichael Lange       }
12464eca1733SMichael Lange     }
12474eca1733SMichael Lange     for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */
12484eca1733SMichael Lange       rowners[p].rank  = -3;
12494eca1733SMichael Lange       rowners[p].index = -3;
12504eca1733SMichael Lange     }
12514eca1733SMichael Lange     ierr = PetscSFReduceBegin(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
12524eca1733SMichael Lange     ierr = PetscSFReduceEnd(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
12534eca1733SMichael Lange     ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12544eca1733SMichael Lange     ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
12554eca1733SMichael Lange     for (p = 0; p < numLeaves; ++p) {
12564eca1733SMichael Lange       if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed");
12574eca1733SMichael Lange       if (lowners[p].rank != rank) ++numGhostPoints;
12584eca1733SMichael Lange     }
12594eca1733SMichael Lange     ierr = PetscMalloc1(numGhostPoints,    &ghostPoints);CHKERRQ(ierr);
12604eca1733SMichael Lange     ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr);
12614eca1733SMichael Lange     for (p = 0, gp = 0; p < numLeaves; ++p) {
12624eca1733SMichael Lange       if (lowners[p].rank != rank) {
12634eca1733SMichael Lange         ghostPoints[gp]        = leaves ? leaves[p] : p;
12644eca1733SMichael Lange         remotePoints[gp].rank  = lowners[p].rank;
12654eca1733SMichael Lange         remotePoints[gp].index = lowners[p].index;
12664eca1733SMichael Lange         ++gp;
12674eca1733SMichael Lange       }
12684eca1733SMichael Lange     }
12694eca1733SMichael Lange     ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr);
12704eca1733SMichael Lange     ierr = PetscSFSetGraph((dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
12714eca1733SMichael Lange     ierr = PetscSFSetFromOptions((dmParallel)->sf);CHKERRQ(ierr);
12724eca1733SMichael Lange   }
12734eca1733SMichael Lange   pmesh->useCone    = mesh->useCone;
12744eca1733SMichael Lange   pmesh->useClosure = mesh->useClosure;
12754eca1733SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
12764eca1733SMichael Lange   PetscFunctionReturn(0);
12774eca1733SMichael Lange }
12784eca1733SMichael Lange 
12794eca1733SMichael Lange #undef __FUNCT__
128070034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistribute"
128170034214SMatthew G. Knepley /*@C
128270034214SMatthew G. Knepley   DMPlexDistribute - Distributes the mesh and any associated sections.
128370034214SMatthew G. Knepley 
128470034214SMatthew G. Knepley   Not Collective
128570034214SMatthew G. Knepley 
128670034214SMatthew G. Knepley   Input Parameter:
128770034214SMatthew G. Knepley + dm  - The original DMPlex object
128870034214SMatthew G. Knepley . partitioner - The partitioning package, or NULL for the default
128970034214SMatthew G. Knepley - overlap - The overlap of partitions, 0 is the default
129070034214SMatthew G. Knepley 
129170034214SMatthew G. Knepley   Output Parameter:
129270034214SMatthew G. Knepley + sf - The PetscSF used for point distribution
129370034214SMatthew G. Knepley - parallelMesh - The distributed DMPlex object, or NULL
129470034214SMatthew G. Knepley 
1295a9c22940SMatthew G. Knepley   Note: If the mesh was not distributed, the return value is NULL.
1296a9c22940SMatthew G. Knepley 
1297a9c22940SMatthew G. Knepley   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
1298a9c22940SMatthew G. Knepley   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
1299a9c22940SMatthew G. Knepley   representation on the mesh.
130070034214SMatthew G. Knepley 
130170034214SMatthew G. Knepley   Level: intermediate
130270034214SMatthew G. Knepley 
130370034214SMatthew G. Knepley .keywords: mesh, elements
1304a9c22940SMatthew G. Knepley .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
130570034214SMatthew G. Knepley @*/
130670034214SMatthew G. Knepley PetscErrorCode DMPlexDistribute(DM dm, const char partitioner[], PetscInt overlap, PetscSF *sf, DM *dmParallel)
130770034214SMatthew G. Knepley {
130870034214SMatthew G. Knepley   DM_Plex               *mesh   = (DM_Plex*) dm->data, *pmesh;
130970034214SMatthew G. Knepley   MPI_Comm               comm;
131043331d4aSMichael Lange   PetscInt               dim, numRemoteRanks, nroots, nleaves;
1311a157612eSMichael Lange   DM                     dmOverlap;
1312a157612eSMichael Lange   IS                     cellPart,        part;
1313a157612eSMichael Lange   PetscSection           cellPartSection, partSection;
131443331d4aSMichael Lange   PetscSFNode           *remoteRanks, *newRemote;
131543331d4aSMichael Lange   const PetscSFNode     *oldRemote;
131643331d4aSMichael Lange   PetscSF                partSF, pointSF, overlapPointSF, overlapSF;
131770034214SMatthew G. Knepley   ISLocalToGlobalMapping renumbering;
131870034214SMatthew G. Knepley   PetscBool              flg;
131970034214SMatthew G. Knepley   PetscMPIInt            rank, numProcs, p;
132070034214SMatthew G. Knepley   PetscErrorCode         ierr;
132170034214SMatthew G. Knepley 
132270034214SMatthew G. Knepley   PetscFunctionBegin;
132370034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
132470034214SMatthew G. Knepley   if (sf) PetscValidPointer(sf,4);
132570034214SMatthew G. Knepley   PetscValidPointer(dmParallel,5);
132670034214SMatthew G. Knepley 
132770034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
132870034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
132970034214SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
133070034214SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
133170034214SMatthew G. Knepley 
133270034214SMatthew G. Knepley   *dmParallel = NULL;
133370034214SMatthew G. Knepley   if (numProcs == 1) PetscFunctionReturn(0);
133470034214SMatthew G. Knepley 
1335c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
133670034214SMatthew G. Knepley   /* Create cell partition - We need to rewrite to use IS, use the MatPartition stuff */
133777623264SMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
133877623264SMatthew G. Knepley   if (overlap > 1) SETERRQ(comm, PETSC_ERR_SUP, "Overlap > 1 not yet implemented");
133977623264SMatthew G. Knepley   ierr = PetscSectionCreate(comm, &cellPartSection);CHKERRQ(ierr);
1340a157612eSMichael Lange   ierr = PetscPartitionerPartition(mesh->partitioner, dm, PETSC_FALSE, cellPartSection, &cellPart, NULL, NULL);CHKERRQ(ierr);
134170034214SMatthew G. Knepley   /* Create SF assuming a serial partition for all processes: Could check for IS length here */
134270034214SMatthew G. Knepley   if (!rank) numRemoteRanks = numProcs;
134370034214SMatthew G. Knepley   else       numRemoteRanks = 0;
134470034214SMatthew G. Knepley   ierr = PetscMalloc1(numRemoteRanks, &remoteRanks);CHKERRQ(ierr);
134570034214SMatthew G. Knepley   for (p = 0; p < numRemoteRanks; ++p) {
134670034214SMatthew G. Knepley     remoteRanks[p].rank  = p;
134770034214SMatthew G. Knepley     remoteRanks[p].index = 0;
134870034214SMatthew G. Knepley   }
134970034214SMatthew G. Knepley   ierr = PetscSFCreate(comm, &partSF);CHKERRQ(ierr);
135070034214SMatthew G. Knepley   ierr = PetscSFSetGraph(partSF, 1, numRemoteRanks, NULL, PETSC_OWN_POINTER, remoteRanks, PETSC_OWN_POINTER);CHKERRQ(ierr);
135170034214SMatthew G. Knepley   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
135270034214SMatthew G. Knepley   if (flg) {
1353a157612eSMichael Lange     ierr = PetscPrintf(comm, "Original Cell Partition:\n");CHKERRQ(ierr);
135470034214SMatthew G. Knepley     ierr = PetscSectionView(cellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
135570034214SMatthew G. Knepley     ierr = ISView(cellPart, NULL);CHKERRQ(ierr);
135670034214SMatthew G. Knepley     ierr = PetscSFView(partSF, NULL);CHKERRQ(ierr);
135770034214SMatthew G. Knepley   }
135870034214SMatthew G. Knepley   /* Close the partition over the mesh */
135970034214SMatthew G. Knepley   ierr = DMPlexCreatePartitionClosure(dm, cellPartSection, cellPart, &partSection, &part);CHKERRQ(ierr);
136070034214SMatthew G. Knepley   /* Create new mesh */
136170034214SMatthew G. Knepley   ierr  = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr);
1362c73cfb54SMatthew G. Knepley   ierr  = DMSetDimension(*dmParallel, dim);CHKERRQ(ierr);
136370034214SMatthew G. Knepley   ierr  = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr);
136470034214SMatthew G. Knepley   pmesh = (DM_Plex*) (*dmParallel)->data;
13654eca1733SMichael Lange   pmesh->useAnchors = mesh->useAnchors;
13664eca1733SMichael Lange 
136770034214SMatthew G. Knepley   /* Distribute sieve points and the global point numbering (replaces creating remote bases) */
136870034214SMatthew G. Knepley   ierr = PetscSFConvertPartition(partSF, partSection, part, &renumbering, &pointSF);CHKERRQ(ierr);
136970034214SMatthew G. Knepley   if (flg) {
137070034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Point Partition:\n");CHKERRQ(ierr);
137170034214SMatthew G. Knepley     ierr = PetscSectionView(partSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
137270034214SMatthew G. Knepley     ierr = ISView(part, NULL);CHKERRQ(ierr);
137370034214SMatthew G. Knepley     ierr = PetscSFView(pointSF, NULL);CHKERRQ(ierr);
137470034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Point Renumbering after partition:\n");CHKERRQ(ierr);
137570034214SMatthew G. Knepley     ierr = ISLocalToGlobalMappingView(renumbering, NULL);CHKERRQ(ierr);
137670034214SMatthew G. Knepley   }
137777623264SMatthew G. Knepley   ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
137870034214SMatthew G. Knepley 
1379a157612eSMichael Lange   /* Migrate data to a non-overlapping parallel DM */
1380cc71bff1SMichael Lange   ierr = DMPlexDistributeCones(dm, pointSF, renumbering, *dmParallel);CHKERRQ(ierr);
13810df0e737SMichael Lange   ierr = DMPlexDistributeCoordinates(dm, pointSF, *dmParallel);CHKERRQ(ierr);
13822eb1fa14SMichael Lange   ierr = DMPlexDistributeLabels(dm, pointSF, *dmParallel);CHKERRQ(ierr);
13839734c634SMichael Lange   ierr = DMPlexDistributeSetupHybrid(dm, pointSF, renumbering, *dmParallel);CHKERRQ(ierr);
1384a6f36705SMichael Lange   ierr = DMPlexDistributeSetupTree(dm, pointSF, renumbering, *dmParallel);CHKERRQ(ierr);
138570034214SMatthew G. Knepley 
1386a157612eSMichael Lange   /* Build the point SF without overlap */
1387a157612eSMichael Lange   ierr = DMPlexDistributeSF(dm, pointSF, partSection, part, NULL, NULL, *dmParallel);CHKERRQ(ierr);
1388a157612eSMichael Lange 
1389a157612eSMichael Lange   if (overlap > 0) {
13903d822a50SMichael Lange     ierr = PetscLogEventBegin(DMPLEX_DistributeOverlap,dm,0,0,0);CHKERRQ(ierr);
1391a157612eSMichael Lange     /* Add the partition overlap to the distributed DM */
1392a157612eSMichael Lange     ierr = DMPlexDistributeOverlap(*dmParallel, overlap, renumbering, &overlapSF, &dmOverlap);CHKERRQ(ierr);
1393a157612eSMichael Lange     ierr = DMDestroy(dmParallel);CHKERRQ(ierr);
1394a157612eSMichael Lange     *dmParallel = dmOverlap;
13953d822a50SMichael Lange     if (flg) {
13963d822a50SMichael Lange       ierr = PetscPrintf(comm, "Overlap Migration SF:\n");CHKERRQ(ierr);
13973d822a50SMichael Lange       ierr = PetscSFView(overlapSF, NULL);CHKERRQ(ierr);
13983d822a50SMichael Lange     }
139943331d4aSMichael Lange 
140043331d4aSMichael Lange     /* Re-map the pointSF to establish the full migration pattern */
140143331d4aSMichael Lange     ierr = PetscSFGetGraph(pointSF, &nroots, NULL, NULL, &oldRemote);CHKERRQ(ierr);
140243331d4aSMichael Lange     ierr = PetscSFGetGraph(overlapSF, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr);
140343331d4aSMichael Lange     ierr = PetscMalloc1(nleaves, &newRemote);CHKERRQ(ierr);
140443331d4aSMichael Lange     ierr = PetscSFBcastBegin(overlapSF, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr);
140543331d4aSMichael Lange     ierr = PetscSFBcastEnd(overlapSF, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr);
140643331d4aSMichael Lange     ierr = PetscSFCreate(comm, &overlapPointSF);CHKERRQ(ierr);
140743331d4aSMichael Lange     ierr = PetscSFSetGraph(overlapPointSF, nroots, nleaves, NULL, PETSC_OWN_POINTER, newRemote, PETSC_OWN_POINTER);CHKERRQ(ierr);
140843331d4aSMichael Lange     ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr);
140943331d4aSMichael Lange     pointSF = overlapPointSF;
14103d822a50SMichael Lange     ierr = PetscLogEventEnd(DMPLEX_DistributeOverlap,dm,0,0,0);CHKERRQ(ierr);
1411a157612eSMichael Lange   }
1412bf5244c7SMatthew G. Knepley   /* Cleanup Partition */
1413bf5244c7SMatthew G. Knepley   ierr = ISLocalToGlobalMappingDestroy(&renumbering);CHKERRQ(ierr);
1414bf5244c7SMatthew G. Knepley   ierr = PetscSFDestroy(&partSF);CHKERRQ(ierr);
1415bf5244c7SMatthew G. Knepley   ierr = PetscSectionDestroy(&partSection);CHKERRQ(ierr);
1416bf5244c7SMatthew G. Knepley   ierr = ISDestroy(&part);CHKERRQ(ierr);
14174eca1733SMichael Lange   ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr);
14184eca1733SMichael Lange   ierr = ISDestroy(&cellPart);CHKERRQ(ierr);
141986719b60SMatthew G. Knepley   /* Copy BC */
142086719b60SMatthew G. Knepley   ierr = DMPlexCopyBoundary(dm, *dmParallel);CHKERRQ(ierr);
142170034214SMatthew G. Knepley   /* Cleanup */
142270034214SMatthew G. Knepley   if (sf) {*sf = pointSF;}
142370034214SMatthew G. Knepley   else    {ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr);}
142470034214SMatthew G. Knepley   ierr = DMSetFromOptions(*dmParallel);CHKERRQ(ierr);
142570034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
142670034214SMatthew G. Knepley   PetscFunctionReturn(0);
142770034214SMatthew G. Knepley }
1428a157612eSMichael Lange 
1429a157612eSMichael Lange #undef __FUNCT__
1430a157612eSMichael Lange #define __FUNCT__ "DMPlexDistributeOverlap"
1431a157612eSMichael Lange /*@C
1432a157612eSMichael Lange   DMPlexDistribute - Add partition overlap to a distributed non-overlapping DM.
1433a157612eSMichael Lange 
1434a157612eSMichael Lange   Not Collective
1435a157612eSMichael Lange 
1436a157612eSMichael Lange   Input Parameter:
1437a157612eSMichael Lange + dm  - The non-overlapping distrbuted DMPlex object
1438a157612eSMichael Lange - overlap - The overlap of partitions, 0 is the default
1439a157612eSMichael Lange 
1440a157612eSMichael Lange   Output Parameter:
1441a157612eSMichael Lange + sf - The PetscSF used for point distribution
1442a157612eSMichael Lange - dmOverlap - The overlapping distributed DMPlex object, or NULL
1443a157612eSMichael Lange 
1444a157612eSMichael Lange   Note: If the mesh was not distributed, the return value is NULL.
1445a157612eSMichael Lange 
1446a157612eSMichael Lange   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
1447a157612eSMichael Lange   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
1448a157612eSMichael Lange   representation on the mesh.
1449a157612eSMichael Lange 
1450a157612eSMichael Lange   Level: intermediate
1451a157612eSMichael Lange 
1452a157612eSMichael Lange .keywords: mesh, elements
1453a157612eSMichael Lange .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
1454a157612eSMichael Lange @*/
1455a157612eSMichael Lange PetscErrorCode DMPlexDistributeOverlap(DM dm, PetscInt overlap, ISLocalToGlobalMapping renumbering, PetscSF *sf, DM *dmOverlap)
1456a157612eSMichael Lange {
1457a157612eSMichael Lange   MPI_Comm               comm;
1458a157612eSMichael Lange   PetscMPIInt            rank;
1459a157612eSMichael Lange   PetscSection           rootSection, leafSection;
1460a157612eSMichael Lange   IS                     rootrank, leafrank;
1461a157612eSMichael Lange   PetscSection           coneSection;
1462a157612eSMichael Lange   PetscSF                overlapSF, migrationSF, pointSF, newPointSF;
1463a157612eSMichael Lange   PetscSFNode           *ghostRemote;
1464a157612eSMichael Lange   const PetscSFNode     *overlapRemote;
1465a157612eSMichael Lange   ISLocalToGlobalMapping overlapRenumbering;
1466a157612eSMichael Lange   const PetscInt        *renumberingArray, *overlapLocal;
1467a157612eSMichael Lange   PetscInt               dim, p, pStart, pEnd, conesSize, idx;
1468a157612eSMichael Lange   PetscInt               numGhostPoints, numOverlapPoints, numSharedPoints, overlapLeaves;
1469a157612eSMichael Lange   PetscInt              *cones, *ghostLocal, *overlapRenumberingArray, *pointIDs, *recvPointIDs;
1470a157612eSMichael Lange   PetscErrorCode         ierr;
1471a157612eSMichael Lange 
1472a157612eSMichael Lange   PetscFunctionBegin;
1473a157612eSMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1474a157612eSMichael Lange   if (sf) PetscValidPointer(sf, 3);
1475a157612eSMichael Lange   PetscValidPointer(dmOverlap, 4);
1476a157612eSMichael Lange 
1477a157612eSMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1478a157612eSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1479a157612eSMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
1480a157612eSMichael Lange 
1481a157612eSMichael Lange   /* Compute point overlap with neighbouring processes on the distributed DM */
1482a157612eSMichael Lange   ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1483a157612eSMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
1484a157612eSMichael Lange   ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr);
1485a157612eSMichael Lange   ierr = DMPlexDistributeOwnership(dm, rootSection, &rootrank, leafSection, &leafrank);CHKERRQ(ierr);
1486a157612eSMichael Lange   ierr = DMPlexCreateOverlap(dm, rootSection, rootrank, leafSection, leafrank, &overlapSF);CHKERRQ(ierr);
1487a157612eSMichael Lange   ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1488a157612eSMichael Lange 
1489a157612eSMichael Lange   /* Build dense migration SF that maps the non-overlapping partition to the overlapping one */
1490a157612eSMichael Lange   ierr = DMPlexCreateOverlapMigrationSF(dm, overlapSF, &migrationSF);CHKERRQ(ierr);
1491a157612eSMichael Lange 
1492a157612eSMichael Lange   /* Convert cones to global numbering before migrating them */
1493a157612eSMichael Lange   ierr = DMPlexGetConeSection(dm, &coneSection);CHKERRQ(ierr);
1494a157612eSMichael Lange   ierr = PetscSectionGetStorageSize(coneSection, &conesSize);CHKERRQ(ierr);
1495a157612eSMichael Lange   ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
1496a157612eSMichael Lange   ierr = ISLocalToGlobalMappingApplyBlock(renumbering, conesSize, cones, cones);CHKERRQ(ierr);
1497a157612eSMichael Lange 
1498a157612eSMichael Lange   /* Derive the new local-to-global mapping from the old one */
1499a157612eSMichael Lange   ierr = PetscSFGetGraph(migrationSF, NULL, &overlapLeaves, &overlapLocal, &overlapRemote);CHKERRQ(ierr);
1500a157612eSMichael Lange   ierr = PetscMalloc1(overlapLeaves, &overlapRenumberingArray);CHKERRQ(ierr);
1501a157612eSMichael Lange   ierr = ISLocalToGlobalMappingGetBlockIndices(renumbering, &renumberingArray);CHKERRQ(ierr);
1502a157612eSMichael Lange   ierr = PetscSFBcastBegin(migrationSF, MPIU_INT, renumberingArray, overlapRenumberingArray);CHKERRQ(ierr);
1503a157612eSMichael Lange   ierr = PetscSFBcastEnd(migrationSF, MPIU_INT, renumberingArray, overlapRenumberingArray);CHKERRQ(ierr);
1504a157612eSMichael Lange   ierr = ISLocalToGlobalMappingCreate(comm, 1, overlapLeaves, (const PetscInt*) overlapRenumberingArray, PETSC_OWN_POINTER, &overlapRenumbering);CHKERRQ(ierr);
1505a157612eSMichael Lange 
1506a157612eSMichael Lange   /* Build the overlapping DM */
1507a157612eSMichael Lange   ierr = DMPlexCreate(comm, dmOverlap);CHKERRQ(ierr);
1508a157612eSMichael Lange   ierr = DMSetDimension(*dmOverlap, dim);CHKERRQ(ierr);
1509a157612eSMichael Lange   ierr = PetscObjectSetName((PetscObject) *dmOverlap, "Parallel Mesh");CHKERRQ(ierr);
1510a157612eSMichael Lange   ierr = DMPlexDistributeCones(dm, migrationSF, overlapRenumbering, *dmOverlap);CHKERRQ(ierr);
1511a157612eSMichael Lange   ierr = DMPlexDistributeCoordinates(dm, migrationSF, *dmOverlap);CHKERRQ(ierr);
1512a157612eSMichael Lange   ierr = DMPlexDistributeLabels(dm, migrationSF, *dmOverlap);CHKERRQ(ierr);
1513a157612eSMichael Lange   ierr = DMPlexDistributeSetupHybrid(dm, migrationSF, overlapRenumbering, *dmOverlap);CHKERRQ(ierr);
1514a157612eSMichael Lange 
1515a157612eSMichael Lange   /* Build the new point SF by propagating the depthShift generate remote root indices */
1516a157612eSMichael Lange   ierr = DMGetPointSF(dm, &pointSF);CHKERRQ(ierr);
1517a157612eSMichael Lange   ierr = PetscSFGetGraph(pointSF, NULL, &numSharedPoints, NULL, NULL);CHKERRQ(ierr);
1518a157612eSMichael Lange   ierr = PetscSFGetGraph(overlapSF, NULL, &numOverlapPoints, NULL, NULL);CHKERRQ(ierr);
1519a157612eSMichael Lange   numGhostPoints = numSharedPoints + numOverlapPoints;
1520a157612eSMichael Lange   ierr = PetscMalloc2(numGhostPoints, &ghostLocal, numGhostPoints, &ghostRemote);CHKERRQ(ierr);
1521a157612eSMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
1522a157612eSMichael Lange   ierr = PetscMalloc2(pEnd-pStart, &pointIDs, overlapLeaves, &recvPointIDs);CHKERRQ(ierr);
1523a157612eSMichael Lange   for (p=0; p<overlapLeaves; p++) {
1524a157612eSMichael Lange     if (overlapRemote[p].rank == rank) pointIDs[overlapRemote[p].index] = overlapLocal[p];
1525a157612eSMichael Lange   }
1526a157612eSMichael Lange   ierr = PetscSFBcastBegin(migrationSF, MPIU_INT, pointIDs, recvPointIDs);CHKERRQ(ierr);
1527a157612eSMichael Lange   ierr = PetscSFBcastEnd(migrationSF, MPIU_INT, pointIDs, recvPointIDs);CHKERRQ(ierr);
1528a157612eSMichael Lange   for (idx=0, p=0; p<overlapLeaves; p++) {
1529a157612eSMichael Lange     if (overlapRemote[p].rank != rank) {
1530a157612eSMichael Lange       ghostLocal[idx] = overlapLocal[p];
1531a157612eSMichael Lange       ghostRemote[idx].index = recvPointIDs[p];
1532a157612eSMichael Lange       ghostRemote[idx].rank = overlapRemote[p].rank;
1533a157612eSMichael Lange       idx++;
1534a157612eSMichael Lange     }
1535a157612eSMichael Lange   }
1536a157612eSMichael Lange   ierr = DMPlexGetChart(*dmOverlap, &pStart, &pEnd);CHKERRQ(ierr);
1537a157612eSMichael Lange   ierr = PetscSFCreate(comm, &newPointSF);;CHKERRQ(ierr);
1538a157612eSMichael Lange   ierr = PetscSFSetGraph(newPointSF, pEnd - pStart, numGhostPoints, ghostLocal, PETSC_OWN_POINTER, ghostRemote, PETSC_OWN_POINTER);CHKERRQ(ierr);
1539a157612eSMichael Lange   ierr = DMSetPointSF(*dmOverlap, newPointSF);CHKERRQ(ierr);
1540a157612eSMichael Lange   /* Cleanup overlap partition */
1541a157612eSMichael Lange   ierr = ISLocalToGlobalMappingDestroy(&overlapRenumbering);CHKERRQ(ierr);
1542a157612eSMichael Lange   ierr = PetscSFDestroy(&overlapSF);CHKERRQ(ierr);
1543a157612eSMichael Lange   ierr = PetscFree2(pointIDs, recvPointIDs);CHKERRQ(ierr);
1544a157612eSMichael Lange   if (sf) *sf = migrationSF;
1545a157612eSMichael Lange   else    {ierr = PetscSFDestroy(&migrationSF);CHKERRQ(ierr);}
1546a157612eSMichael Lange   ierr = DMSetFromOptions(*dmOverlap);CHKERRQ(ierr);
1547a157612eSMichael Lange   PetscFunctionReturn(0);
1548a157612eSMichael Lange }
1549