xref: /petsc/src/dm/impls/plex/plexdistribute.c (revision b9f40539eb16c55bb65092e9c620e0c7f5cb6588)
170034214SMatthew G. Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
270034214SMatthew G. Knepley 
370034214SMatthew G. Knepley #undef __FUNCT__
470034214SMatthew G. Knepley #define __FUNCT__ "DMPlexSetAdjacencyUseCone"
570034214SMatthew G. Knepley /*@
670034214SMatthew G. Knepley   DMPlexSetAdjacencyUseCone - Define adjacency in the mesh using either the cone or the support first
770034214SMatthew G. Knepley 
870034214SMatthew G. Knepley   Input Parameters:
970034214SMatthew G. Knepley + dm      - The DM object
1070034214SMatthew G. Knepley - useCone - Flag to use the cone first
1170034214SMatthew G. Knepley 
1270034214SMatthew G. Knepley   Level: intermediate
1370034214SMatthew G. Knepley 
1470034214SMatthew G. Knepley   Notes:
1570034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
164b6b44bdSMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
1770034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
1870034214SMatthew G. Knepley 
1970034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
2070034214SMatthew G. Knepley @*/
2170034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseCone(DM dm, PetscBool useCone)
2270034214SMatthew G. Knepley {
2370034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
2470034214SMatthew G. Knepley 
2570034214SMatthew G. Knepley   PetscFunctionBegin;
2670034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2770034214SMatthew G. Knepley   mesh->useCone = useCone;
2870034214SMatthew G. Knepley   PetscFunctionReturn(0);
2970034214SMatthew G. Knepley }
3070034214SMatthew G. Knepley 
3170034214SMatthew G. Knepley #undef __FUNCT__
3270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacencyUseCone"
3370034214SMatthew G. Knepley /*@
3470034214SMatthew G. Knepley   DMPlexGetAdjacencyUseCone - Query whether adjacency in the mesh uses the cone or the support first
3570034214SMatthew G. Knepley 
3670034214SMatthew G. Knepley   Input Parameter:
3770034214SMatthew G. Knepley . dm      - The DM object
3870034214SMatthew G. Knepley 
3970034214SMatthew G. Knepley   Output Parameter:
4070034214SMatthew G. Knepley . useCone - Flag to use the cone first
4170034214SMatthew G. Knepley 
4270034214SMatthew G. Knepley   Level: intermediate
4370034214SMatthew G. Knepley 
4470034214SMatthew G. Knepley   Notes:
4570034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
464b6b44bdSMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
4770034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
4870034214SMatthew G. Knepley 
4970034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
5070034214SMatthew G. Knepley @*/
5170034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseCone(DM dm, PetscBool *useCone)
5270034214SMatthew G. Knepley {
5370034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
5470034214SMatthew G. Knepley 
5570034214SMatthew G. Knepley   PetscFunctionBegin;
5670034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5770034214SMatthew G. Knepley   PetscValidIntPointer(useCone, 2);
5870034214SMatthew G. Knepley   *useCone = mesh->useCone;
5970034214SMatthew G. Knepley   PetscFunctionReturn(0);
6070034214SMatthew G. Knepley }
6170034214SMatthew G. Knepley 
6270034214SMatthew G. Knepley #undef __FUNCT__
6370034214SMatthew G. Knepley #define __FUNCT__ "DMPlexSetAdjacencyUseClosure"
6470034214SMatthew G. Knepley /*@
6570034214SMatthew G. Knepley   DMPlexSetAdjacencyUseClosure - Define adjacency in the mesh using the transitive closure
6670034214SMatthew G. Knepley 
6770034214SMatthew G. Knepley   Input Parameters:
6870034214SMatthew G. Knepley + dm      - The DM object
6970034214SMatthew G. Knepley - useClosure - Flag to use the closure
7070034214SMatthew G. Knepley 
7170034214SMatthew G. Knepley   Level: intermediate
7270034214SMatthew G. Knepley 
7370034214SMatthew G. Knepley   Notes:
7470034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
754b6b44bdSMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
7670034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
7770034214SMatthew G. Knepley 
7870034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
7970034214SMatthew G. Knepley @*/
8070034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseClosure(DM dm, PetscBool useClosure)
8170034214SMatthew G. Knepley {
8270034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
8370034214SMatthew G. Knepley 
8470034214SMatthew G. Knepley   PetscFunctionBegin;
8570034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8670034214SMatthew G. Knepley   mesh->useClosure = useClosure;
8770034214SMatthew G. Knepley   PetscFunctionReturn(0);
8870034214SMatthew G. Knepley }
8970034214SMatthew G. Knepley 
9070034214SMatthew G. Knepley #undef __FUNCT__
9170034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacencyUseClosure"
9270034214SMatthew G. Knepley /*@
9370034214SMatthew G. Knepley   DMPlexGetAdjacencyUseClosure - Query whether adjacency in the mesh uses the transitive closure
9470034214SMatthew G. Knepley 
9570034214SMatthew G. Knepley   Input Parameter:
9670034214SMatthew G. Knepley . dm      - The DM object
9770034214SMatthew G. Knepley 
9870034214SMatthew G. Knepley   Output Parameter:
9970034214SMatthew G. Knepley . useClosure - Flag to use the closure
10070034214SMatthew G. Knepley 
10170034214SMatthew G. Knepley   Level: intermediate
10270034214SMatthew G. Knepley 
10370034214SMatthew G. Knepley   Notes:
10470034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)),   useCone = PETSC_FALSE, useClosure = PETSC_TRUE
1054b6b44bdSMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
10670034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)),   useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
10770034214SMatthew G. Knepley 
10870034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
10970034214SMatthew G. Knepley @*/
11070034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseClosure(DM dm, PetscBool *useClosure)
11170034214SMatthew G. Knepley {
11270034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
11370034214SMatthew G. Knepley 
11470034214SMatthew G. Knepley   PetscFunctionBegin;
11570034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11670034214SMatthew G. Knepley   PetscValidIntPointer(useClosure, 2);
11770034214SMatthew G. Knepley   *useClosure = mesh->useClosure;
11870034214SMatthew G. Knepley   PetscFunctionReturn(0);
11970034214SMatthew G. Knepley }
12070034214SMatthew G. Knepley 
12170034214SMatthew G. Knepley #undef __FUNCT__
122a17985deSToby Isaac #define __FUNCT__ "DMPlexSetAdjacencyUseAnchors"
1238b0b4c70SToby Isaac /*@
124a17985deSToby Isaac   DMPlexSetAdjacencyUseAnchors - Define adjacency in the mesh using the point-to-point constraints.
1258b0b4c70SToby Isaac 
1268b0b4c70SToby Isaac   Input Parameters:
1278b0b4c70SToby Isaac + dm      - The DM object
1285b317d89SToby Isaac - useAnchors - Flag to use the constraints.  If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place.
1298b0b4c70SToby Isaac 
1308b0b4c70SToby Isaac   Level: intermediate
1318b0b4c70SToby Isaac 
132a17985deSToby Isaac .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
1338b0b4c70SToby Isaac @*/
1345b317d89SToby Isaac PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM dm, PetscBool useAnchors)
1358b0b4c70SToby Isaac {
1368b0b4c70SToby Isaac   DM_Plex *mesh = (DM_Plex *) dm->data;
1378b0b4c70SToby Isaac 
1388b0b4c70SToby Isaac   PetscFunctionBegin;
1398b0b4c70SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1405b317d89SToby Isaac   mesh->useAnchors = useAnchors;
1418b0b4c70SToby Isaac   PetscFunctionReturn(0);
1428b0b4c70SToby Isaac }
1438b0b4c70SToby Isaac 
1448b0b4c70SToby Isaac #undef __FUNCT__
145a17985deSToby Isaac #define __FUNCT__ "DMPlexGetAdjacencyUseAnchors"
1468b0b4c70SToby Isaac /*@
147a17985deSToby Isaac   DMPlexGetAdjacencyUseAnchors - Query whether adjacency in the mesh uses the point-to-point constraints.
1488b0b4c70SToby Isaac 
1498b0b4c70SToby Isaac   Input Parameter:
1508b0b4c70SToby Isaac . dm      - The DM object
1518b0b4c70SToby Isaac 
1528b0b4c70SToby Isaac   Output Parameter:
1535b317d89SToby Isaac . useAnchors - Flag to use the closure.  If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place.
1548b0b4c70SToby Isaac 
1558b0b4c70SToby Isaac   Level: intermediate
1568b0b4c70SToby Isaac 
157a17985deSToby Isaac .seealso: DMPlexSetAdjacencyUseAnchors(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors()
1588b0b4c70SToby Isaac @*/
1595b317d89SToby Isaac PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM dm, PetscBool *useAnchors)
1608b0b4c70SToby Isaac {
1618b0b4c70SToby Isaac   DM_Plex *mesh = (DM_Plex *) dm->data;
1628b0b4c70SToby Isaac 
1638b0b4c70SToby Isaac   PetscFunctionBegin;
1648b0b4c70SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1655b317d89SToby Isaac   PetscValidIntPointer(useAnchors, 2);
1665b317d89SToby Isaac   *useAnchors = mesh->useAnchors;
1678b0b4c70SToby Isaac   PetscFunctionReturn(0);
1688b0b4c70SToby Isaac }
1698b0b4c70SToby Isaac 
1708b0b4c70SToby Isaac #undef __FUNCT__
17170034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Cone_Internal"
17270034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Cone_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
17370034214SMatthew G. Knepley {
17470034214SMatthew G. Knepley   const PetscInt *cone = NULL;
17570034214SMatthew G. Knepley   PetscInt        numAdj = 0, maxAdjSize = *adjSize, coneSize, c;
17670034214SMatthew G. Knepley   PetscErrorCode  ierr;
17770034214SMatthew G. Knepley 
17870034214SMatthew G. Knepley   PetscFunctionBeginHot;
17970034214SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
18070034214SMatthew G. Knepley   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
1814b6b44bdSMatthew G. Knepley   for (c = 0; c <= coneSize; ++c) {
1824b6b44bdSMatthew G. Knepley     const PetscInt  point   = !c ? p : cone[c-1];
18370034214SMatthew G. Knepley     const PetscInt *support = NULL;
18470034214SMatthew G. Knepley     PetscInt        supportSize, s, q;
18570034214SMatthew G. Knepley 
1864b6b44bdSMatthew G. Knepley     ierr = DMPlexGetSupportSize(dm, point, &supportSize);CHKERRQ(ierr);
1874b6b44bdSMatthew G. Knepley     ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr);
18870034214SMatthew G. Knepley     for (s = 0; s < supportSize; ++s) {
18970034214SMatthew G. Knepley       for (q = 0; q < numAdj || (adj[numAdj++] = support[s],0); ++q) {
19070034214SMatthew G. Knepley         if (support[s] == adj[q]) break;
19170034214SMatthew G. Knepley       }
19270034214SMatthew G. Knepley       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
19370034214SMatthew G. Knepley     }
19470034214SMatthew G. Knepley   }
19570034214SMatthew G. Knepley   *adjSize = numAdj;
19670034214SMatthew G. Knepley   PetscFunctionReturn(0);
19770034214SMatthew G. Knepley }
19870034214SMatthew G. Knepley 
19970034214SMatthew G. Knepley #undef __FUNCT__
20070034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Support_Internal"
20170034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Support_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
20270034214SMatthew G. Knepley {
20370034214SMatthew G. Knepley   const PetscInt *support = NULL;
20470034214SMatthew G. Knepley   PetscInt        numAdj   = 0, maxAdjSize = *adjSize, supportSize, s;
20570034214SMatthew G. Knepley   PetscErrorCode  ierr;
20670034214SMatthew G. Knepley 
20770034214SMatthew G. Knepley   PetscFunctionBeginHot;
20870034214SMatthew G. Knepley   ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr);
20970034214SMatthew G. Knepley   ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr);
2104b6b44bdSMatthew G. Knepley   for (s = 0; s <= supportSize; ++s) {
2114b6b44bdSMatthew G. Knepley     const PetscInt  point = !s ? p : support[s-1];
21270034214SMatthew G. Knepley     const PetscInt *cone  = NULL;
21370034214SMatthew G. Knepley     PetscInt        coneSize, c, q;
21470034214SMatthew G. Knepley 
2154b6b44bdSMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr);
2164b6b44bdSMatthew G. Knepley     ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr);
21770034214SMatthew G. Knepley     for (c = 0; c < coneSize; ++c) {
21870034214SMatthew G. Knepley       for (q = 0; q < numAdj || (adj[numAdj++] = cone[c],0); ++q) {
21970034214SMatthew G. Knepley         if (cone[c] == adj[q]) break;
22070034214SMatthew G. Knepley       }
22170034214SMatthew G. Knepley       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
22270034214SMatthew G. Knepley     }
22370034214SMatthew G. Knepley   }
22470034214SMatthew G. Knepley   *adjSize = numAdj;
22570034214SMatthew G. Knepley   PetscFunctionReturn(0);
22670034214SMatthew G. Knepley }
22770034214SMatthew G. Knepley 
22870034214SMatthew G. Knepley #undef __FUNCT__
22970034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Transitive_Internal"
23070034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Transitive_Internal(DM dm, PetscInt p, PetscBool useClosure, PetscInt *adjSize, PetscInt adj[])
23170034214SMatthew G. Knepley {
23270034214SMatthew G. Knepley   PetscInt      *star = NULL;
23370034214SMatthew G. Knepley   PetscInt       numAdj = 0, maxAdjSize = *adjSize, starSize, s;
23470034214SMatthew G. Knepley   PetscErrorCode ierr;
23570034214SMatthew G. Knepley 
23670034214SMatthew G. Knepley   PetscFunctionBeginHot;
23770034214SMatthew G. Knepley   ierr = DMPlexGetTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr);
23870034214SMatthew G. Knepley   for (s = 0; s < starSize*2; s += 2) {
23970034214SMatthew G. Knepley     const PetscInt *closure = NULL;
24070034214SMatthew G. Knepley     PetscInt        closureSize, c, q;
24170034214SMatthew G. Knepley 
24270034214SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
24370034214SMatthew G. Knepley     for (c = 0; c < closureSize*2; c += 2) {
24470034214SMatthew G. Knepley       for (q = 0; q < numAdj || (adj[numAdj++] = closure[c],0); ++q) {
24570034214SMatthew G. Knepley         if (closure[c] == adj[q]) break;
24670034214SMatthew G. Knepley       }
24770034214SMatthew G. Knepley       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
24870034214SMatthew G. Knepley     }
24970034214SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
25070034214SMatthew G. Knepley   }
25170034214SMatthew G. Knepley   ierr = DMPlexRestoreTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr);
25270034214SMatthew G. Knepley   *adjSize = numAdj;
25370034214SMatthew G. Knepley   PetscFunctionReturn(0);
25470034214SMatthew G. Knepley }
25570034214SMatthew G. Knepley 
25670034214SMatthew G. Knepley #undef __FUNCT__
25770034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Internal"
2585b317d89SToby Isaac PetscErrorCode DMPlexGetAdjacency_Internal(DM dm, PetscInt p, PetscBool useCone, PetscBool useTransitiveClosure, PetscBool useAnchors, PetscInt *adjSize, PetscInt *adj[])
25970034214SMatthew G. Knepley {
26079bad331SMatthew G. Knepley   static PetscInt asiz = 0;
2618b0b4c70SToby Isaac   PetscInt maxAnchors = 1;
2628b0b4c70SToby Isaac   PetscInt aStart = -1, aEnd = -1;
2638b0b4c70SToby Isaac   PetscInt maxAdjSize;
2648b0b4c70SToby Isaac   PetscSection aSec = NULL;
2658b0b4c70SToby Isaac   IS aIS = NULL;
2668b0b4c70SToby Isaac   const PetscInt *anchors;
26770034214SMatthew G. Knepley   PetscErrorCode  ierr;
26870034214SMatthew G. Knepley 
26970034214SMatthew G. Knepley   PetscFunctionBeginHot;
2705b317d89SToby Isaac   if (useAnchors) {
271a17985deSToby Isaac     ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr);
2728b0b4c70SToby Isaac     if (aSec) {
2738b0b4c70SToby Isaac       ierr = PetscSectionGetMaxDof(aSec,&maxAnchors);CHKERRQ(ierr);
27424c766afSToby Isaac       maxAnchors = PetscMax(1,maxAnchors);
2758b0b4c70SToby Isaac       ierr = PetscSectionGetChart(aSec,&aStart,&aEnd);CHKERRQ(ierr);
2768b0b4c70SToby Isaac       ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr);
2778b0b4c70SToby Isaac     }
2788b0b4c70SToby Isaac   }
27979bad331SMatthew G. Knepley   if (!*adj) {
28024c766afSToby Isaac     PetscInt depth, coneSeries, supportSeries, maxC, maxS, pStart, pEnd;
28179bad331SMatthew G. Knepley 
28224c766afSToby Isaac     ierr  = DMPlexGetChart(dm, &pStart,&pEnd);CHKERRQ(ierr);
28379bad331SMatthew G. Knepley     ierr  = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
28424c766afSToby Isaac     ierr  = DMPlexGetMaxSizes(dm, &maxC, &maxS);CHKERRQ(ierr);
28524c766afSToby Isaac     coneSeries    = (maxC > 1) ? ((PetscPowInt(maxC,depth+1)-1)/(maxC-1)) : depth+1;
28624c766afSToby Isaac     supportSeries = (maxS > 1) ? ((PetscPowInt(maxS,depth+1)-1)/(maxS-1)) : depth+1;
28724c766afSToby Isaac     asiz  = PetscMax(PetscPowInt(maxS,depth)*coneSeries,PetscPowInt(maxC,depth)*supportSeries);
2888b0b4c70SToby Isaac     asiz *= maxAnchors;
28924c766afSToby Isaac     asiz  = PetscMin(asiz,pEnd-pStart);
29079bad331SMatthew G. Knepley     ierr  = PetscMalloc1(asiz,adj);CHKERRQ(ierr);
29179bad331SMatthew G. Knepley   }
29279bad331SMatthew G. Knepley   if (*adjSize < 0) *adjSize = asiz;
2938b0b4c70SToby Isaac   maxAdjSize = *adjSize;
29470034214SMatthew G. Knepley   if (useTransitiveClosure) {
29579bad331SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Transitive_Internal(dm, p, useCone, adjSize, *adj);CHKERRQ(ierr);
29670034214SMatthew G. Knepley   } else if (useCone) {
29779bad331SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Cone_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
29870034214SMatthew G. Knepley   } else {
29979bad331SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Support_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
30070034214SMatthew G. Knepley   }
3015b317d89SToby Isaac   if (useAnchors && aSec) {
3028b0b4c70SToby Isaac     PetscInt origSize = *adjSize;
3038b0b4c70SToby Isaac     PetscInt numAdj = origSize;
3048b0b4c70SToby Isaac     PetscInt i = 0, j;
3058b0b4c70SToby Isaac     PetscInt *orig = *adj;
3068b0b4c70SToby Isaac 
3078b0b4c70SToby Isaac     while (i < origSize) {
3088b0b4c70SToby Isaac       PetscInt p = orig[i];
3098b0b4c70SToby Isaac       PetscInt aDof = 0;
3108b0b4c70SToby Isaac 
3118b0b4c70SToby Isaac       if (p >= aStart && p < aEnd) {
3128b0b4c70SToby Isaac         ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr);
3138b0b4c70SToby Isaac       }
3148b0b4c70SToby Isaac       if (aDof) {
3158b0b4c70SToby Isaac         PetscInt aOff;
3168b0b4c70SToby Isaac         PetscInt s, q;
3178b0b4c70SToby Isaac 
3188b0b4c70SToby Isaac         for (j = i + 1; j < numAdj; j++) {
3198b0b4c70SToby Isaac           orig[j - 1] = orig[j];
3208b0b4c70SToby Isaac         }
3218b0b4c70SToby Isaac         origSize--;
3228b0b4c70SToby Isaac         numAdj--;
3238b0b4c70SToby Isaac         ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr);
3248b0b4c70SToby Isaac         for (s = 0; s < aDof; ++s) {
3258b0b4c70SToby Isaac           for (q = 0; q < numAdj || (orig[numAdj++] = anchors[aOff+s],0); ++q) {
3268b0b4c70SToby Isaac             if (anchors[aOff+s] == orig[q]) break;
3278b0b4c70SToby Isaac           }
3288b0b4c70SToby Isaac           if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
3298b0b4c70SToby Isaac         }
3308b0b4c70SToby Isaac       }
3318b0b4c70SToby Isaac       else {
3328b0b4c70SToby Isaac         i++;
3338b0b4c70SToby Isaac       }
3348b0b4c70SToby Isaac     }
3358b0b4c70SToby Isaac     *adjSize = numAdj;
3368b0b4c70SToby Isaac     ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr);
3378b0b4c70SToby Isaac   }
33870034214SMatthew G. Knepley   PetscFunctionReturn(0);
33970034214SMatthew G. Knepley }
34070034214SMatthew G. Knepley 
34170034214SMatthew G. Knepley #undef __FUNCT__
34270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency"
34370034214SMatthew G. Knepley /*@
34470034214SMatthew G. Knepley   DMPlexGetAdjacency - Return all points adjacent to the given point
34570034214SMatthew G. Knepley 
34670034214SMatthew G. Knepley   Input Parameters:
34770034214SMatthew G. Knepley + dm - The DM object
34870034214SMatthew G. Knepley . p  - The point
34970034214SMatthew G. Knepley . adjSize - The maximum size of adj if it is non-NULL, or PETSC_DETERMINE
35070034214SMatthew G. Knepley - adj - Either NULL so that the array is allocated, or an existing array with size adjSize
35170034214SMatthew G. Knepley 
35270034214SMatthew G. Knepley   Output Parameters:
35370034214SMatthew G. Knepley + adjSize - The number of adjacent points
35470034214SMatthew G. Knepley - adj - The adjacent points
35570034214SMatthew G. Knepley 
35670034214SMatthew G. Knepley   Level: advanced
35770034214SMatthew G. Knepley 
35870034214SMatthew G. Knepley   Notes: The user must PetscFree the adj array if it was not passed in.
35970034214SMatthew G. Knepley 
36070034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexDistribute(), DMCreateMatrix(), DMPlexPreallocateOperator()
36170034214SMatthew G. Knepley @*/
36270034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacency(DM dm, PetscInt p, PetscInt *adjSize, PetscInt *adj[])
36370034214SMatthew G. Knepley {
36470034214SMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
36570034214SMatthew G. Knepley   PetscErrorCode ierr;
36670034214SMatthew G. Knepley 
36770034214SMatthew G. Knepley   PetscFunctionBeginHot;
36870034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
36970034214SMatthew G. Knepley   PetscValidPointer(adjSize,3);
37070034214SMatthew G. Knepley   PetscValidPointer(adj,4);
3715b317d89SToby Isaac   ierr = DMPlexGetAdjacency_Internal(dm, p, mesh->useCone, mesh->useClosure, mesh->useAnchors, adjSize, adj);CHKERRQ(ierr);
37270034214SMatthew G. Knepley   PetscFunctionReturn(0);
37370034214SMatthew G. Knepley }
374b0a623aaSMatthew G. Knepley #undef __FUNCT__
375b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexCreateTwoSidedProcessSF"
376b0a623aaSMatthew G. Knepley /*@
377b0a623aaSMatthew G. Knepley   DMPlexCreateTwoSidedProcessSF - Create an SF which just has process connectivity
378b0a623aaSMatthew G. Knepley 
379b0a623aaSMatthew G. Knepley   Collective on DM
380b0a623aaSMatthew G. Knepley 
381b0a623aaSMatthew G. Knepley   Input Parameters:
382b0a623aaSMatthew G. Knepley + dm      - The DM
383b0a623aaSMatthew G. Knepley - sfPoint - The PetscSF which encodes point connectivity
384b0a623aaSMatthew G. Knepley 
385b0a623aaSMatthew G. Knepley   Output Parameters:
386b0a623aaSMatthew G. Knepley + processRanks - A list of process neighbors, or NULL
387b0a623aaSMatthew G. Knepley - sfProcess    - An SF encoding the two-sided process connectivity, or NULL
388b0a623aaSMatthew G. Knepley 
389b0a623aaSMatthew G. Knepley   Level: developer
390b0a623aaSMatthew G. Knepley 
391b0a623aaSMatthew G. Knepley .seealso: PetscSFCreate(), DMPlexCreateProcessSF()
392b0a623aaSMatthew G. Knepley @*/
393b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM dm, PetscSF sfPoint, PetscSection rootRankSection, IS rootRanks, PetscSection leafRankSection, IS leafRanks, IS *processRanks, PetscSF *sfProcess)
394b0a623aaSMatthew G. Knepley {
395b0a623aaSMatthew G. Knepley   const PetscSFNode *remotePoints;
396b0a623aaSMatthew G. Knepley   PetscInt          *localPointsNew;
397b0a623aaSMatthew G. Knepley   PetscSFNode       *remotePointsNew;
398b0a623aaSMatthew G. Knepley   const PetscInt    *nranks;
399b0a623aaSMatthew G. Knepley   PetscInt          *ranksNew;
400b0a623aaSMatthew G. Knepley   PetscBT            neighbors;
401b0a623aaSMatthew G. Knepley   PetscInt           pStart, pEnd, p, numLeaves, l, numNeighbors, n;
402b0a623aaSMatthew G. Knepley   PetscMPIInt        numProcs, proc, rank;
403b0a623aaSMatthew G. Knepley   PetscErrorCode     ierr;
404b0a623aaSMatthew G. Knepley 
405b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
406b0a623aaSMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
407b0a623aaSMatthew G. Knepley   PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2);
408b0a623aaSMatthew G. Knepley   if (processRanks) {PetscValidPointer(processRanks, 3);}
409b0a623aaSMatthew G. Knepley   if (sfProcess)    {PetscValidPointer(sfProcess, 4);}
410b0a623aaSMatthew G. Knepley   ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr);
411b0a623aaSMatthew G. Knepley   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr);
412b0a623aaSMatthew G. Knepley   ierr = PetscSFGetGraph(sfPoint, NULL, &numLeaves, NULL, &remotePoints);CHKERRQ(ierr);
413b0a623aaSMatthew G. Knepley   ierr = PetscBTCreate(numProcs, &neighbors);CHKERRQ(ierr);
414b0a623aaSMatthew G. Knepley   ierr = PetscBTMemzero(numProcs, neighbors);CHKERRQ(ierr);
415b0a623aaSMatthew G. Knepley   /* Compute root-to-leaf process connectivity */
416b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetChart(rootRankSection, &pStart, &pEnd);CHKERRQ(ierr);
417b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(rootRanks, &nranks);CHKERRQ(ierr);
418b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
419b0a623aaSMatthew G. Knepley     PetscInt ndof, noff, n;
420b0a623aaSMatthew G. Knepley 
421b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetDof(rootRankSection, p, &ndof);CHKERRQ(ierr);
422b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(rootRankSection, p, &noff);CHKERRQ(ierr);
423b0a623aaSMatthew G. Knepley     for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);}
424b0a623aaSMatthew G. Knepley   }
425b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(rootRanks, &nranks);CHKERRQ(ierr);
426b0a623aaSMatthew G. Knepley   /* Compute leaf-to-neighbor process connectivity */
427b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetChart(leafRankSection, &pStart, &pEnd);CHKERRQ(ierr);
428b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(leafRanks, &nranks);CHKERRQ(ierr);
429b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
430b0a623aaSMatthew G. Knepley     PetscInt ndof, noff, n;
431b0a623aaSMatthew G. Knepley 
432b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetDof(leafRankSection, p, &ndof);CHKERRQ(ierr);
433b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(leafRankSection, p, &noff);CHKERRQ(ierr);
434b0a623aaSMatthew G. Knepley     for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);}
435b0a623aaSMatthew G. Knepley   }
436b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(leafRanks, &nranks);CHKERRQ(ierr);
437b0a623aaSMatthew G. Knepley   /* Compute leaf-to-root process connectivity */
438b0a623aaSMatthew G. Knepley   for (l = 0; l < numLeaves; ++l) {PetscBTSet(neighbors, remotePoints[l].rank);}
439b0a623aaSMatthew G. Knepley   /* Calculate edges */
440b0a623aaSMatthew G. Knepley   PetscBTClear(neighbors, rank);
441b0a623aaSMatthew G. Knepley   for(proc = 0, numNeighbors = 0; proc < numProcs; ++proc) {if (PetscBTLookup(neighbors, proc)) ++numNeighbors;}
442b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(numNeighbors, &ranksNew);CHKERRQ(ierr);
443b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(numNeighbors, &localPointsNew);CHKERRQ(ierr);
444b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(numNeighbors, &remotePointsNew);CHKERRQ(ierr);
445b0a623aaSMatthew G. Knepley   for(proc = 0, n = 0; proc < numProcs; ++proc) {
446b0a623aaSMatthew G. Knepley     if (PetscBTLookup(neighbors, proc)) {
447b0a623aaSMatthew G. Knepley       ranksNew[n]              = proc;
448b0a623aaSMatthew G. Knepley       localPointsNew[n]        = proc;
449b0a623aaSMatthew G. Knepley       remotePointsNew[n].index = rank;
450b0a623aaSMatthew G. Knepley       remotePointsNew[n].rank  = proc;
451b0a623aaSMatthew G. Knepley       ++n;
452b0a623aaSMatthew G. Knepley     }
453b0a623aaSMatthew G. Knepley   }
454b0a623aaSMatthew G. Knepley   ierr = PetscBTDestroy(&neighbors);CHKERRQ(ierr);
455b0a623aaSMatthew G. Knepley   if (processRanks) {ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), numNeighbors, ranksNew, PETSC_OWN_POINTER, processRanks);CHKERRQ(ierr);}
456b0a623aaSMatthew G. Knepley   else              {ierr = PetscFree(ranksNew);CHKERRQ(ierr);}
457b0a623aaSMatthew G. Knepley   if (sfProcess) {
458b0a623aaSMatthew G. Knepley     ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess);CHKERRQ(ierr);
459b0a623aaSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) *sfProcess, "Two-Sided Process SF");CHKERRQ(ierr);
460b0a623aaSMatthew G. Knepley     ierr = PetscSFSetFromOptions(*sfProcess);CHKERRQ(ierr);
461b0a623aaSMatthew G. Knepley     ierr = PetscSFSetGraph(*sfProcess, numProcs, numNeighbors, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr);
462b0a623aaSMatthew G. Knepley   }
463b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
464b0a623aaSMatthew G. Knepley }
465b0a623aaSMatthew G. Knepley 
466b0a623aaSMatthew G. Knepley #undef __FUNCT__
467b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeOwnership"
468b0a623aaSMatthew G. Knepley /*@
469b0a623aaSMatthew G. Knepley   DMPlexDistributeOwnership - Compute owner information for shared points. This basically gets two-sided for an SF.
470b0a623aaSMatthew G. Knepley 
471b0a623aaSMatthew G. Knepley   Collective on DM
472b0a623aaSMatthew G. Knepley 
473b0a623aaSMatthew G. Knepley   Input Parameter:
474b0a623aaSMatthew G. Knepley . dm - The DM
475b0a623aaSMatthew G. Knepley 
476b0a623aaSMatthew G. Knepley   Output Parameters:
477b0a623aaSMatthew G. Knepley + rootSection - The number of leaves for a given root point
478b0a623aaSMatthew G. Knepley . rootrank    - The rank of each edge into the root point
479b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point
480b0a623aaSMatthew G. Knepley - leafrank    - The rank of each process sharing a leaf point
481b0a623aaSMatthew G. Knepley 
482b0a623aaSMatthew G. Knepley   Level: developer
483b0a623aaSMatthew G. Knepley 
484b0a623aaSMatthew G. Knepley .seealso: DMPlexCreateOverlap()
485b0a623aaSMatthew G. Knepley @*/
486b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexDistributeOwnership(DM dm, PetscSection rootSection, IS *rootrank, PetscSection leafSection, IS *leafrank)
487b0a623aaSMatthew G. Knepley {
488b0a623aaSMatthew G. Knepley   MPI_Comm        comm;
489b0a623aaSMatthew G. Knepley   PetscSF         sfPoint;
490b0a623aaSMatthew G. Knepley   const PetscInt *rootdegree;
491b0a623aaSMatthew G. Knepley   PetscInt       *myrank, *remoterank;
492b0a623aaSMatthew G. Knepley   PetscInt        pStart, pEnd, p, nedges;
493b0a623aaSMatthew G. Knepley   PetscMPIInt     rank;
494b0a623aaSMatthew G. Knepley   PetscErrorCode  ierr;
495b0a623aaSMatthew G. Knepley 
496b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
497b0a623aaSMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
498b0a623aaSMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
499b0a623aaSMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
500b0a623aaSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
501b0a623aaSMatthew G. Knepley   /* Compute number of leaves for each root */
502b0a623aaSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) rootSection, "Root Section");CHKERRQ(ierr);
503b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetChart(rootSection, pStart, pEnd);CHKERRQ(ierr);
504b0a623aaSMatthew G. Knepley   ierr = PetscSFComputeDegreeBegin(sfPoint, &rootdegree);CHKERRQ(ierr);
505b0a623aaSMatthew G. Knepley   ierr = PetscSFComputeDegreeEnd(sfPoint, &rootdegree);CHKERRQ(ierr);
506b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {ierr = PetscSectionSetDof(rootSection, p, rootdegree[p-pStart]);CHKERRQ(ierr);}
507b0a623aaSMatthew G. Knepley   ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr);
508b0a623aaSMatthew G. Knepley   /* Gather rank of each leaf to root */
509b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(rootSection, &nedges);CHKERRQ(ierr);
510b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(pEnd-pStart, &myrank);CHKERRQ(ierr);
511b0a623aaSMatthew G. Knepley   ierr = PetscMalloc1(nedges,  &remoterank);CHKERRQ(ierr);
512b0a623aaSMatthew G. Knepley   for (p = 0; p < pEnd-pStart; ++p) myrank[p] = rank;
513b0a623aaSMatthew G. Knepley   ierr = PetscSFGatherBegin(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr);
514b0a623aaSMatthew G. Knepley   ierr = PetscSFGatherEnd(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr);
515b0a623aaSMatthew G. Knepley   ierr = PetscFree(myrank);CHKERRQ(ierr);
516b0a623aaSMatthew G. Knepley   ierr = ISCreateGeneral(comm, nedges, remoterank, PETSC_OWN_POINTER, rootrank);CHKERRQ(ierr);
517b0a623aaSMatthew G. Knepley   /* Distribute remote ranks to leaves */
518b0a623aaSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) leafSection, "Leaf Section");CHKERRQ(ierr);
519b0a623aaSMatthew G. Knepley   ierr = DMPlexDistributeFieldIS(dm, sfPoint, rootSection, *rootrank, leafSection, leafrank);CHKERRQ(ierr);
520b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
521b0a623aaSMatthew G. Knepley }
522b0a623aaSMatthew G. Knepley 
523b0a623aaSMatthew G. Knepley #undef __FUNCT__
524b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexCreateOverlap"
525278397a0SMatthew G. Knepley /*@C
526b0a623aaSMatthew G. Knepley   DMPlexCreateOverlap - Compute owner information for shared points. This basically gets two-sided for an SF.
527b0a623aaSMatthew G. Knepley 
528b0a623aaSMatthew G. Knepley   Collective on DM
529b0a623aaSMatthew G. Knepley 
530b0a623aaSMatthew G. Knepley   Input Parameters:
531b0a623aaSMatthew G. Knepley + dm          - The DM
53224d039d7SMichael Lange . levels      - Number of overlap levels
533b0a623aaSMatthew G. Knepley . rootSection - The number of leaves for a given root point
534b0a623aaSMatthew G. Knepley . rootrank    - The rank of each edge into the root point
535b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point
536b0a623aaSMatthew G. Knepley - leafrank    - The rank of each process sharing a leaf point
537b0a623aaSMatthew G. Knepley 
538b0a623aaSMatthew G. Knepley   Output Parameters:
5391fd9873aSMichael Lange + overlapSF   - The star forest mapping remote overlap points into a contiguous leaf space
540b0a623aaSMatthew G. Knepley 
541b0a623aaSMatthew G. Knepley   Level: developer
542b0a623aaSMatthew G. Knepley 
5431fd9873aSMichael Lange .seealso: DMPlexDistributeOwnership(), DMPlexDistribute()
544b0a623aaSMatthew G. Knepley @*/
54524d039d7SMichael Lange PetscErrorCode DMPlexCreateOverlap(DM dm, PetscInt levels, PetscSection rootSection, IS rootrank, PetscSection leafSection, IS leafrank, PetscSF *overlapSF)
546b0a623aaSMatthew G. Knepley {
547e540f424SMichael Lange   MPI_Comm           comm;
548b0a623aaSMatthew G. Knepley   DMLabel            ovAdjByRank; /* A DMLabel containing all points adjacent to shared points, separated by rank (value in label) */
549b0a623aaSMatthew G. Knepley   PetscSF            sfPoint, sfProc;
550e540f424SMichael Lange   DMLabel            ovLeafLabel;
551b0a623aaSMatthew G. Knepley   const PetscSFNode *remote;
552b0a623aaSMatthew G. Knepley   const PetscInt    *local;
5531fd9873aSMichael Lange   const PetscInt    *nrank, *rrank;
554b0a623aaSMatthew G. Knepley   PetscInt          *adj = NULL;
5551fd9873aSMichael Lange   PetscInt           pStart, pEnd, p, sStart, sEnd, nleaves, l;
556b0a623aaSMatthew G. Knepley   PetscMPIInt        rank, numProcs;
55726a7d390SMatthew G. Knepley   PetscBool          useCone, useClosure, flg;
558b0a623aaSMatthew G. Knepley   PetscErrorCode     ierr;
559b0a623aaSMatthew G. Knepley 
560b0a623aaSMatthew G. Knepley   PetscFunctionBegin;
561e540f424SMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
562e540f424SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
563e540f424SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
564b0a623aaSMatthew G. Knepley   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
565b0a623aaSMatthew G. Knepley   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
566b0a623aaSMatthew G. Knepley   ierr = PetscSectionGetChart(leafSection, &sStart, &sEnd);CHKERRQ(ierr);
567b0a623aaSMatthew G. Knepley   ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr);
568b0a623aaSMatthew G. Knepley   ierr = DMLabelCreate("Overlap adjacency", &ovAdjByRank);CHKERRQ(ierr);
569b0a623aaSMatthew G. Knepley   /* Handle leaves: shared with the root point */
570b0a623aaSMatthew G. Knepley   for (l = 0; l < nleaves; ++l) {
571b0a623aaSMatthew G. Knepley     PetscInt adjSize = PETSC_DETERMINE, a;
572b0a623aaSMatthew G. Knepley 
573b0a623aaSMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, local[l], &adjSize, &adj);CHKERRQ(ierr);
574b0a623aaSMatthew G. Knepley     for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remote[l].rank);CHKERRQ(ierr);}
575b0a623aaSMatthew G. Knepley   }
576b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(rootrank, &rrank);CHKERRQ(ierr);
577b0a623aaSMatthew G. Knepley   ierr = ISGetIndices(leafrank, &nrank);CHKERRQ(ierr);
578b0a623aaSMatthew G. Knepley   /* Handle roots */
579b0a623aaSMatthew G. Knepley   for (p = pStart; p < pEnd; ++p) {
580b0a623aaSMatthew G. Knepley     PetscInt adjSize = PETSC_DETERMINE, neighbors = 0, noff, n, a;
581b0a623aaSMatthew G. Knepley 
582b0a623aaSMatthew G. Knepley     if ((p >= sStart) && (p < sEnd)) {
583b0a623aaSMatthew G. Knepley       /* Some leaves share a root with other leaves on different processes */
584b0a623aaSMatthew G. Knepley       ierr = PetscSectionGetDof(leafSection, p, &neighbors);CHKERRQ(ierr);
585b0a623aaSMatthew G. Knepley       if (neighbors) {
586b0a623aaSMatthew G. Knepley         ierr = PetscSectionGetOffset(leafSection, p, &noff);CHKERRQ(ierr);
587b0a623aaSMatthew G. Knepley         ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
588b0a623aaSMatthew G. Knepley         for (n = 0; n < neighbors; ++n) {
589b0a623aaSMatthew G. Knepley           const PetscInt remoteRank = nrank[noff+n];
590b0a623aaSMatthew G. Knepley 
591b0a623aaSMatthew G. Knepley           if (remoteRank == rank) continue;
592b0a623aaSMatthew G. Knepley           for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);}
593b0a623aaSMatthew G. Knepley         }
594b0a623aaSMatthew G. Knepley       }
595b0a623aaSMatthew G. Knepley     }
596b0a623aaSMatthew G. Knepley     /* Roots are shared with leaves */
597b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetDof(rootSection, p, &neighbors);CHKERRQ(ierr);
598b0a623aaSMatthew G. Knepley     if (!neighbors) continue;
599b0a623aaSMatthew G. Knepley     ierr = PetscSectionGetOffset(rootSection, p, &noff);CHKERRQ(ierr);
600b0a623aaSMatthew G. Knepley     ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr);
601b0a623aaSMatthew G. Knepley     for (n = 0; n < neighbors; ++n) {
602b0a623aaSMatthew G. Knepley       const PetscInt remoteRank = rrank[noff+n];
603b0a623aaSMatthew G. Knepley 
604b0a623aaSMatthew G. Knepley       if (remoteRank == rank) continue;
605b0a623aaSMatthew G. Knepley       for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);}
606b0a623aaSMatthew G. Knepley     }
607b0a623aaSMatthew G. Knepley   }
608b0a623aaSMatthew G. Knepley   ierr = PetscFree(adj);CHKERRQ(ierr);
609b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(rootrank, &rrank);CHKERRQ(ierr);
610b0a623aaSMatthew G. Knepley   ierr = ISRestoreIndices(leafrank, &nrank);CHKERRQ(ierr);
61124d039d7SMichael Lange   /* Add additional overlap levels */
61224d039d7SMichael Lange   for (l = 1; l < levels; l++) {ierr = DMPlexPartitionLabelAdjacency(dm, ovAdjByRank);CHKERRQ(ierr);}
61326a7d390SMatthew G. Knepley   /* We require the closure in the overlap */
61426a7d390SMatthew G. Knepley   ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr);
61526a7d390SMatthew G. Knepley   ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr);
61626a7d390SMatthew G. Knepley   if (useCone || !useClosure) {
6175abbe4feSMichael Lange     ierr = DMPlexPartitionLabelClosure(dm, ovAdjByRank);CHKERRQ(ierr);
61826a7d390SMatthew G. Knepley   }
619e540f424SMichael Lange   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-overlap_view", &flg);CHKERRQ(ierr);
620e540f424SMichael Lange   if (flg) {
621b0a623aaSMatthew G. Knepley     ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr);
622b0a623aaSMatthew G. Knepley     ierr = DMLabelView(ovAdjByRank, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
623b0a623aaSMatthew G. Knepley   }
6241fd9873aSMichael Lange   /* Make process SF and invert sender to receiver label */
625b0a623aaSMatthew G. Knepley   ierr = DMPlexCreateTwoSidedProcessSF(dm, sfPoint, rootSection, rootrank, leafSection, leafrank, NULL, &sfProc);CHKERRQ(ierr);
626e540f424SMichael Lange   ierr = DMLabelCreate("ovLeafs", &ovLeafLabel);CHKERRQ(ierr);
6271fd9873aSMichael Lange   ierr = DMPlexPartitionLabelInvert(dm, ovAdjByRank, sfProc, ovLeafLabel);CHKERRQ(ierr);
628e540f424SMichael Lange   /* Don't import points from yourself */
6291fd9873aSMichael Lange   ierr = DMLabelClearStratum(ovLeafLabel, rank);CHKERRQ(ierr);
630e540f424SMichael Lange   /* Don't import points already in the pointSF */
631e540f424SMichael Lange   for (l = 0; l < nleaves; ++l) {
632e540f424SMichael Lange     ierr = DMLabelClearValue(ovLeafLabel, remote[l].index, remote[l].rank);CHKERRQ(ierr);
633e540f424SMichael Lange   }
634aa3148a8SMichael Lange   /* Convert receiver label to SF */
635aa3148a8SMichael Lange   ierr = DMPlexPartitionLabelCreateSF(dm, ovLeafLabel, overlapSF);CHKERRQ(ierr);
636e540f424SMichael Lange   ierr = PetscObjectSetName((PetscObject) *overlapSF, "Overlap SF");CHKERRQ(ierr);
637e540f424SMichael Lange   ierr = PetscSFSetFromOptions(*overlapSF);CHKERRQ(ierr);
638e540f424SMichael Lange   if (flg) {
639e540f424SMichael Lange     ierr = PetscPrintf(comm, "Overlap SF\n");CHKERRQ(ierr);
640e540f424SMichael Lange     ierr = PetscSFView(*overlapSF, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
641e540f424SMichael Lange   }
642e540f424SMichael Lange   /* Clean up */
6431fd9873aSMichael Lange   ierr = DMLabelDestroy(&ovAdjByRank);CHKERRQ(ierr);
644aa3148a8SMichael Lange   ierr = DMLabelDestroy(&ovLeafLabel);CHKERRQ(ierr);
645b0a623aaSMatthew G. Knepley   ierr = PetscSFDestroy(&sfProc);CHKERRQ(ierr);
646b0a623aaSMatthew G. Knepley   PetscFunctionReturn(0);
647b0a623aaSMatthew G. Knepley }
64870034214SMatthew G. Knepley 
64970034214SMatthew G. Knepley #undef __FUNCT__
65046f9b1c3SMichael Lange #define __FUNCT__ "DMPlexCreateOverlapMigrationSF"
65146f9b1c3SMichael Lange PetscErrorCode DMPlexCreateOverlapMigrationSF(DM dm, PetscSF overlapSF, PetscSF *migrationSF)
65246f9b1c3SMichael Lange {
65346f9b1c3SMichael Lange   MPI_Comm           comm;
65446f9b1c3SMichael Lange   PetscMPIInt        rank, numProcs;
65546f9b1c3SMichael Lange   PetscInt           d, dim, p, pStart, pEnd, nroots, nleaves, newLeaves, point, numSharedPoints;
65646f9b1c3SMichael Lange   PetscInt          *pointDepths, *remoteDepths, *ilocal;
65746f9b1c3SMichael Lange   PetscInt          *depthRecv, *depthShift, *depthIdx;
65846f9b1c3SMichael Lange   PetscSFNode       *iremote;
65946f9b1c3SMichael Lange   PetscSF            pointSF;
66046f9b1c3SMichael Lange   const PetscInt    *sharedLocal;
66146f9b1c3SMichael Lange   const PetscSFNode *overlapRemote, *sharedRemote;
66246f9b1c3SMichael Lange   PetscErrorCode     ierr;
66346f9b1c3SMichael Lange 
66446f9b1c3SMichael Lange   PetscFunctionBegin;
66546f9b1c3SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
66646f9b1c3SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
66746f9b1c3SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
66846f9b1c3SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
66946f9b1c3SMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
67046f9b1c3SMichael Lange 
67146f9b1c3SMichael Lange   /* Before building the migration SF we need to know the new stratum offsets */
67246f9b1c3SMichael Lange   ierr = PetscSFGetGraph(overlapSF, &nroots, &nleaves, NULL, &overlapRemote);CHKERRQ(ierr);
67346f9b1c3SMichael Lange   ierr = PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);CHKERRQ(ierr);
67446f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
67546f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
67646f9b1c3SMichael Lange     for (p=pStart; p<pEnd; p++) pointDepths[p] = d;
67746f9b1c3SMichael Lange   }
67846f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) remoteDepths[p] = -1;
67946f9b1c3SMichael Lange   ierr = PetscSFBcastBegin(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
68046f9b1c3SMichael Lange   ierr = PetscSFBcastEnd(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr);
68146f9b1c3SMichael Lange 
68246f9b1c3SMichael Lange   /* Count recevied points in each stratum and compute the internal strata shift */
68346f9b1c3SMichael Lange   ierr = PetscMalloc3(dim+1, &depthRecv, dim+1, &depthShift, dim+1, &depthIdx);CHKERRQ(ierr);
68446f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) depthRecv[d]=0;
68546f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) depthRecv[remoteDepths[p]]++;
68646f9b1c3SMichael Lange   depthShift[dim] = 0;
68746f9b1c3SMichael Lange   for (d=0; d<dim; d++) depthShift[d] = depthRecv[dim];
68846f9b1c3SMichael Lange   for (d=1; d<dim; d++) depthShift[d] += depthRecv[0];
68946f9b1c3SMichael Lange   for (d=dim-2; d>0; d--) depthShift[d] += depthRecv[d+1];
69046f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
69146f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
69246f9b1c3SMichael Lange     depthIdx[d] = pStart + depthShift[d];
69346f9b1c3SMichael Lange   }
69446f9b1c3SMichael Lange 
69546f9b1c3SMichael Lange   /* Form the overlap SF build an SF that describes the full overlap migration SF */
69646f9b1c3SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
69746f9b1c3SMichael Lange   newLeaves = pEnd - pStart + nleaves;
69809b7985cSMatthew G. Knepley   ierr = PetscMalloc1(newLeaves, &ilocal);CHKERRQ(ierr);
69909b7985cSMatthew G. Knepley   ierr = PetscMalloc1(newLeaves, &iremote);CHKERRQ(ierr);
70046f9b1c3SMichael Lange   /* First map local points to themselves */
70146f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
70246f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
70346f9b1c3SMichael Lange     for (p=pStart; p<pEnd; p++) {
70446f9b1c3SMichael Lange       point = p + depthShift[d];
70546f9b1c3SMichael Lange       ilocal[point] = point;
70646f9b1c3SMichael Lange       iremote[point].index = p;
70746f9b1c3SMichael Lange       iremote[point].rank = rank;
70846f9b1c3SMichael Lange       depthIdx[d]++;
70946f9b1c3SMichael Lange     }
71046f9b1c3SMichael Lange   }
71146f9b1c3SMichael Lange 
71246f9b1c3SMichael Lange   /* Add in the remote roots for currently shared points */
71346f9b1c3SMichael Lange   ierr = DMGetPointSF(dm, &pointSF);CHKERRQ(ierr);
71446f9b1c3SMichael Lange   ierr = PetscSFGetGraph(pointSF, NULL, &numSharedPoints, &sharedLocal, &sharedRemote);CHKERRQ(ierr);
71546f9b1c3SMichael Lange   for (d=0; d<dim+1; d++) {
71646f9b1c3SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr);
71746f9b1c3SMichael Lange     for (p=0; p<numSharedPoints; p++) {
71846f9b1c3SMichael Lange       if (pStart <= sharedLocal[p] && sharedLocal[p] < pEnd) {
71946f9b1c3SMichael Lange         point = sharedLocal[p] + depthShift[d];
72046f9b1c3SMichael Lange         iremote[point].index = sharedRemote[p].index;
72146f9b1c3SMichael Lange         iremote[point].rank = sharedRemote[p].rank;
72246f9b1c3SMichael Lange       }
72346f9b1c3SMichael Lange     }
72446f9b1c3SMichael Lange   }
72546f9b1c3SMichael Lange 
72646f9b1c3SMichael Lange   /* Now add the incoming overlap points */
72746f9b1c3SMichael Lange   for (p=0; p<nleaves; p++) {
72846f9b1c3SMichael Lange     point = depthIdx[remoteDepths[p]];
72946f9b1c3SMichael Lange     ilocal[point] = point;
73046f9b1c3SMichael Lange     iremote[point].index = overlapRemote[p].index;
73146f9b1c3SMichael Lange     iremote[point].rank = overlapRemote[p].rank;
73246f9b1c3SMichael Lange     depthIdx[remoteDepths[p]]++;
73346f9b1c3SMichael Lange   }
73415fff7beSMatthew G. Knepley   ierr = PetscFree2(pointDepths,remoteDepths);CHKERRQ(ierr);
73546f9b1c3SMichael Lange 
73646f9b1c3SMichael Lange   ierr = PetscSFCreate(comm, migrationSF);CHKERRQ(ierr);
73746f9b1c3SMichael Lange   ierr = PetscObjectSetName((PetscObject) *migrationSF, "Overlap Migration SF");CHKERRQ(ierr);
73846f9b1c3SMichael Lange   ierr = PetscSFSetFromOptions(*migrationSF);CHKERRQ(ierr);
73946f9b1c3SMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
74046f9b1c3SMichael Lange   ierr = PetscSFSetGraph(*migrationSF, pEnd-pStart, newLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr);
74146f9b1c3SMichael Lange 
74246f9b1c3SMichael Lange   ierr = PetscFree3(depthRecv, depthShift, depthIdx);CHKERRQ(ierr);
74370034214SMatthew G. Knepley   PetscFunctionReturn(0);
74470034214SMatthew G. Knepley }
74570034214SMatthew G. Knepley 
74670034214SMatthew G. Knepley #undef __FUNCT__
74770034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeField"
74870034214SMatthew G. Knepley /*@
74970034214SMatthew G. Knepley   DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
75070034214SMatthew G. Knepley 
75170034214SMatthew G. Knepley   Collective on DM
75270034214SMatthew G. Knepley 
75370034214SMatthew G. Knepley   Input Parameters:
75470034214SMatthew G. Knepley + dm - The DMPlex object
75570034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
75670034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
75770034214SMatthew G. Knepley - originalVec - The existing data
75870034214SMatthew G. Knepley 
75970034214SMatthew G. Knepley   Output Parameters:
76070034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout
76170034214SMatthew G. Knepley - newVec - The new data
76270034214SMatthew G. Knepley 
76370034214SMatthew G. Knepley   Level: developer
76470034214SMatthew G. Knepley 
76570034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeFieldIS(), DMPlexDistributeData()
76670034214SMatthew G. Knepley @*/
76770034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec)
76870034214SMatthew G. Knepley {
76970034214SMatthew G. Knepley   PetscSF        fieldSF;
77070034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
77170034214SMatthew G. Knepley   PetscScalar   *originalValues, *newValues;
77270034214SMatthew G. Knepley   PetscErrorCode ierr;
77370034214SMatthew G. Knepley 
77470034214SMatthew G. Knepley   PetscFunctionBegin;
77570034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
77670034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
77770034214SMatthew G. Knepley 
77870034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
77970034214SMatthew G. Knepley   ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr);
78070034214SMatthew G. Knepley   ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr);
78170034214SMatthew G. Knepley 
78270034214SMatthew G. Knepley   ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr);
78370034214SMatthew G. Knepley   ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr);
78470034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
78570034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
78670034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
78770034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
78870034214SMatthew G. Knepley   ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr);
78970034214SMatthew G. Knepley   ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr);
79070034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
79170034214SMatthew G. Knepley   PetscFunctionReturn(0);
79270034214SMatthew G. Knepley }
79370034214SMatthew G. Knepley 
79470034214SMatthew G. Knepley #undef __FUNCT__
79570034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeFieldIS"
79670034214SMatthew G. Knepley /*@
79770034214SMatthew G. Knepley   DMPlexDistributeFieldIS - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
79870034214SMatthew G. Knepley 
79970034214SMatthew G. Knepley   Collective on DM
80070034214SMatthew G. Knepley 
80170034214SMatthew G. Knepley   Input Parameters:
80270034214SMatthew G. Knepley + dm - The DMPlex object
80370034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
80470034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
80570034214SMatthew G. Knepley - originalIS - The existing data
80670034214SMatthew G. Knepley 
80770034214SMatthew G. Knepley   Output Parameters:
80870034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout
80970034214SMatthew G. Knepley - newIS - The new data
81070034214SMatthew G. Knepley 
81170034214SMatthew G. Knepley   Level: developer
81270034214SMatthew G. Knepley 
81370034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexDistributeData()
81470034214SMatthew G. Knepley @*/
81570034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeFieldIS(DM dm, PetscSF pointSF, PetscSection originalSection, IS originalIS, PetscSection newSection, IS *newIS)
81670034214SMatthew G. Knepley {
81770034214SMatthew G. Knepley   PetscSF         fieldSF;
81870034214SMatthew G. Knepley   PetscInt       *newValues, *remoteOffsets, fieldSize;
81970034214SMatthew G. Knepley   const PetscInt *originalValues;
82070034214SMatthew G. Knepley   PetscErrorCode  ierr;
82170034214SMatthew G. Knepley 
82270034214SMatthew G. Knepley   PetscFunctionBegin;
82370034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
82470034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
82570034214SMatthew G. Knepley 
82670034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
82770034214SMatthew G. Knepley   ierr = PetscMalloc(fieldSize * sizeof(PetscInt), &newValues);CHKERRQ(ierr);
82870034214SMatthew G. Knepley 
82970034214SMatthew G. Knepley   ierr = ISGetIndices(originalIS, &originalValues);CHKERRQ(ierr);
83070034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
831aaf8c182SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr);
832aaf8c182SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr);
83370034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
83470034214SMatthew G. Knepley   ierr = ISRestoreIndices(originalIS, &originalValues);CHKERRQ(ierr);
83570034214SMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) pointSF), fieldSize, newValues, PETSC_OWN_POINTER, newIS);CHKERRQ(ierr);
83670034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
83770034214SMatthew G. Knepley   PetscFunctionReturn(0);
83870034214SMatthew G. Knepley }
83970034214SMatthew G. Knepley 
84070034214SMatthew G. Knepley #undef __FUNCT__
84170034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeData"
84270034214SMatthew G. Knepley /*@
84370034214SMatthew G. Knepley   DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
84470034214SMatthew G. Knepley 
84570034214SMatthew G. Knepley   Collective on DM
84670034214SMatthew G. Knepley 
84770034214SMatthew G. Knepley   Input Parameters:
84870034214SMatthew G. Knepley + dm - The DMPlex object
84970034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
85070034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
85170034214SMatthew G. Knepley . datatype - The type of data
85270034214SMatthew G. Knepley - originalData - The existing data
85370034214SMatthew G. Knepley 
85470034214SMatthew G. Knepley   Output Parameters:
85570034214SMatthew G. Knepley + newSection - The PetscSection describing the new data layout
85670034214SMatthew G. Knepley - newData - The new data
85770034214SMatthew G. Knepley 
85870034214SMatthew G. Knepley   Level: developer
85970034214SMatthew G. Knepley 
86070034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField()
86170034214SMatthew G. Knepley @*/
86270034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData)
86370034214SMatthew G. Knepley {
86470034214SMatthew G. Knepley   PetscSF        fieldSF;
86570034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
86670034214SMatthew G. Knepley   PetscMPIInt    dataSize;
86770034214SMatthew G. Knepley   PetscErrorCode ierr;
86870034214SMatthew G. Knepley 
86970034214SMatthew G. Knepley   PetscFunctionBegin;
87070034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
87170034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
87270034214SMatthew G. Knepley 
87370034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
87470034214SMatthew G. Knepley   ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr);
87570034214SMatthew G. Knepley   ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr);
87670034214SMatthew G. Knepley 
87770034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
87870034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
87970034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
88070034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
88170034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
88270034214SMatthew G. Knepley   PetscFunctionReturn(0);
88370034214SMatthew G. Knepley }
88470034214SMatthew G. Knepley 
88570034214SMatthew G. Knepley #undef __FUNCT__
886cc71bff1SMichael Lange #define __FUNCT__ "DMPlexDistributeCones"
887cc71bff1SMichael Lange PetscErrorCode DMPlexDistributeCones(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
888cc71bff1SMichael Lange {
889cc71bff1SMichael Lange   DM_Plex               *pmesh   = (DM_Plex*) (dmParallel)->data;
890cc71bff1SMichael Lange   MPI_Comm               comm;
891cc71bff1SMichael Lange   PetscSF                coneSF;
892cc71bff1SMichael Lange   PetscSection           originalConeSection, newConeSection;
893cc71bff1SMichael Lange   PetscInt              *remoteOffsets, *cones, *newCones, newConesSize;
894cc71bff1SMichael Lange   PetscBool              flg;
895cc71bff1SMichael Lange   PetscErrorCode         ierr;
896cc71bff1SMichael Lange 
897cc71bff1SMichael Lange   PetscFunctionBegin;
898cc71bff1SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
899cc71bff1SMichael Lange   PetscValidPointer(dmParallel,4);
900cc71bff1SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
901cc71bff1SMichael Lange 
902cc71bff1SMichael Lange   /* Distribute cone section */
903cc71bff1SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
904cc71bff1SMichael Lange   ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr);
905cc71bff1SMichael Lange   ierr = DMPlexGetConeSection(dmParallel, &newConeSection);CHKERRQ(ierr);
906cc71bff1SMichael Lange   ierr = PetscSFDistributeSection(migrationSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr);
907cc71bff1SMichael Lange   ierr = DMSetUp(dmParallel);CHKERRQ(ierr);
908cc71bff1SMichael Lange   {
909cc71bff1SMichael Lange     PetscInt pStart, pEnd, p;
910cc71bff1SMichael Lange 
911cc71bff1SMichael Lange     ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr);
912cc71bff1SMichael Lange     for (p = pStart; p < pEnd; ++p) {
913cc71bff1SMichael Lange       PetscInt coneSize;
914cc71bff1SMichael Lange       ierr               = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr);
915cc71bff1SMichael Lange       pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize);
916cc71bff1SMichael Lange     }
917cc71bff1SMichael Lange   }
918cc71bff1SMichael Lange   /* Communicate and renumber cones */
919cc71bff1SMichael Lange   ierr = PetscSFCreateSectionSF(migrationSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr);
920cc71bff1SMichael Lange   ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
921cc71bff1SMichael Lange   ierr = DMPlexGetCones(dmParallel, &newCones);CHKERRQ(ierr);
922cc71bff1SMichael Lange   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
923cc71bff1SMichael Lange   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
924cc71bff1SMichael Lange   ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr);
925cc71bff1SMichael Lange   ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr);
9263533c52bSMatthew G. Knepley #if PETSC_USE_DEBUG
9273533c52bSMatthew G. Knepley   {
9283533c52bSMatthew G. Knepley     PetscInt  p;
9293533c52bSMatthew G. Knepley     PetscBool valid = PETSC_TRUE;
9303533c52bSMatthew G. Knepley     for (p = 0; p < newConesSize; ++p) {
9313533c52bSMatthew G. Knepley       if (newCones[p] < 0) {valid = PETSC_FALSE; ierr = PetscPrintf(PETSC_COMM_SELF, "Point %d not in overlap SF\n", p);CHKERRQ(ierr);}
9323533c52bSMatthew G. Knepley     }
9333533c52bSMatthew G. Knepley     if (!valid) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map");
9343533c52bSMatthew G. Knepley   }
9353533c52bSMatthew G. Knepley #endif
936cc71bff1SMichael Lange   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr);
937cc71bff1SMichael Lange   if (flg) {
938cc71bff1SMichael Lange     ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr);
939cc71bff1SMichael Lange     ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
940cc71bff1SMichael Lange     ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr);
941cc71bff1SMichael Lange     ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
942cc71bff1SMichael Lange     ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr);
943cc71bff1SMichael Lange   }
944cc71bff1SMichael Lange   ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr);
945cc71bff1SMichael Lange   ierr = DMPlexGetConeOrientations(dmParallel, &newCones);CHKERRQ(ierr);
946cc71bff1SMichael Lange   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
947cc71bff1SMichael Lange   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
948cc71bff1SMichael Lange   ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr);
949cc71bff1SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
950cc71bff1SMichael Lange   /* Create supports and stratify sieve */
951cc71bff1SMichael Lange   {
952cc71bff1SMichael Lange     PetscInt pStart, pEnd;
953cc71bff1SMichael Lange 
954cc71bff1SMichael Lange     ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr);
955cc71bff1SMichael Lange     ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr);
956cc71bff1SMichael Lange   }
957cc71bff1SMichael Lange   ierr = DMPlexSymmetrize(dmParallel);CHKERRQ(ierr);
958cc71bff1SMichael Lange   ierr = DMPlexStratify(dmParallel);CHKERRQ(ierr);
959cc71bff1SMichael Lange   PetscFunctionReturn(0);
960cc71bff1SMichael Lange }
961cc71bff1SMichael Lange 
962cc71bff1SMichael Lange #undef __FUNCT__
9630df0e737SMichael Lange #define __FUNCT__ "DMPlexDistributeCoordinates"
9640df0e737SMichael Lange PetscErrorCode DMPlexDistributeCoordinates(DM dm, PetscSF migrationSF, DM dmParallel)
9650df0e737SMichael Lange {
9660df0e737SMichael Lange   MPI_Comm         comm;
9670df0e737SMichael Lange   PetscSection     originalCoordSection, newCoordSection;
9680df0e737SMichael Lange   Vec              originalCoordinates, newCoordinates;
9690df0e737SMichael Lange   PetscInt         bs;
9700df0e737SMichael Lange   const char      *name;
9710df0e737SMichael Lange   const PetscReal *maxCell, *L;
9720df0e737SMichael Lange   PetscErrorCode   ierr;
9730df0e737SMichael Lange 
9740df0e737SMichael Lange   PetscFunctionBegin;
9750df0e737SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9760df0e737SMichael Lange   PetscValidPointer(dmParallel, 3);
9770df0e737SMichael Lange 
9780df0e737SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
9790df0e737SMichael Lange   ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr);
9800df0e737SMichael Lange   ierr = DMGetCoordinateSection(dmParallel, &newCoordSection);CHKERRQ(ierr);
9810df0e737SMichael Lange   ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr);
9820df0e737SMichael Lange   if (originalCoordinates) {
9830df0e737SMichael Lange     ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr);
9840df0e737SMichael Lange     ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr);
9850df0e737SMichael Lange     ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr);
9860df0e737SMichael Lange 
9870df0e737SMichael Lange     ierr = DMPlexDistributeField(dm, migrationSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr);
9880df0e737SMichael Lange     ierr = DMSetCoordinatesLocal(dmParallel, newCoordinates);CHKERRQ(ierr);
9890df0e737SMichael Lange     ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr);
9900df0e737SMichael Lange     ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr);
9910df0e737SMichael Lange     ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr);
9920df0e737SMichael Lange   }
9930df0e737SMichael Lange   ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr);
9940df0e737SMichael Lange   if (L) {ierr = DMSetPeriodicity(dmParallel, maxCell, L);CHKERRQ(ierr);}
9950df0e737SMichael Lange   PetscFunctionReturn(0);
9960df0e737SMichael Lange }
9970df0e737SMichael Lange 
9980df0e737SMichael Lange #undef __FUNCT__
9990df0e737SMichael Lange #define __FUNCT__ "DMPlexDistributeLabels"
1000d995df53SMatthew G. Knepley /* Here we are assuming that process 0 always has everything */
10012eb1fa14SMichael Lange PetscErrorCode DMPlexDistributeLabels(DM dm, PetscSF migrationSF, DM dmParallel)
10020df0e737SMichael Lange {
10030df0e737SMichael Lange   MPI_Comm       comm;
10040df0e737SMichael Lange   PetscMPIInt    rank;
1005d995df53SMatthew G. Knepley   PetscInt       numLabels, numLocalLabels, l;
1006d995df53SMatthew G. Knepley   PetscBool      hasLabels = PETSC_FALSE;
10070df0e737SMichael Lange   PetscErrorCode ierr;
10080df0e737SMichael Lange 
10090df0e737SMichael Lange   PetscFunctionBegin;
10100df0e737SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10112eb1fa14SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 3);
10120df0e737SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
10130df0e737SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
10140df0e737SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
10150df0e737SMichael Lange 
1016d995df53SMatthew G. Knepley   /* Everyone must have either the same number of labels, or none */
1017d995df53SMatthew G. Knepley   ierr = DMPlexGetNumLabels(dm, &numLocalLabels);CHKERRQ(ierr);
1018d995df53SMatthew G. Knepley   numLabels = numLocalLabels;
10190df0e737SMichael Lange   ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
1020627847f0SMatthew G. Knepley   if (numLabels == numLocalLabels) hasLabels = PETSC_TRUE;
1021bdd2d751SMichael Lange   for (l = numLabels-1; l >= 0; --l) {
1022bdd2d751SMichael Lange     DMLabel     label = NULL, labelNew = NULL;
10230df0e737SMichael Lange     PetscBool   isdepth;
10240df0e737SMichael Lange 
1025d995df53SMatthew G. Knepley     if (hasLabels) {
1026bdd2d751SMichael Lange       ierr = DMPlexGetLabelByNum(dm, l, &label);CHKERRQ(ierr);
10270df0e737SMichael Lange       /* Skip "depth" because it is recreated */
1028bdd2d751SMichael Lange       ierr = PetscStrcmp(label->name, "depth", &isdepth);CHKERRQ(ierr);
1029bdd2d751SMichael Lange     }
10300df0e737SMichael Lange     ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
1031bdd2d751SMichael Lange     if (isdepth) continue;
1032bdd2d751SMichael Lange     ierr = DMLabelDistribute(label, migrationSF, &labelNew);CHKERRQ(ierr);
1033bdd2d751SMichael Lange     ierr = DMPlexAddLabel(dmParallel, labelNew);CHKERRQ(ierr);
10340df0e737SMichael Lange   }
10350df0e737SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
10360df0e737SMichael Lange   PetscFunctionReturn(0);
10370df0e737SMichael Lange }
10380df0e737SMichael Lange 
10399734c634SMichael Lange #undef __FUNCT__
10409734c634SMichael Lange #define __FUNCT__ "DMPlexDistributeSetupHybrid"
10419734c634SMichael Lange PetscErrorCode DMPlexDistributeSetupHybrid(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
10429734c634SMichael Lange {
10439734c634SMichael Lange   DM_Plex        *mesh  = (DM_Plex*) dm->data;
10449734c634SMichael Lange   DM_Plex        *pmesh = (DM_Plex*) (dmParallel)->data;
10459734c634SMichael Lange   MPI_Comm        comm;
10469734c634SMichael Lange   const PetscInt *gpoints;
10479734c634SMichael Lange   PetscInt        dim, depth, n, d;
10489734c634SMichael Lange   PetscErrorCode  ierr;
10499734c634SMichael Lange 
10509734c634SMichael Lange   PetscFunctionBegin;
10519734c634SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
10529734c634SMichael Lange   PetscValidPointer(dmParallel, 4);
10539734c634SMichael Lange 
10549734c634SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
10559734c634SMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
10569734c634SMichael Lange 
10579734c634SMichael Lange   /* Setup hybrid structure */
10589734c634SMichael Lange   for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];}
10599734c634SMichael Lange   ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr);
10609734c634SMichael Lange   ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr);
10619734c634SMichael Lange   ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr);
10629734c634SMichael Lange   ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
10639734c634SMichael Lange   for (d = 0; d <= dim; ++d) {
10649734c634SMichael Lange     PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p;
10659734c634SMichael Lange 
10669734c634SMichael Lange     if (pmax < 0) continue;
10679734c634SMichael Lange     ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr);
10689734c634SMichael Lange     ierr = DMPlexGetDepthStratum(dmParallel, d, NULL, &pEnd);CHKERRQ(ierr);
10699734c634SMichael Lange     ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr);
10709734c634SMichael Lange     for (p = 0; p < n; ++p) {
10719734c634SMichael Lange       const PetscInt point = gpoints[p];
10729734c634SMichael Lange 
10739734c634SMichael Lange       if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax;
10749734c634SMichael Lange     }
10759734c634SMichael Lange     if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax;
10769734c634SMichael Lange     else            pmesh->hybridPointMax[d] = -1;
10779734c634SMichael Lange   }
10789734c634SMichael Lange   ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr);
10799734c634SMichael Lange   PetscFunctionReturn(0);
10809734c634SMichael Lange }
10810df0e737SMichael Lange 
1082a6f36705SMichael Lange #undef __FUNCT__
1083a6f36705SMichael Lange #define __FUNCT__ "DMPlexDistributeSetupTree"
1084a6f36705SMichael Lange PetscErrorCode DMPlexDistributeSetupTree(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel)
1085a6f36705SMichael Lange {
108615078cd4SMichael Lange   DM_Plex        *mesh  = (DM_Plex*) dm->data;
108715078cd4SMichael Lange   DM_Plex        *pmesh = (DM_Plex*) (dmParallel)->data;
1088a6f36705SMichael Lange   MPI_Comm        comm;
1089a6f36705SMichael Lange   DM              refTree;
1090a6f36705SMichael Lange   PetscSection    origParentSection, newParentSection;
1091a6f36705SMichael Lange   PetscInt        *origParents, *origChildIDs;
1092a6f36705SMichael Lange   PetscBool       flg;
1093a6f36705SMichael Lange   PetscErrorCode  ierr;
1094a6f36705SMichael Lange 
1095a6f36705SMichael Lange   PetscFunctionBegin;
1096a6f36705SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1097a6f36705SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 4);
1098a6f36705SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
1099a6f36705SMichael Lange 
1100a6f36705SMichael Lange   /* Set up tree */
1101a6f36705SMichael Lange   ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr);
1102a6f36705SMichael Lange   ierr = DMPlexSetReferenceTree(dmParallel,refTree);CHKERRQ(ierr);
1103a6f36705SMichael Lange   ierr = DMPlexGetTree(dm,&origParentSection,&origParents,&origChildIDs,NULL,NULL);CHKERRQ(ierr);
1104a6f36705SMichael Lange   if (origParentSection) {
1105a6f36705SMichael Lange     PetscInt        pStart, pEnd;
1106a6f36705SMichael Lange     PetscInt        *newParents, *newChildIDs;
1107a6f36705SMichael Lange     PetscInt        *remoteOffsetsParents, newParentSize;
1108a6f36705SMichael Lange     PetscSF         parentSF;
1109a6f36705SMichael Lange 
1110a6f36705SMichael Lange     ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
1111a6f36705SMichael Lange     ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dmParallel),&newParentSection);CHKERRQ(ierr);
1112a6f36705SMichael Lange     ierr = PetscSectionSetChart(newParentSection,pStart,pEnd);CHKERRQ(ierr);
1113a6f36705SMichael Lange     ierr = PetscSFDistributeSection(migrationSF, origParentSection, &remoteOffsetsParents, newParentSection);CHKERRQ(ierr);
1114a6f36705SMichael Lange     ierr = PetscSFCreateSectionSF(migrationSF, origParentSection, remoteOffsetsParents, newParentSection, &parentSF);CHKERRQ(ierr);
1115a6f36705SMichael Lange     ierr = PetscSectionGetStorageSize(newParentSection,&newParentSize);CHKERRQ(ierr);
1116a6f36705SMichael Lange     ierr = PetscMalloc2(newParentSize,&newParents,newParentSize,&newChildIDs);CHKERRQ(ierr);
1117a6f36705SMichael Lange     ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
1118a6f36705SMichael Lange     ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr);
1119a6f36705SMichael Lange     ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1120a6f36705SMichael Lange     ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr);
1121a6f36705SMichael Lange     ierr = ISGlobalToLocalMappingApplyBlock(renumbering,IS_GTOLM_MASK, newParentSize, newParents, NULL, newParents);CHKERRQ(ierr);
1122a6f36705SMichael Lange     ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-parents_view", &flg);CHKERRQ(ierr);
1123a6f36705SMichael Lange     if (flg) {
1124a6f36705SMichael Lange       ierr = PetscPrintf(comm, "Serial Parent Section: \n");CHKERRQ(ierr);
1125a6f36705SMichael Lange       ierr = PetscSectionView(origParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1126a6f36705SMichael Lange       ierr = PetscPrintf(comm, "Parallel Parent Section: \n");CHKERRQ(ierr);
1127a6f36705SMichael Lange       ierr = PetscSectionView(newParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1128a6f36705SMichael Lange       ierr = PetscSFView(parentSF, NULL);CHKERRQ(ierr);
1129a6f36705SMichael Lange     }
1130a6f36705SMichael Lange     ierr = DMPlexSetTree(dmParallel,newParentSection,newParents,newChildIDs);CHKERRQ(ierr);
1131a6f36705SMichael Lange     ierr = PetscSectionDestroy(&newParentSection);CHKERRQ(ierr);
1132a6f36705SMichael Lange     ierr = PetscFree2(newParents,newChildIDs);CHKERRQ(ierr);
1133a6f36705SMichael Lange     ierr = PetscSFDestroy(&parentSF);CHKERRQ(ierr);
1134a6f36705SMichael Lange   }
113515078cd4SMichael Lange   pmesh->useAnchors = mesh->useAnchors;
1136a6f36705SMichael Lange   PetscFunctionReturn(0);
1137a6f36705SMichael Lange }
11380df0e737SMichael Lange 
11390df0e737SMichael Lange #undef __FUNCT__
11404eca1733SMichael Lange #define __FUNCT__ "DMPlexDistributeSF"
1141f8987ae8SMichael Lange PetscErrorCode DMPlexDistributeSF(DM dm, PetscSF migrationSF, DM dmParallel)
11424eca1733SMichael Lange {
11434eca1733SMichael Lange   DM_Plex               *mesh  = (DM_Plex*) dm->data;
11444eca1733SMichael Lange   DM_Plex               *pmesh = (DM_Plex*) (dmParallel)->data;
11454eca1733SMichael Lange   PetscMPIInt            rank, numProcs;
11464eca1733SMichael Lange   MPI_Comm               comm;
11474eca1733SMichael Lange   PetscErrorCode         ierr;
11484eca1733SMichael Lange 
11494eca1733SMichael Lange   PetscFunctionBegin;
11504eca1733SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11514eca1733SMichael Lange   PetscValidPointer(dmParallel,7);
11524eca1733SMichael Lange 
11534eca1733SMichael Lange   /* Create point SF for parallel mesh */
11544eca1733SMichael Lange   ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
11554eca1733SMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr);
11564eca1733SMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
11574eca1733SMichael Lange   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
11584eca1733SMichael Lange   {
11594eca1733SMichael Lange     const PetscInt *leaves;
11604eca1733SMichael Lange     PetscSFNode    *remotePoints, *rowners, *lowners;
11614eca1733SMichael Lange     PetscInt        numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints;
11624eca1733SMichael Lange     PetscInt        pStart, pEnd;
11634eca1733SMichael Lange 
11644eca1733SMichael Lange     ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
11654eca1733SMichael Lange     ierr = PetscSFGetGraph(migrationSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr);
11664eca1733SMichael Lange     ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr);
11674eca1733SMichael Lange     for (p=0; p<numRoots; p++) {
11684eca1733SMichael Lange       rowners[p].rank  = -1;
11694eca1733SMichael Lange       rowners[p].index = -1;
11704eca1733SMichael Lange     }
11714eca1733SMichael Lange     ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
11724eca1733SMichael Lange     ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
11734eca1733SMichael Lange     for (p = 0; p < numLeaves; ++p) {
11744eca1733SMichael Lange       if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */
11754eca1733SMichael Lange         lowners[p].rank  = rank;
11764eca1733SMichael Lange         lowners[p].index = leaves ? leaves[p] : p;
11774eca1733SMichael Lange       } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */
11784eca1733SMichael Lange         lowners[p].rank  = -2;
11794eca1733SMichael Lange         lowners[p].index = -2;
11804eca1733SMichael Lange       }
11814eca1733SMichael Lange     }
11824eca1733SMichael Lange     for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */
11834eca1733SMichael Lange       rowners[p].rank  = -3;
11844eca1733SMichael Lange       rowners[p].index = -3;
11854eca1733SMichael Lange     }
11864eca1733SMichael Lange     ierr = PetscSFReduceBegin(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
11874eca1733SMichael Lange     ierr = PetscSFReduceEnd(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
11884eca1733SMichael Lange     ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
11894eca1733SMichael Lange     ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
11904eca1733SMichael Lange     for (p = 0; p < numLeaves; ++p) {
11914eca1733SMichael Lange       if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed");
11924eca1733SMichael Lange       if (lowners[p].rank != rank) ++numGhostPoints;
11934eca1733SMichael Lange     }
11944eca1733SMichael Lange     ierr = PetscMalloc1(numGhostPoints, &ghostPoints);CHKERRQ(ierr);
11954eca1733SMichael Lange     ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr);
11964eca1733SMichael Lange     for (p = 0, gp = 0; p < numLeaves; ++p) {
11974eca1733SMichael Lange       if (lowners[p].rank != rank) {
11984eca1733SMichael Lange         ghostPoints[gp]        = leaves ? leaves[p] : p;
11994eca1733SMichael Lange         remotePoints[gp].rank  = lowners[p].rank;
12004eca1733SMichael Lange         remotePoints[gp].index = lowners[p].index;
12014eca1733SMichael Lange         ++gp;
12024eca1733SMichael Lange       }
12034eca1733SMichael Lange     }
12044eca1733SMichael Lange     ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr);
12054eca1733SMichael Lange     ierr = PetscSFSetGraph((dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
12064eca1733SMichael Lange     ierr = PetscSFSetFromOptions((dmParallel)->sf);CHKERRQ(ierr);
12074eca1733SMichael Lange   }
12084eca1733SMichael Lange   pmesh->useCone    = mesh->useCone;
12094eca1733SMichael Lange   pmesh->useClosure = mesh->useClosure;
12104eca1733SMichael Lange   ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
12114eca1733SMichael Lange   PetscFunctionReturn(0);
12124eca1733SMichael Lange }
12134eca1733SMichael Lange 
12144eca1733SMichael Lange #undef __FUNCT__
121515078cd4SMichael Lange #define __FUNCT__ "DMPlexMigrate"
121615078cd4SMichael Lange /*@C
121715078cd4SMichael Lange   DMPlexMigrate  - Migrates internal DM data over the supplied star forest
121815078cd4SMichael Lange 
121915078cd4SMichael Lange   Input Parameter:
122015078cd4SMichael Lange + dm       - The source DMPlex object
122115078cd4SMichael Lange . sf       - The start forest communication context describing the migration pattern
122215078cd4SMichael Lange 
122315078cd4SMichael Lange   Output Parameter:
122415078cd4SMichael Lange - targetDM - The target DMPlex object
122515078cd4SMichael Lange 
1226*b9f40539SMichael Lange   Level: intermediate
122715078cd4SMichael Lange 
122815078cd4SMichael Lange .seealso: DMPlexDistribute(), DMPlexDistributeOverlap()
122915078cd4SMichael Lange @*/
1230*b9f40539SMichael Lange PetscErrorCode DMPlexMigrate(DM dm, PetscSF sf, DM targetDM)
123115078cd4SMichael Lange {
1232*b9f40539SMichael Lange   MPI_Comm               comm;
1233*b9f40539SMichael Lange   PetscInt               dim, nroots;
1234*b9f40539SMichael Lange   PetscSF                sfPoint;
123515078cd4SMichael Lange   ISLocalToGlobalMapping ltogMigration;
1236*b9f40539SMichael Lange   PetscBool              flg;
123715078cd4SMichael Lange   PetscErrorCode         ierr;
123815078cd4SMichael Lange 
123915078cd4SMichael Lange   PetscFunctionBegin;
124015078cd4SMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1241*b9f40539SMichael Lange   ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr);
124215078cd4SMichael Lange   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
124315078cd4SMichael Lange   ierr = DMSetDimension(targetDM, dim);CHKERRQ(ierr);
124415078cd4SMichael Lange 
1245*b9f40539SMichael Lange   ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr);
1246*b9f40539SMichael Lange   ierr = PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL);CHKERRQ(ierr);
1247*b9f40539SMichael Lange   if (nroots > 0) {
1248*b9f40539SMichael Lange     IS                     isOriginal;
1249*b9f40539SMichael Lange     ISLocalToGlobalMapping ltogOriginal;
1250*b9f40539SMichael Lange     PetscInt               n, size, nleaves, conesSize;
1251*b9f40539SMichael Lange     PetscInt              *numbering_orig, *numbering_new, *cones;
125215078cd4SMichael Lange     PetscSection           coneSection;
1253*b9f40539SMichael Lange     /* Get the original point numbering */
1254*b9f40539SMichael Lange     ierr = DMPlexCreatePointNumbering(dm, &isOriginal);CHKERRQ(ierr);
1255*b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingCreateIS(isOriginal, &ltogOriginal);CHKERRQ(ierr);
1256*b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingGetSize(ltogOriginal, &size);CHKERRQ(ierr);
1257*b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingGetBlockIndices(ltogOriginal, (const PetscInt**)&numbering_orig);CHKERRQ(ierr);
1258*b9f40539SMichael Lange     /* Convert to positive global numbers */
1259*b9f40539SMichael Lange     for (n=0; n<size; n++) {if (numbering_orig[n] < 0) numbering_orig[n] = -(numbering_orig[n]+1);}
1260*b9f40539SMichael Lange     /* Derive the new local-to-global mapping from the old one */
1261*b9f40539SMichael Lange     ierr = PetscSFGetGraph(sf, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr);
1262*b9f40539SMichael Lange     ierr = PetscMalloc1(nleaves, &numbering_new);CHKERRQ(ierr);
1263*b9f40539SMichael Lange     ierr = PetscSFBcastBegin(sf, MPIU_INT, (PetscInt *) numbering_orig, numbering_new);CHKERRQ(ierr);
1264*b9f40539SMichael Lange     ierr = PetscSFBcastEnd(sf, MPIU_INT, (PetscInt *) numbering_orig, numbering_new);CHKERRQ(ierr);
1265*b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingCreate(comm, 1, nleaves, (const PetscInt*) numbering_new, PETSC_OWN_POINTER, &ltogMigration);CHKERRQ(ierr);
1266*b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingRestoreIndices(ltogOriginal, (const PetscInt**)&numbering_orig);CHKERRQ(ierr);
126715078cd4SMichael Lange     /* Convert cones to global numbering before migrating them */
126815078cd4SMichael Lange     ierr = DMPlexGetConeSection(dm, &coneSection);CHKERRQ(ierr);
126915078cd4SMichael Lange     ierr = PetscSectionGetStorageSize(coneSection, &conesSize);CHKERRQ(ierr);
127015078cd4SMichael Lange     ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
1271*b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingApplyBlock(ltogOriginal, conesSize, cones, cones);CHKERRQ(ierr);
1272*b9f40539SMichael Lange     ierr = ISDestroy(&isOriginal);CHKERRQ(ierr);
1273*b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingDestroy(&ltogOriginal);
127415078cd4SMichael Lange   } else {
127515078cd4SMichael Lange     /* Assume zero-to-all-ranks pattern and derive LToG from SF */
127615078cd4SMichael Lange     ierr = ISLocalToGlobalMappingCreateSF(sf, 0, &ltogMigration);CHKERRQ(ierr);
127715078cd4SMichael Lange   }
1278*b9f40539SMichael Lange   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
1279*b9f40539SMichael Lange   if (flg) {
1280*b9f40539SMichael Lange     ierr = PetscPrintf(comm, "Point renumbering for DM migration:\n");CHKERRQ(ierr);
1281*b9f40539SMichael Lange     ierr = ISLocalToGlobalMappingView(ltogMigration, NULL);CHKERRQ(ierr);
1282*b9f40539SMichael Lange   }
128315078cd4SMichael Lange   /* Migrate DM data to target DM */
128415078cd4SMichael Lange   ierr = DMPlexDistributeCones(dm, sf, ltogMigration, targetDM);CHKERRQ(ierr);
128515078cd4SMichael Lange   ierr = DMPlexDistributeCoordinates(dm, sf, targetDM);CHKERRQ(ierr);
128615078cd4SMichael Lange   ierr = DMPlexDistributeLabels(dm, sf, targetDM);CHKERRQ(ierr);
128715078cd4SMichael Lange   ierr = DMPlexDistributeSetupHybrid(dm, sf, ltogMigration, targetDM);CHKERRQ(ierr);
128815078cd4SMichael Lange   ierr = DMPlexDistributeSetupTree(dm, sf, ltogMigration, targetDM);CHKERRQ(ierr);
128915078cd4SMichael Lange   ierr = DMPlexCopyBoundary(dm, targetDM);CHKERRQ(ierr);
1290*b9f40539SMichael Lange   ierr = ISLocalToGlobalMappingDestroy(&ltogMigration);
129115078cd4SMichael Lange   PetscFunctionReturn(0);
129215078cd4SMichael Lange }
129315078cd4SMichael Lange 
129415078cd4SMichael Lange #undef __FUNCT__
129570034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistribute"
129670034214SMatthew G. Knepley /*@C
129770034214SMatthew G. Knepley   DMPlexDistribute - Distributes the mesh and any associated sections.
129870034214SMatthew G. Knepley 
129970034214SMatthew G. Knepley   Not Collective
130070034214SMatthew G. Knepley 
130170034214SMatthew G. Knepley   Input Parameter:
130270034214SMatthew G. Knepley + dm  - The original DMPlex object
130370034214SMatthew G. Knepley - overlap - The overlap of partitions, 0 is the default
130470034214SMatthew G. Knepley 
130570034214SMatthew G. Knepley   Output Parameter:
130670034214SMatthew G. Knepley + sf - The PetscSF used for point distribution
130770034214SMatthew G. Knepley - parallelMesh - The distributed DMPlex object, or NULL
130870034214SMatthew G. Knepley 
130970034214SMatthew G. Knepley   Note: If the mesh was not distributed, the return value is NULL.
131070034214SMatthew G. Knepley 
131170034214SMatthew G. Knepley   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
131270034214SMatthew G. Knepley   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
131370034214SMatthew G. Knepley   representation on the mesh.
131470034214SMatthew G. Knepley 
131570034214SMatthew G. Knepley   Level: intermediate
131670034214SMatthew G. Knepley 
131770034214SMatthew G. Knepley .keywords: mesh, elements
131870034214SMatthew G. Knepley .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
131970034214SMatthew G. Knepley @*/
132080cf41d5SMatthew G. Knepley PetscErrorCode DMPlexDistribute(DM dm, PetscInt overlap, PetscSF *sf, DM *dmParallel)
132170034214SMatthew G. Knepley {
132270034214SMatthew G. Knepley   MPI_Comm               comm;
132315078cd4SMichael Lange   PetscPartitioner       partitioner;
1324f8987ae8SMichael Lange   IS                     cellPart;
1325f8987ae8SMichael Lange   PetscSection           cellPartSection;
1326f8987ae8SMichael Lange   DMLabel                lblPartition, lblMigration;
132715078cd4SMichael Lange   PetscSF                sfProcess, sfMigration;
132870034214SMatthew G. Knepley   PetscBool              flg;
132970034214SMatthew G. Knepley   PetscMPIInt            rank, numProcs, p;
133070034214SMatthew G. Knepley   PetscErrorCode         ierr;
133170034214SMatthew G. Knepley 
133270034214SMatthew G. Knepley   PetscFunctionBegin;
133370034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
133470034214SMatthew G. Knepley   if (sf) PetscValidPointer(sf,4);
133570034214SMatthew G. Knepley   PetscValidPointer(dmParallel,5);
133670034214SMatthew G. Knepley 
133770034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
133870034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
133970034214SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
134070034214SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
134170034214SMatthew G. Knepley 
134270034214SMatthew G. Knepley   *dmParallel = NULL;
134370034214SMatthew G. Knepley   if (numProcs == 1) PetscFunctionReturn(0);
134470034214SMatthew G. Knepley 
134515078cd4SMichael Lange   /* Create cell partition */
134677623264SMatthew G. Knepley   ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
134777623264SMatthew G. Knepley   ierr = PetscSectionCreate(comm, &cellPartSection);CHKERRQ(ierr);
134815078cd4SMichael Lange   ierr = DMPlexGetPartitioner(dm, &partitioner);CHKERRQ(ierr);
134915078cd4SMichael Lange   ierr = PetscPartitionerPartition(partitioner, dm, cellPartSection, &cellPart);CHKERRQ(ierr);
1350f8987ae8SMichael Lange   {
1351f8987ae8SMichael Lange     /* Convert partition to DMLabel */
1352f8987ae8SMichael Lange     PetscInt proc, pStart, pEnd, npoints, poffset;
1353f8987ae8SMichael Lange     const PetscInt *points;
1354f8987ae8SMichael Lange     ierr = DMLabelCreate("Point Partition", &lblPartition);CHKERRQ(ierr);
1355f8987ae8SMichael Lange     ierr = ISGetIndices(cellPart, &points);CHKERRQ(ierr);
1356f8987ae8SMichael Lange     ierr = PetscSectionGetChart(cellPartSection, &pStart, &pEnd);CHKERRQ(ierr);
1357f8987ae8SMichael Lange     for (proc = pStart; proc < pEnd; proc++) {
1358f8987ae8SMichael Lange       ierr = PetscSectionGetDof(cellPartSection, proc, &npoints);CHKERRQ(ierr);
1359f8987ae8SMichael Lange       ierr = PetscSectionGetOffset(cellPartSection, proc, &poffset);CHKERRQ(ierr);
1360f8987ae8SMichael Lange       for (p = poffset; p < poffset+npoints; p++) {
1361f8987ae8SMichael Lange         ierr = DMLabelSetValue(lblPartition, points[p], proc);CHKERRQ(ierr);
136270034214SMatthew G. Knepley       }
1363f8987ae8SMichael Lange     }
1364f8987ae8SMichael Lange     ierr = ISRestoreIndices(cellPart, &points);CHKERRQ(ierr);
1365f8987ae8SMichael Lange   }
1366f8987ae8SMichael Lange   ierr = DMPlexPartitionLabelClosure(dm, lblPartition);CHKERRQ(ierr);
1367f8987ae8SMichael Lange   {
1368f8987ae8SMichael Lange     /* Build a global process SF */
1369f8987ae8SMichael Lange     PetscSFNode *remoteProc;
1370f8987ae8SMichael Lange     ierr = PetscMalloc1(numProcs, &remoteProc);CHKERRQ(ierr);
1371f8987ae8SMichael Lange     for (p = 0; p < numProcs; ++p) {
1372f8987ae8SMichael Lange       remoteProc[p].rank  = p;
1373f8987ae8SMichael Lange       remoteProc[p].index = rank;
1374f8987ae8SMichael Lange     }
1375f8987ae8SMichael Lange     ierr = PetscSFCreate(comm, &sfProcess);CHKERRQ(ierr);
1376f8987ae8SMichael Lange     ierr = PetscObjectSetName((PetscObject) sfProcess, "Process SF");CHKERRQ(ierr);
1377f8987ae8SMichael Lange     ierr = PetscSFSetGraph(sfProcess, numProcs, numProcs, NULL, PETSC_OWN_POINTER, remoteProc, PETSC_OWN_POINTER);CHKERRQ(ierr);
1378f8987ae8SMichael Lange   }
1379f8987ae8SMichael Lange   ierr = DMLabelCreate("Point migration", &lblMigration);CHKERRQ(ierr);
1380f8987ae8SMichael Lange   ierr = DMPlexPartitionLabelInvert(dm, lblPartition, sfProcess, lblMigration);CHKERRQ(ierr);
1381f8987ae8SMichael Lange   ierr = DMPlexPartitionLabelCreateSF(dm, lblMigration, &sfMigration);
1382f8987ae8SMichael Lange   ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
138370034214SMatthew G. Knepley   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
138470034214SMatthew G. Knepley   if (flg) {
1385f8987ae8SMichael Lange     ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr);
1386f8987ae8SMichael Lange     ierr = DMLabelView(lblPartition, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1387f8987ae8SMichael Lange     ierr = PetscSFView(sfMigration, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
138870034214SMatthew G. Knepley   }
1389f8987ae8SMichael Lange 
139015078cd4SMichael Lange   /* Create non-overlapping parallel DM and migrate internal data */
139170034214SMatthew G. Knepley   ierr = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr);
139270034214SMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr);
1393*b9f40539SMichael Lange   ierr = DMPlexMigrate(dm, sfMigration, *dmParallel);CHKERRQ(ierr);
139470034214SMatthew G. Knepley 
1395a157612eSMichael Lange   /* Build the point SF without overlap */
1396f8987ae8SMichael Lange   ierr = DMPlexDistributeSF(dm, sfMigration, *dmParallel);CHKERRQ(ierr);
139770034214SMatthew G. Knepley   if (flg) {
139815fff7beSMatthew G. Knepley     ierr = PetscPrintf(comm, "Point SF:\n");CHKERRQ(ierr);
139915fff7beSMatthew G. Knepley     ierr = PetscSFView((*dmParallel)->sf, NULL);CHKERRQ(ierr);
140070034214SMatthew G. Knepley   }
140170034214SMatthew G. Knepley 
1402a157612eSMichael Lange   if (overlap > 0) {
140315078cd4SMichael Lange     DM                 dmOverlap;
140415078cd4SMichael Lange     PetscInt           nroots, nleaves;
140515078cd4SMichael Lange     PetscSFNode       *newRemote;
140615078cd4SMichael Lange     const PetscSFNode *oldRemote;
140715078cd4SMichael Lange     PetscSF            sfOverlap, sfOverlapPoint;
14083d822a50SMichael Lange     ierr = PetscLogEventBegin(DMPLEX_DistributeOverlap,dm,0,0,0);CHKERRQ(ierr);
1409a157612eSMichael Lange     /* Add the partition overlap to the distributed DM */
1410*b9f40539SMichael Lange     ierr = DMPlexDistributeOverlap(*dmParallel, overlap, &sfOverlap, &dmOverlap);CHKERRQ(ierr);
1411a157612eSMichael Lange     ierr = DMDestroy(dmParallel);CHKERRQ(ierr);
1412a157612eSMichael Lange     *dmParallel = dmOverlap;
1413c389ea9fSToby Isaac     if (flg) {
14143d822a50SMichael Lange       ierr = PetscPrintf(comm, "Overlap Migration SF:\n");CHKERRQ(ierr);
141515078cd4SMichael Lange       ierr = PetscSFView(sfOverlap, NULL);CHKERRQ(ierr);
1416c389ea9fSToby Isaac     }
141743331d4aSMichael Lange 
1418f8987ae8SMichael Lange     /* Re-map the migration SF to establish the full migration pattern */
1419f8987ae8SMichael Lange     ierr = PetscSFGetGraph(sfMigration, &nroots, NULL, NULL, &oldRemote);CHKERRQ(ierr);
142015078cd4SMichael Lange     ierr = PetscSFGetGraph(sfOverlap, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr);
142143331d4aSMichael Lange     ierr = PetscMalloc1(nleaves, &newRemote);CHKERRQ(ierr);
142215078cd4SMichael Lange     ierr = PetscSFBcastBegin(sfOverlap, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr);
142315078cd4SMichael Lange     ierr = PetscSFBcastEnd(sfOverlap, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr);
142415078cd4SMichael Lange     ierr = PetscSFCreate(comm, &sfOverlapPoint);CHKERRQ(ierr);
142515078cd4SMichael Lange     ierr = PetscSFSetGraph(sfOverlapPoint, nroots, nleaves, NULL, PETSC_OWN_POINTER, newRemote, PETSC_OWN_POINTER);CHKERRQ(ierr);
142615078cd4SMichael Lange     ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr);
1427f8987ae8SMichael Lange     ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr);
142815078cd4SMichael Lange     sfMigration = sfOverlapPoint;
14293d822a50SMichael Lange     ierr = PetscLogEventEnd(DMPLEX_DistributeOverlap,dm,0,0,0);CHKERRQ(ierr);
1430c389ea9fSToby Isaac   }
1431bf5244c7SMatthew G. Knepley   /* Cleanup Partition */
1432f8987ae8SMichael Lange   ierr = PetscSFDestroy(&sfProcess);CHKERRQ(ierr);
1433f8987ae8SMichael Lange   ierr = DMLabelDestroy(&lblPartition);CHKERRQ(ierr);
1434f8987ae8SMichael Lange   ierr = DMLabelDestroy(&lblMigration);CHKERRQ(ierr);
14354eca1733SMichael Lange   ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr);
14364eca1733SMichael Lange   ierr = ISDestroy(&cellPart);CHKERRQ(ierr);
1437f8987ae8SMichael Lange   if (sf) {*sf = sfMigration;}
1438f8987ae8SMichael Lange   else    {ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr);}
143970034214SMatthew G. Knepley   ierr = DMSetFromOptions(*dmParallel);CHKERRQ(ierr);
144070034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
144170034214SMatthew G. Knepley   PetscFunctionReturn(0);
144270034214SMatthew G. Knepley }
1443a157612eSMichael Lange 
1444a157612eSMichael Lange #undef __FUNCT__
1445a157612eSMichael Lange #define __FUNCT__ "DMPlexDistributeOverlap"
1446a157612eSMichael Lange /*@C
1447a157612eSMichael Lange   DMPlexDistribute - Add partition overlap to a distributed non-overlapping DM.
1448a157612eSMichael Lange 
1449a157612eSMichael Lange   Not Collective
1450a157612eSMichael Lange 
1451a157612eSMichael Lange   Input Parameter:
1452a157612eSMichael Lange + dm  - The non-overlapping distrbuted DMPlex object
1453a157612eSMichael Lange - overlap - The overlap of partitions, 0 is the default
1454a157612eSMichael Lange 
1455a157612eSMichael Lange   Output Parameter:
1456a157612eSMichael Lange + sf - The PetscSF used for point distribution
1457a157612eSMichael Lange - dmOverlap - The overlapping distributed DMPlex object, or NULL
1458a157612eSMichael Lange 
1459a157612eSMichael Lange   Note: If the mesh was not distributed, the return value is NULL.
1460a157612eSMichael Lange 
1461a157612eSMichael Lange   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
1462a157612eSMichael Lange   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
1463a157612eSMichael Lange   representation on the mesh.
1464a157612eSMichael Lange 
1465a157612eSMichael Lange   Level: intermediate
1466a157612eSMichael Lange 
1467a157612eSMichael Lange .keywords: mesh, elements
1468a157612eSMichael Lange .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
1469a157612eSMichael Lange @*/
1470*b9f40539SMichael Lange PetscErrorCode DMPlexDistributeOverlap(DM dm, PetscInt overlap, PetscSF *sf, DM *dmOverlap)
1471a157612eSMichael Lange {
1472a157612eSMichael Lange   MPI_Comm               comm;
1473a157612eSMichael Lange   PetscMPIInt            rank;
1474a157612eSMichael Lange   PetscSection           rootSection, leafSection;
1475a157612eSMichael Lange   IS                     rootrank, leafrank;
1476a157612eSMichael Lange   PetscSF                overlapSF, migrationSF, pointSF, newPointSF;
1477a157612eSMichael Lange   PetscSFNode           *ghostRemote;
1478a157612eSMichael Lange   const PetscSFNode     *overlapRemote;
147915078cd4SMichael Lange   const PetscInt        *overlapLocal;
148015078cd4SMichael Lange   PetscInt               p, pStart, pEnd, idx;
1481a157612eSMichael Lange   PetscInt               numGhostPoints, numOverlapPoints, numSharedPoints, overlapLeaves;
148215078cd4SMichael Lange   PetscInt              *ghostLocal, *pointIDs, *recvPointIDs;
1483a157612eSMichael Lange   PetscErrorCode         ierr;
1484a157612eSMichael Lange 
1485a157612eSMichael Lange   PetscFunctionBegin;
1486a157612eSMichael Lange   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1487a157612eSMichael Lange   if (sf) PetscValidPointer(sf, 3);
1488a157612eSMichael Lange   PetscValidPointer(dmOverlap, 4);
1489a157612eSMichael Lange 
1490a157612eSMichael Lange   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1491a157612eSMichael Lange   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1492a157612eSMichael Lange 
1493a157612eSMichael Lange   /* Compute point overlap with neighbouring processes on the distributed DM */
1494a157612eSMichael Lange   ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1495a157612eSMichael Lange   ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr);
1496a157612eSMichael Lange   ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr);
1497a157612eSMichael Lange   ierr = DMPlexDistributeOwnership(dm, rootSection, &rootrank, leafSection, &leafrank);CHKERRQ(ierr);
149824d039d7SMichael Lange   ierr = DMPlexCreateOverlap(dm, overlap, rootSection, rootrank, leafSection, leafrank, &overlapSF);CHKERRQ(ierr);
149915fff7beSMatthew G. Knepley   ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr);
150015fff7beSMatthew G. Knepley   ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr);
150115fff7beSMatthew G. Knepley   ierr = ISDestroy(&rootrank);CHKERRQ(ierr);
150215fff7beSMatthew G. Knepley   ierr = ISDestroy(&leafrank);CHKERRQ(ierr);
1503a157612eSMichael Lange   ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr);
1504a157612eSMichael Lange 
1505a157612eSMichael Lange   /* Build dense migration SF that maps the non-overlapping partition to the overlapping one */
1506a157612eSMichael Lange   ierr = DMPlexCreateOverlapMigrationSF(dm, overlapSF, &migrationSF);CHKERRQ(ierr);
1507a157612eSMichael Lange 
1508a157612eSMichael Lange   /* Build the overlapping DM */
1509a157612eSMichael Lange   ierr = DMPlexCreate(comm, dmOverlap);CHKERRQ(ierr);
1510a157612eSMichael Lange   ierr = PetscObjectSetName((PetscObject) *dmOverlap, "Parallel Mesh");CHKERRQ(ierr);
1511*b9f40539SMichael Lange   ierr = DMPlexMigrate(dm, migrationSF, *dmOverlap);CHKERRQ(ierr);
1512a157612eSMichael Lange 
1513a157612eSMichael Lange   /* Build the new point SF by propagating the depthShift generate remote root indices */
1514a157612eSMichael Lange   ierr = DMGetPointSF(dm, &pointSF);CHKERRQ(ierr);
1515a157612eSMichael Lange   ierr = PetscSFGetGraph(pointSF, NULL, &numSharedPoints, NULL, NULL);CHKERRQ(ierr);
1516a157612eSMichael Lange   ierr = PetscSFGetGraph(overlapSF, NULL, &numOverlapPoints, NULL, NULL);CHKERRQ(ierr);
1517a157612eSMichael Lange   numGhostPoints = numSharedPoints + numOverlapPoints;
151809b7985cSMatthew G. Knepley   ierr = PetscMalloc1(numGhostPoints, &ghostLocal);CHKERRQ(ierr);
151909b7985cSMatthew G. Knepley   ierr = PetscMalloc1(numGhostPoints, &ghostRemote);CHKERRQ(ierr);
1520a157612eSMichael Lange   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
152115078cd4SMichael Lange   ierr = PetscSFGetGraph(migrationSF, NULL, &overlapLeaves, &overlapLocal, &overlapRemote);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);
154015fff7beSMatthew G. Knepley   ierr = PetscSFDestroy(&newPointSF);CHKERRQ(ierr);
1541a157612eSMichael Lange   /* Cleanup overlap partition */
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   PetscFunctionReturn(0);
1547a157612eSMichael Lange }
1548