xref: /petsc/src/dm/impls/plex/plexdistribute.c (revision 1e8fc25d9fb60c227919fdb72934256ec06026ef)
170034214SMatthew G. Knepley #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
270034214SMatthew G. Knepley #include <petscsf.h>
370034214SMatthew G. Knepley 
470034214SMatthew G. Knepley #undef __FUNCT__
570034214SMatthew G. Knepley #define __FUNCT__ "DMPlexSetAdjacencyUseCone"
670034214SMatthew G. Knepley /*@
770034214SMatthew G. Knepley   DMPlexSetAdjacencyUseCone - Define adjacency in the mesh using either the cone or the support first
870034214SMatthew G. Knepley 
970034214SMatthew G. Knepley   Input Parameters:
1070034214SMatthew G. Knepley + dm      - The DM object
1170034214SMatthew G. Knepley - useCone - Flag to use the cone first
1270034214SMatthew G. Knepley 
1370034214SMatthew G. Knepley   Level: intermediate
1470034214SMatthew G. Knepley 
1570034214SMatthew G. Knepley   Notes:
1670034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE
1770034214SMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in star(cone(p)),    useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
1870034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
1970034214SMatthew G. Knepley 
2070034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
2170034214SMatthew G. Knepley @*/
2270034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseCone(DM dm, PetscBool useCone)
2370034214SMatthew G. Knepley {
2470034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
2570034214SMatthew G. Knepley 
2670034214SMatthew G. Knepley   PetscFunctionBegin;
2770034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
2870034214SMatthew G. Knepley   mesh->useCone = useCone;
2970034214SMatthew G. Knepley   PetscFunctionReturn(0);
3070034214SMatthew G. Knepley }
3170034214SMatthew G. Knepley 
3270034214SMatthew G. Knepley #undef __FUNCT__
3370034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacencyUseCone"
3470034214SMatthew G. Knepley /*@
3570034214SMatthew G. Knepley   DMPlexGetAdjacencyUseCone - Query whether adjacency in the mesh uses the cone or the support first
3670034214SMatthew G. Knepley 
3770034214SMatthew G. Knepley   Input Parameter:
3870034214SMatthew G. Knepley . dm      - The DM object
3970034214SMatthew G. Knepley 
4070034214SMatthew G. Knepley   Output Parameter:
4170034214SMatthew G. Knepley . useCone - Flag to use the cone first
4270034214SMatthew G. Knepley 
4370034214SMatthew G. Knepley   Level: intermediate
4470034214SMatthew G. Knepley 
4570034214SMatthew G. Knepley   Notes:
4670034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE
4770034214SMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in star(cone(p)),    useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
4870034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
4970034214SMatthew G. Knepley 
5070034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator()
5170034214SMatthew G. Knepley @*/
5270034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseCone(DM dm, PetscBool *useCone)
5370034214SMatthew G. Knepley {
5470034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
5570034214SMatthew G. Knepley 
5670034214SMatthew G. Knepley   PetscFunctionBegin;
5770034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
5870034214SMatthew G. Knepley   PetscValidIntPointer(useCone, 2);
5970034214SMatthew G. Knepley   *useCone = mesh->useCone;
6070034214SMatthew G. Knepley   PetscFunctionReturn(0);
6170034214SMatthew G. Knepley }
6270034214SMatthew G. Knepley 
6370034214SMatthew G. Knepley #undef __FUNCT__
6470034214SMatthew G. Knepley #define __FUNCT__ "DMPlexSetAdjacencyUseClosure"
6570034214SMatthew G. Knepley /*@
6670034214SMatthew G. Knepley   DMPlexSetAdjacencyUseClosure - Define adjacency in the mesh using the transitive closure
6770034214SMatthew G. Knepley 
6870034214SMatthew G. Knepley   Input Parameters:
6970034214SMatthew G. Knepley + dm      - The DM object
7070034214SMatthew G. Knepley - useClosure - Flag to use the closure
7170034214SMatthew G. Knepley 
7270034214SMatthew G. Knepley   Level: intermediate
7370034214SMatthew G. Knepley 
7470034214SMatthew G. Knepley   Notes:
7570034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE
7670034214SMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in star(cone(p)),    useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
7770034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
7870034214SMatthew G. Knepley 
7970034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
8070034214SMatthew G. Knepley @*/
8170034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseClosure(DM dm, PetscBool useClosure)
8270034214SMatthew G. Knepley {
8370034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
8470034214SMatthew G. Knepley 
8570034214SMatthew G. Knepley   PetscFunctionBegin;
8670034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
8770034214SMatthew G. Knepley   mesh->useClosure = useClosure;
8870034214SMatthew G. Knepley   PetscFunctionReturn(0);
8970034214SMatthew G. Knepley }
9070034214SMatthew G. Knepley 
9170034214SMatthew G. Knepley #undef __FUNCT__
9270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacencyUseClosure"
9370034214SMatthew G. Knepley /*@
9470034214SMatthew G. Knepley   DMPlexGetAdjacencyUseClosure - Query whether adjacency in the mesh uses the transitive closure
9570034214SMatthew G. Knepley 
9670034214SMatthew G. Knepley   Input Parameter:
9770034214SMatthew G. Knepley . dm      - The DM object
9870034214SMatthew G. Knepley 
9970034214SMatthew G. Knepley   Output Parameter:
10070034214SMatthew G. Knepley . useClosure - Flag to use the closure
10170034214SMatthew G. Knepley 
10270034214SMatthew G. Knepley   Level: intermediate
10370034214SMatthew G. Knepley 
10470034214SMatthew G. Knepley   Notes:
10570034214SMatthew G. Knepley $     FEM:   Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE
10670034214SMatthew G. Knepley $     FVM:   Two points p and q are adjacent if q \in star(cone(p)),    useCone = PETSC_TRUE,  useClosure = PETSC_FALSE
10770034214SMatthew G. Knepley $     FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE,  useClosure = PETSC_TRUE
10870034214SMatthew G. Knepley 
10970034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator()
11070034214SMatthew G. Knepley @*/
11170034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseClosure(DM dm, PetscBool *useClosure)
11270034214SMatthew G. Knepley {
11370034214SMatthew G. Knepley   DM_Plex *mesh = (DM_Plex *) dm->data;
11470034214SMatthew G. Knepley 
11570034214SMatthew G. Knepley   PetscFunctionBegin;
11670034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
11770034214SMatthew G. Knepley   PetscValidIntPointer(useClosure, 2);
11870034214SMatthew G. Knepley   *useClosure = mesh->useClosure;
11970034214SMatthew G. Knepley   PetscFunctionReturn(0);
12070034214SMatthew G. Knepley }
12170034214SMatthew G. Knepley 
12270034214SMatthew G. Knepley #undef __FUNCT__
12370034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Cone_Internal"
12470034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Cone_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
12570034214SMatthew G. Knepley {
12670034214SMatthew G. Knepley   const PetscInt *cone = NULL;
12770034214SMatthew G. Knepley   PetscInt        numAdj = 0, maxAdjSize = *adjSize, coneSize, c;
12870034214SMatthew G. Knepley   PetscErrorCode  ierr;
12970034214SMatthew G. Knepley 
13070034214SMatthew G. Knepley   PetscFunctionBeginHot;
13170034214SMatthew G. Knepley   ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr);
13270034214SMatthew G. Knepley   ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr);
13370034214SMatthew G. Knepley   for (c = 0; c < coneSize; ++c) {
13470034214SMatthew G. Knepley     const PetscInt *support = NULL;
13570034214SMatthew G. Knepley     PetscInt        supportSize, s, q;
13670034214SMatthew G. Knepley 
13770034214SMatthew G. Knepley     ierr = DMPlexGetSupportSize(dm, cone[c], &supportSize);CHKERRQ(ierr);
13870034214SMatthew G. Knepley     ierr = DMPlexGetSupport(dm, cone[c], &support);CHKERRQ(ierr);
13970034214SMatthew G. Knepley     for (s = 0; s < supportSize; ++s) {
14070034214SMatthew G. Knepley       for (q = 0; q < numAdj || (adj[numAdj++] = support[s],0); ++q) {
14170034214SMatthew G. Knepley         if (support[s] == adj[q]) break;
14270034214SMatthew G. Knepley       }
14370034214SMatthew G. Knepley       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
14470034214SMatthew G. Knepley     }
14570034214SMatthew G. Knepley   }
14670034214SMatthew G. Knepley   *adjSize = numAdj;
14770034214SMatthew G. Knepley   PetscFunctionReturn(0);
14870034214SMatthew G. Knepley }
14970034214SMatthew G. Knepley 
15070034214SMatthew G. Knepley #undef __FUNCT__
15170034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Support_Internal"
15270034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Support_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[])
15370034214SMatthew G. Knepley {
15470034214SMatthew G. Knepley   const PetscInt *support = NULL;
15570034214SMatthew G. Knepley   PetscInt        numAdj   = 0, maxAdjSize = *adjSize, supportSize, s;
15670034214SMatthew G. Knepley   PetscErrorCode  ierr;
15770034214SMatthew G. Knepley 
15870034214SMatthew G. Knepley   PetscFunctionBeginHot;
15970034214SMatthew G. Knepley   ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr);
16070034214SMatthew G. Knepley   ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr);
16170034214SMatthew G. Knepley   for (s = 0; s < supportSize; ++s) {
16270034214SMatthew G. Knepley     const PetscInt *cone = NULL;
16370034214SMatthew G. Knepley     PetscInt        coneSize, c, q;
16470034214SMatthew G. Knepley 
16570034214SMatthew G. Knepley     ierr = DMPlexGetConeSize(dm, support[s], &coneSize);CHKERRQ(ierr);
16670034214SMatthew G. Knepley     ierr = DMPlexGetCone(dm, support[s], &cone);CHKERRQ(ierr);
16770034214SMatthew G. Knepley     for (c = 0; c < coneSize; ++c) {
16870034214SMatthew G. Knepley       for (q = 0; q < numAdj || (adj[numAdj++] = cone[c],0); ++q) {
16970034214SMatthew G. Knepley         if (cone[c] == adj[q]) break;
17070034214SMatthew G. Knepley       }
17170034214SMatthew G. Knepley       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
17270034214SMatthew G. Knepley     }
17370034214SMatthew G. Knepley   }
17470034214SMatthew G. Knepley   *adjSize = numAdj;
17570034214SMatthew G. Knepley   PetscFunctionReturn(0);
17670034214SMatthew G. Knepley }
17770034214SMatthew G. Knepley 
17870034214SMatthew G. Knepley #undef __FUNCT__
17970034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Transitive_Internal"
18070034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Transitive_Internal(DM dm, PetscInt p, PetscBool useClosure, PetscInt *adjSize, PetscInt adj[])
18170034214SMatthew G. Knepley {
18270034214SMatthew G. Knepley   PetscInt      *star = NULL;
18370034214SMatthew G. Knepley   PetscInt       numAdj = 0, maxAdjSize = *adjSize, starSize, s;
18470034214SMatthew G. Knepley   PetscErrorCode ierr;
18570034214SMatthew G. Knepley 
18670034214SMatthew G. Knepley   PetscFunctionBeginHot;
18770034214SMatthew G. Knepley   ierr = DMPlexGetTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr);
18870034214SMatthew G. Knepley   for (s = 0; s < starSize*2; s += 2) {
18970034214SMatthew G. Knepley     const PetscInt *closure = NULL;
19070034214SMatthew G. Knepley     PetscInt        closureSize, c, q;
19170034214SMatthew G. Knepley 
19270034214SMatthew G. Knepley     ierr = DMPlexGetTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
19370034214SMatthew G. Knepley     for (c = 0; c < closureSize*2; c += 2) {
19470034214SMatthew G. Knepley       for (q = 0; q < numAdj || (adj[numAdj++] = closure[c],0); ++q) {
19570034214SMatthew G. Knepley         if (closure[c] == adj[q]) break;
19670034214SMatthew G. Knepley       }
19770034214SMatthew G. Knepley       if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize);
19870034214SMatthew G. Knepley     }
19970034214SMatthew G. Knepley     ierr = DMPlexRestoreTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr);
20070034214SMatthew G. Knepley   }
20170034214SMatthew G. Knepley   ierr = DMPlexRestoreTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr);
20270034214SMatthew G. Knepley   *adjSize = numAdj;
20370034214SMatthew G. Knepley   PetscFunctionReturn(0);
20470034214SMatthew G. Knepley }
20570034214SMatthew G. Knepley 
20670034214SMatthew G. Knepley #undef __FUNCT__
20770034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Internal"
20879bad331SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacency_Internal(DM dm, PetscInt p, PetscBool useCone, PetscBool useTransitiveClosure, PetscInt *adjSize, PetscInt *adj[])
20970034214SMatthew G. Knepley {
21079bad331SMatthew G. Knepley   static PetscInt asiz = 0;
21170034214SMatthew G. Knepley   PetscErrorCode  ierr;
21270034214SMatthew G. Knepley 
21370034214SMatthew G. Knepley   PetscFunctionBeginHot;
21479bad331SMatthew G. Knepley   if (!*adj) {
21579bad331SMatthew G. Knepley     PetscInt depth, maxConeSize, maxSupportSize;
21679bad331SMatthew G. Knepley 
21779bad331SMatthew G. Knepley     ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
21879bad331SMatthew G. Knepley     ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr);
21979bad331SMatthew G. Knepley     asiz = PetscPowInt(maxConeSize, depth+1) * PetscPowInt(maxSupportSize, depth+1) + 1;
22079bad331SMatthew G. Knepley     ierr = PetscMalloc1(asiz,adj);CHKERRQ(ierr);
22179bad331SMatthew G. Knepley   }
22279bad331SMatthew G. Knepley   if (*adjSize < 0) *adjSize = asiz;
22370034214SMatthew G. Knepley   if (useTransitiveClosure) {
22479bad331SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Transitive_Internal(dm, p, useCone, adjSize, *adj);CHKERRQ(ierr);
22570034214SMatthew G. Knepley   } else if (useCone) {
22679bad331SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Cone_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
22770034214SMatthew G. Knepley   } else {
22879bad331SMatthew G. Knepley     ierr = DMPlexGetAdjacency_Support_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr);
22970034214SMatthew G. Knepley   }
23070034214SMatthew G. Knepley   PetscFunctionReturn(0);
23170034214SMatthew G. Knepley }
23270034214SMatthew G. Knepley 
23370034214SMatthew G. Knepley #undef __FUNCT__
23470034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency"
23570034214SMatthew G. Knepley /*@
23670034214SMatthew G. Knepley   DMPlexGetAdjacency - Return all points adjacent to the given point
23770034214SMatthew G. Knepley 
23870034214SMatthew G. Knepley   Input Parameters:
23970034214SMatthew G. Knepley + dm - The DM object
24070034214SMatthew G. Knepley . p  - The point
24170034214SMatthew G. Knepley . adjSize - The maximum size of adj if it is non-NULL, or PETSC_DETERMINE
24270034214SMatthew G. Knepley - adj - Either NULL so that the array is allocated, or an existing array with size adjSize
24370034214SMatthew G. Knepley 
24470034214SMatthew G. Knepley   Output Parameters:
24570034214SMatthew G. Knepley + adjSize - The number of adjacent points
24670034214SMatthew G. Knepley - adj - The adjacent points
24770034214SMatthew G. Knepley 
24870034214SMatthew G. Knepley   Level: advanced
24970034214SMatthew G. Knepley 
25070034214SMatthew G. Knepley   Notes: The user must PetscFree the adj array if it was not passed in.
25170034214SMatthew G. Knepley 
25270034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexDistribute(), DMCreateMatrix(), DMPlexPreallocateOperator()
25370034214SMatthew G. Knepley @*/
25470034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacency(DM dm, PetscInt p, PetscInt *adjSize, PetscInt *adj[])
25570034214SMatthew G. Knepley {
25670034214SMatthew G. Knepley   DM_Plex       *mesh = (DM_Plex *) dm->data;
25770034214SMatthew G. Knepley   PetscErrorCode ierr;
25870034214SMatthew G. Knepley 
25970034214SMatthew G. Knepley   PetscFunctionBeginHot;
26070034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
26170034214SMatthew G. Knepley   PetscValidPointer(adjSize,3);
26270034214SMatthew G. Knepley   PetscValidPointer(adj,4);
26379bad331SMatthew G. Knepley   ierr = DMPlexGetAdjacency_Internal(dm, p, mesh->useCone, mesh->useClosure, adjSize, adj);CHKERRQ(ierr);
26470034214SMatthew G. Knepley   PetscFunctionReturn(0);
26570034214SMatthew G. Knepley }
26670034214SMatthew G. Knepley 
26770034214SMatthew G. Knepley #undef __FUNCT__
26870034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeField"
26970034214SMatthew G. Knepley /*@
27070034214SMatthew G. Knepley   DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
27170034214SMatthew G. Knepley 
27270034214SMatthew G. Knepley   Collective on DM
27370034214SMatthew G. Knepley 
27470034214SMatthew G. Knepley   Input Parameters:
27570034214SMatthew G. Knepley + dm - The DMPlex object
27670034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
27770034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
27870034214SMatthew G. Knepley - originalVec - The existing data
27970034214SMatthew G. Knepley 
28070034214SMatthew G. Knepley   Output Parameters:
28170034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout
28270034214SMatthew G. Knepley - newVec - The new data
28370034214SMatthew G. Knepley 
28470034214SMatthew G. Knepley   Level: developer
28570034214SMatthew G. Knepley 
286*1e8fc25dSMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeFieldIS(), DMPlexDistributeData()
28770034214SMatthew G. Knepley @*/
28870034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec)
28970034214SMatthew G. Knepley {
29070034214SMatthew G. Knepley   PetscSF        fieldSF;
29170034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
29270034214SMatthew G. Knepley   PetscScalar   *originalValues, *newValues;
29370034214SMatthew G. Knepley   PetscErrorCode ierr;
29470034214SMatthew G. Knepley 
29570034214SMatthew G. Knepley   PetscFunctionBegin;
29670034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
29770034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
29870034214SMatthew G. Knepley 
29970034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
30070034214SMatthew G. Knepley   ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr);
30170034214SMatthew G. Knepley   ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr);
30270034214SMatthew G. Knepley 
30370034214SMatthew G. Knepley   ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr);
30470034214SMatthew G. Knepley   ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr);
30570034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
30670034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
30770034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr);
30870034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
30970034214SMatthew G. Knepley   ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr);
31070034214SMatthew G. Knepley   ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr);
31170034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
31270034214SMatthew G. Knepley   PetscFunctionReturn(0);
31370034214SMatthew G. Knepley }
31470034214SMatthew G. Knepley 
31570034214SMatthew G. Knepley #undef __FUNCT__
316*1e8fc25dSMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeFieldIS"
317*1e8fc25dSMatthew G. Knepley /*@
318*1e8fc25dSMatthew G. Knepley   DMPlexDistributeFieldIS - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
319*1e8fc25dSMatthew G. Knepley 
320*1e8fc25dSMatthew G. Knepley   Collective on DM
321*1e8fc25dSMatthew G. Knepley 
322*1e8fc25dSMatthew G. Knepley   Input Parameters:
323*1e8fc25dSMatthew G. Knepley + dm - The DMPlex object
324*1e8fc25dSMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
325*1e8fc25dSMatthew G. Knepley . originalSection - The PetscSection for existing data layout
326*1e8fc25dSMatthew G. Knepley - originalIS - The existing data
327*1e8fc25dSMatthew G. Knepley 
328*1e8fc25dSMatthew G. Knepley   Output Parameters:
329*1e8fc25dSMatthew G. Knepley + newSection - The PetscSF describing the new data layout
330*1e8fc25dSMatthew G. Knepley - newIS - The new data
331*1e8fc25dSMatthew G. Knepley 
332*1e8fc25dSMatthew G. Knepley   Level: developer
333*1e8fc25dSMatthew G. Knepley 
334*1e8fc25dSMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexDistributeData()
335*1e8fc25dSMatthew G. Knepley @*/
336*1e8fc25dSMatthew G. Knepley PetscErrorCode DMPlexDistributeFieldIS(DM dm, PetscSF pointSF, PetscSection originalSection, IS originalIS, PetscSection newSection, IS *newIS)
337*1e8fc25dSMatthew G. Knepley {
338*1e8fc25dSMatthew G. Knepley   PetscSF         fieldSF;
339*1e8fc25dSMatthew G. Knepley   PetscInt       *newValues, *remoteOffsets, fieldSize;
340*1e8fc25dSMatthew G. Knepley   const PetscInt *originalValues;
341*1e8fc25dSMatthew G. Knepley   PetscErrorCode  ierr;
342*1e8fc25dSMatthew G. Knepley 
343*1e8fc25dSMatthew G. Knepley   PetscFunctionBegin;
344*1e8fc25dSMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
345*1e8fc25dSMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
346*1e8fc25dSMatthew G. Knepley 
347*1e8fc25dSMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
348*1e8fc25dSMatthew G. Knepley   ierr = PetscMalloc(fieldSize * sizeof(PetscInt), &newValues);CHKERRQ(ierr);
349*1e8fc25dSMatthew G. Knepley 
350*1e8fc25dSMatthew G. Knepley   ierr = ISGetIndices(originalIS, &originalValues);CHKERRQ(ierr);
351*1e8fc25dSMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
352*1e8fc25dSMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, MPIU_INT, originalValues, newValues);CHKERRQ(ierr);
353*1e8fc25dSMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, MPIU_INT, originalValues, newValues);CHKERRQ(ierr);
354*1e8fc25dSMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
355*1e8fc25dSMatthew G. Knepley   ierr = ISRestoreIndices(originalIS, &originalValues);CHKERRQ(ierr);
356*1e8fc25dSMatthew G. Knepley   ierr = ISCreateGeneral(PetscObjectComm((PetscObject) pointSF), fieldSize, newValues, PETSC_OWN_POINTER, newIS);CHKERRQ(ierr);
357*1e8fc25dSMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr);
358*1e8fc25dSMatthew G. Knepley   PetscFunctionReturn(0);
359*1e8fc25dSMatthew G. Knepley }
360*1e8fc25dSMatthew G. Knepley 
361*1e8fc25dSMatthew G. Knepley #undef __FUNCT__
36270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeData"
36370034214SMatthew G. Knepley /*@
36470034214SMatthew G. Knepley   DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution
36570034214SMatthew G. Knepley 
36670034214SMatthew G. Knepley   Collective on DM
36770034214SMatthew G. Knepley 
36870034214SMatthew G. Knepley   Input Parameters:
36970034214SMatthew G. Knepley + dm - The DMPlex object
37070034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern
37170034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout
37270034214SMatthew G. Knepley . datatype - The type of data
37370034214SMatthew G. Knepley - originalData - The existing data
37470034214SMatthew G. Knepley 
37570034214SMatthew G. Knepley   Output Parameters:
376*1e8fc25dSMatthew G. Knepley + newSection - The PetscSection describing the new data layout
37770034214SMatthew G. Knepley - newData - The new data
37870034214SMatthew G. Knepley 
37970034214SMatthew G. Knepley   Level: developer
38070034214SMatthew G. Knepley 
38170034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField()
38270034214SMatthew G. Knepley @*/
38370034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData)
38470034214SMatthew G. Knepley {
38570034214SMatthew G. Knepley   PetscSF        fieldSF;
38670034214SMatthew G. Knepley   PetscInt      *remoteOffsets, fieldSize;
38770034214SMatthew G. Knepley   PetscMPIInt    dataSize;
38870034214SMatthew G. Knepley   PetscErrorCode ierr;
38970034214SMatthew G. Knepley 
39070034214SMatthew G. Knepley   PetscFunctionBegin;
39170034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
39270034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr);
39370034214SMatthew G. Knepley 
39470034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr);
39570034214SMatthew G. Knepley   ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr);
39670034214SMatthew G. Knepley   ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr);
39770034214SMatthew G. Knepley 
39870034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr);
39970034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
40070034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr);
40170034214SMatthew G. Knepley   ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr);
40270034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr);
40370034214SMatthew G. Knepley   PetscFunctionReturn(0);
40470034214SMatthew G. Knepley }
40570034214SMatthew G. Knepley 
40670034214SMatthew G. Knepley #undef __FUNCT__
40770034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistribute"
40870034214SMatthew G. Knepley /*@C
40970034214SMatthew G. Knepley   DMPlexDistribute - Distributes the mesh and any associated sections.
41070034214SMatthew G. Knepley 
41170034214SMatthew G. Knepley   Not Collective
41270034214SMatthew G. Knepley 
41370034214SMatthew G. Knepley   Input Parameter:
41470034214SMatthew G. Knepley + dm  - The original DMPlex object
41570034214SMatthew G. Knepley . partitioner - The partitioning package, or NULL for the default
41670034214SMatthew G. Knepley - overlap - The overlap of partitions, 0 is the default
41770034214SMatthew G. Knepley 
41870034214SMatthew G. Knepley   Output Parameter:
41970034214SMatthew G. Knepley + sf - The PetscSF used for point distribution
42070034214SMatthew G. Knepley - parallelMesh - The distributed DMPlex object, or NULL
42170034214SMatthew G. Knepley 
422a9c22940SMatthew G. Knepley   Note: If the mesh was not distributed, the return value is NULL.
423a9c22940SMatthew G. Knepley 
424a9c22940SMatthew G. Knepley   The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and
425a9c22940SMatthew G. Knepley   DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function
426a9c22940SMatthew G. Knepley   representation on the mesh.
42770034214SMatthew G. Knepley 
42870034214SMatthew G. Knepley   Level: intermediate
42970034214SMatthew G. Knepley 
43070034214SMatthew G. Knepley .keywords: mesh, elements
431a9c22940SMatthew G. Knepley .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure()
43270034214SMatthew G. Knepley @*/
43370034214SMatthew G. Knepley PetscErrorCode DMPlexDistribute(DM dm, const char partitioner[], PetscInt overlap, PetscSF *sf, DM *dmParallel)
43470034214SMatthew G. Knepley {
43570034214SMatthew G. Knepley   DM_Plex               *mesh   = (DM_Plex*) dm->data, *pmesh;
43670034214SMatthew G. Knepley   MPI_Comm               comm;
43770034214SMatthew G. Knepley   const PetscInt         height = 0;
43870034214SMatthew G. Knepley   PetscInt               dim, numRemoteRanks;
43970034214SMatthew G. Knepley   IS                     origCellPart,        origPart,        cellPart,        part;
44070034214SMatthew G. Knepley   PetscSection           origCellPartSection, origPartSection, cellPartSection, partSection;
44170034214SMatthew G. Knepley   PetscSFNode           *remoteRanks;
44270034214SMatthew G. Knepley   PetscSF                partSF, pointSF, coneSF;
44370034214SMatthew G. Knepley   ISLocalToGlobalMapping renumbering;
44470034214SMatthew G. Knepley   PetscSection           originalConeSection, newConeSection;
44570034214SMatthew G. Knepley   PetscInt              *remoteOffsets;
44670034214SMatthew G. Knepley   PetscInt              *cones, *newCones, newConesSize;
44770034214SMatthew G. Knepley   PetscBool              flg;
44870034214SMatthew G. Knepley   PetscMPIInt            rank, numProcs, p;
44970034214SMatthew G. Knepley   PetscErrorCode         ierr;
45070034214SMatthew G. Knepley 
45170034214SMatthew G. Knepley   PetscFunctionBegin;
45270034214SMatthew G. Knepley   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
45370034214SMatthew G. Knepley   if (sf) PetscValidPointer(sf,4);
45470034214SMatthew G. Knepley   PetscValidPointer(dmParallel,5);
45570034214SMatthew G. Knepley 
45670034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
45770034214SMatthew G. Knepley   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
45870034214SMatthew G. Knepley   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
45970034214SMatthew G. Knepley   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
46070034214SMatthew G. Knepley 
46170034214SMatthew G. Knepley   *dmParallel = NULL;
46270034214SMatthew G. Knepley   if (numProcs == 1) PetscFunctionReturn(0);
46370034214SMatthew G. Knepley 
464c73cfb54SMatthew G. Knepley   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
46570034214SMatthew G. Knepley   /* Create cell partition - We need to rewrite to use IS, use the MatPartition stuff */
46670034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_Partition,dm,0,0,0);CHKERRQ(ierr);
46770034214SMatthew G. Knepley   if (overlap > 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Overlap > 1 not yet implemented");
46870034214SMatthew G. Knepley   ierr = DMPlexCreatePartition(dm, partitioner, height, overlap > 0 ? PETSC_TRUE : PETSC_FALSE, &cellPartSection, &cellPart, &origCellPartSection, &origCellPart);CHKERRQ(ierr);
46970034214SMatthew G. Knepley   /* Create SF assuming a serial partition for all processes: Could check for IS length here */
47070034214SMatthew G. Knepley   if (!rank) numRemoteRanks = numProcs;
47170034214SMatthew G. Knepley   else       numRemoteRanks = 0;
47270034214SMatthew G. Knepley   ierr = PetscMalloc1(numRemoteRanks, &remoteRanks);CHKERRQ(ierr);
47370034214SMatthew G. Knepley   for (p = 0; p < numRemoteRanks; ++p) {
47470034214SMatthew G. Knepley     remoteRanks[p].rank  = p;
47570034214SMatthew G. Knepley     remoteRanks[p].index = 0;
47670034214SMatthew G. Knepley   }
47770034214SMatthew G. Knepley   ierr = PetscSFCreate(comm, &partSF);CHKERRQ(ierr);
47870034214SMatthew G. Knepley   ierr = PetscSFSetGraph(partSF, 1, numRemoteRanks, NULL, PETSC_OWN_POINTER, remoteRanks, PETSC_OWN_POINTER);CHKERRQ(ierr);
47970034214SMatthew G. Knepley   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr);
48070034214SMatthew G. Knepley   if (flg) {
48170034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Cell Partition:\n");CHKERRQ(ierr);
48270034214SMatthew G. Knepley     ierr = PetscSectionView(cellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
48370034214SMatthew G. Knepley     ierr = ISView(cellPart, NULL);CHKERRQ(ierr);
48470034214SMatthew G. Knepley     if (origCellPart) {
48570034214SMatthew G. Knepley       ierr = PetscPrintf(comm, "Original Cell Partition:\n");CHKERRQ(ierr);
48670034214SMatthew G. Knepley       ierr = PetscSectionView(origCellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
48770034214SMatthew G. Knepley       ierr = ISView(origCellPart, NULL);CHKERRQ(ierr);
48870034214SMatthew G. Knepley     }
48970034214SMatthew G. Knepley     ierr = PetscSFView(partSF, NULL);CHKERRQ(ierr);
49070034214SMatthew G. Knepley   }
49170034214SMatthew G. Knepley   /* Close the partition over the mesh */
49270034214SMatthew G. Knepley   ierr = DMPlexCreatePartitionClosure(dm, cellPartSection, cellPart, &partSection, &part);CHKERRQ(ierr);
49370034214SMatthew G. Knepley   ierr = ISDestroy(&cellPart);CHKERRQ(ierr);
49470034214SMatthew G. Knepley   ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr);
49570034214SMatthew G. Knepley   /* Create new mesh */
49670034214SMatthew G. Knepley   ierr  = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr);
497c73cfb54SMatthew G. Knepley   ierr  = DMSetDimension(*dmParallel, dim);CHKERRQ(ierr);
49870034214SMatthew G. Knepley   ierr  = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr);
49970034214SMatthew G. Knepley   pmesh = (DM_Plex*) (*dmParallel)->data;
50070034214SMatthew G. Knepley   /* Distribute sieve points and the global point numbering (replaces creating remote bases) */
50170034214SMatthew G. Knepley   ierr = PetscSFConvertPartition(partSF, partSection, part, &renumbering, &pointSF);CHKERRQ(ierr);
50270034214SMatthew G. Knepley   if (flg) {
50370034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Point Partition:\n");CHKERRQ(ierr);
50470034214SMatthew G. Knepley     ierr = PetscSectionView(partSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
50570034214SMatthew G. Knepley     ierr = ISView(part, NULL);CHKERRQ(ierr);
50670034214SMatthew G. Knepley     ierr = PetscSFView(pointSF, NULL);CHKERRQ(ierr);
50770034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Point Renumbering after partition:\n");CHKERRQ(ierr);
50870034214SMatthew G. Knepley     ierr = ISLocalToGlobalMappingView(renumbering, NULL);CHKERRQ(ierr);
50970034214SMatthew G. Knepley   }
51070034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Partition,dm,0,0,0);CHKERRQ(ierr);
51170034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
51270034214SMatthew G. Knepley   /* Distribute cone section */
51370034214SMatthew G. Knepley   ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr);
51470034214SMatthew G. Knepley   ierr = DMPlexGetConeSection(*dmParallel, &newConeSection);CHKERRQ(ierr);
51570034214SMatthew G. Knepley   ierr = PetscSFDistributeSection(pointSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr);
51670034214SMatthew G. Knepley   ierr = DMSetUp(*dmParallel);CHKERRQ(ierr);
51770034214SMatthew G. Knepley   {
51870034214SMatthew G. Knepley     PetscInt pStart, pEnd, p;
51970034214SMatthew G. Knepley 
52070034214SMatthew G. Knepley     ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr);
52170034214SMatthew G. Knepley     for (p = pStart; p < pEnd; ++p) {
52270034214SMatthew G. Knepley       PetscInt coneSize;
52370034214SMatthew G. Knepley       ierr               = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr);
52470034214SMatthew G. Knepley       pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize);
52570034214SMatthew G. Knepley     }
52670034214SMatthew G. Knepley   }
52770034214SMatthew G. Knepley   /* Communicate and renumber cones */
52870034214SMatthew G. Knepley   ierr = PetscSFCreateSectionSF(pointSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr);
52970034214SMatthew G. Knepley   ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
53070034214SMatthew G. Knepley   ierr = DMPlexGetCones(*dmParallel, &newCones);CHKERRQ(ierr);
53170034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
53270034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
53370034214SMatthew G. Knepley   ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr);
5349d90f715SBarry Smith   ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr);
53570034214SMatthew G. Knepley   ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr);
53670034214SMatthew G. Knepley   if (flg) {
53770034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr);
53870034214SMatthew G. Knepley     ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
53970034214SMatthew G. Knepley     ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr);
54070034214SMatthew G. Knepley     ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
54170034214SMatthew G. Knepley     ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr);
54270034214SMatthew G. Knepley   }
54370034214SMatthew G. Knepley   ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr);
54470034214SMatthew G. Knepley   ierr = DMPlexGetConeOrientations(*dmParallel, &newCones);CHKERRQ(ierr);
54570034214SMatthew G. Knepley   ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
54670034214SMatthew G. Knepley   ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr);
54770034214SMatthew G. Knepley   ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr);
54870034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr);
54970034214SMatthew G. Knepley   /* Create supports and stratify sieve */
55070034214SMatthew G. Knepley   {
55170034214SMatthew G. Knepley     PetscInt pStart, pEnd;
55270034214SMatthew G. Knepley 
55370034214SMatthew G. Knepley     ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr);
55470034214SMatthew G. Knepley     ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr);
55570034214SMatthew G. Knepley   }
55670034214SMatthew G. Knepley   ierr = DMPlexSymmetrize(*dmParallel);CHKERRQ(ierr);
55770034214SMatthew G. Knepley   ierr = DMPlexStratify(*dmParallel);CHKERRQ(ierr);
55870034214SMatthew G. Knepley   /* Create point SF for parallel mesh */
55970034214SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
56070034214SMatthew G. Knepley   {
56170034214SMatthew G. Knepley     const PetscInt *leaves;
56270034214SMatthew G. Knepley     PetscSFNode    *remotePoints, *rowners, *lowners;
56370034214SMatthew G. Knepley     PetscInt        numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints;
56470034214SMatthew G. Knepley     PetscInt        pStart, pEnd;
56570034214SMatthew G. Knepley 
56670034214SMatthew G. Knepley     ierr = DMPlexGetChart(*dmParallel, &pStart, &pEnd);CHKERRQ(ierr);
56770034214SMatthew G. Knepley     ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr);
56870034214SMatthew G. Knepley     ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr);
56970034214SMatthew G. Knepley     for (p=0; p<numRoots; p++) {
57070034214SMatthew G. Knepley       rowners[p].rank  = -1;
57170034214SMatthew G. Knepley       rowners[p].index = -1;
57270034214SMatthew G. Knepley     }
57370034214SMatthew G. Knepley     if (origCellPart) {
57470034214SMatthew G. Knepley       /* Make sure points in the original partition are not assigned to other procs */
57570034214SMatthew G. Knepley       const PetscInt *origPoints;
57670034214SMatthew G. Knepley 
57770034214SMatthew G. Knepley       ierr = DMPlexCreatePartitionClosure(dm, origCellPartSection, origCellPart, &origPartSection, &origPart);CHKERRQ(ierr);
57870034214SMatthew G. Knepley       ierr = ISGetIndices(origPart, &origPoints);CHKERRQ(ierr);
57970034214SMatthew G. Knepley       for (p = 0; p < numProcs; ++p) {
58070034214SMatthew G. Knepley         PetscInt dof, off, d;
58170034214SMatthew G. Knepley 
58270034214SMatthew G. Knepley         ierr = PetscSectionGetDof(origPartSection, p, &dof);CHKERRQ(ierr);
58370034214SMatthew G. Knepley         ierr = PetscSectionGetOffset(origPartSection, p, &off);CHKERRQ(ierr);
58470034214SMatthew G. Knepley         for (d = off; d < off+dof; ++d) {
58570034214SMatthew G. Knepley           rowners[origPoints[d]].rank = p;
58670034214SMatthew G. Knepley         }
58770034214SMatthew G. Knepley       }
58870034214SMatthew G. Knepley       ierr = ISRestoreIndices(origPart, &origPoints);CHKERRQ(ierr);
58970034214SMatthew G. Knepley       ierr = ISDestroy(&origPart);CHKERRQ(ierr);
59070034214SMatthew G. Knepley       ierr = PetscSectionDestroy(&origPartSection);CHKERRQ(ierr);
59170034214SMatthew G. Knepley     }
59270034214SMatthew G. Knepley     ierr = ISDestroy(&origCellPart);CHKERRQ(ierr);
59370034214SMatthew G. Knepley     ierr = PetscSectionDestroy(&origCellPartSection);CHKERRQ(ierr);
59470034214SMatthew G. Knepley 
59570034214SMatthew G. Knepley     ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
59670034214SMatthew G. Knepley     ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
59770034214SMatthew G. Knepley     for (p = 0; p < numLeaves; ++p) {
59870034214SMatthew G. Knepley       if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */
59970034214SMatthew G. Knepley         lowners[p].rank  = rank;
60070034214SMatthew G. Knepley         lowners[p].index = leaves ? leaves[p] : p;
60170034214SMatthew G. Knepley       } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */
60270034214SMatthew G. Knepley         lowners[p].rank  = -2;
60370034214SMatthew G. Knepley         lowners[p].index = -2;
60470034214SMatthew G. Knepley       }
60570034214SMatthew G. Knepley     }
60670034214SMatthew G. Knepley     for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */
60770034214SMatthew G. Knepley       rowners[p].rank  = -3;
60870034214SMatthew G. Knepley       rowners[p].index = -3;
60970034214SMatthew G. Knepley     }
61070034214SMatthew G. Knepley     ierr = PetscSFReduceBegin(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
61170034214SMatthew G. Knepley     ierr = PetscSFReduceEnd(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr);
61270034214SMatthew G. Knepley     ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
61370034214SMatthew G. Knepley     ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr);
61470034214SMatthew G. Knepley     for (p = 0; p < numLeaves; ++p) {
61570034214SMatthew G. Knepley       if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed");
61670034214SMatthew G. Knepley       if (lowners[p].rank != rank) ++numGhostPoints;
61770034214SMatthew G. Knepley     }
61870034214SMatthew G. Knepley     ierr = PetscMalloc1(numGhostPoints,    &ghostPoints);CHKERRQ(ierr);
61970034214SMatthew G. Knepley     ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr);
62070034214SMatthew G. Knepley     for (p = 0, gp = 0; p < numLeaves; ++p) {
62170034214SMatthew G. Knepley       if (lowners[p].rank != rank) {
62270034214SMatthew G. Knepley         ghostPoints[gp]        = leaves ? leaves[p] : p;
62370034214SMatthew G. Knepley         remotePoints[gp].rank  = lowners[p].rank;
62470034214SMatthew G. Knepley         remotePoints[gp].index = lowners[p].index;
62570034214SMatthew G. Knepley         ++gp;
62670034214SMatthew G. Knepley       }
62770034214SMatthew G. Knepley     }
62870034214SMatthew G. Knepley     ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr);
62970034214SMatthew G. Knepley     ierr = PetscSFSetGraph((*dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr);
63070034214SMatthew G. Knepley     ierr = PetscSFSetFromOptions((*dmParallel)->sf);CHKERRQ(ierr);
63170034214SMatthew G. Knepley   }
632a42b08eeSMatthew G. Knepley   pmesh->useCone    = mesh->useCone;
633a42b08eeSMatthew G. Knepley   pmesh->useClosure = mesh->useClosure;
63470034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr);
635bf5244c7SMatthew G. Knepley   /* Distribute Coordinates */
636bf5244c7SMatthew G. Knepley   {
637bf5244c7SMatthew G. Knepley     PetscSection     originalCoordSection, newCoordSection;
638bf5244c7SMatthew G. Knepley     Vec              originalCoordinates, newCoordinates;
639bf5244c7SMatthew G. Knepley     PetscInt         bs;
640bf5244c7SMatthew G. Knepley     const char      *name;
641bf5244c7SMatthew G. Knepley     const PetscReal *maxCell, *L;
642bf5244c7SMatthew G. Knepley 
643bf5244c7SMatthew G. Knepley     ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr);
644bf5244c7SMatthew G. Knepley     ierr = DMGetCoordinateSection(*dmParallel, &newCoordSection);CHKERRQ(ierr);
645bf5244c7SMatthew G. Knepley     ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr);
646bf5244c7SMatthew G. Knepley     ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr);
647bf5244c7SMatthew G. Knepley     ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr);
648bf5244c7SMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr);
649bf5244c7SMatthew G. Knepley 
650bf5244c7SMatthew G. Knepley     ierr = DMPlexDistributeField(dm, pointSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr);
651bf5244c7SMatthew G. Knepley     ierr = DMSetCoordinatesLocal(*dmParallel, newCoordinates);CHKERRQ(ierr);
652bf5244c7SMatthew G. Knepley     ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr);
653bf5244c7SMatthew G. Knepley     ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr);
654bf5244c7SMatthew G. Knepley     ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr);
655bf5244c7SMatthew G. Knepley     ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr);
656bf5244c7SMatthew G. Knepley     if (L) {ierr = DMSetPeriodicity(*dmParallel, maxCell, L);CHKERRQ(ierr);}
657bf5244c7SMatthew G. Knepley   }
658bf5244c7SMatthew G. Knepley   /* Distribute labels */
659bf5244c7SMatthew G. Knepley   ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
660bf5244c7SMatthew G. Knepley   {
661bf5244c7SMatthew G. Knepley     DMLabel  next      = mesh->labels, newNext = pmesh->labels;
662bf5244c7SMatthew G. Knepley     PetscInt numLabels = 0, l;
663bf5244c7SMatthew G. Knepley 
664bf5244c7SMatthew G. Knepley     /* Bcast number of labels */
665bf5244c7SMatthew G. Knepley     while (next) {++numLabels; next = next->next;}
666bf5244c7SMatthew G. Knepley     ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr);
667bf5244c7SMatthew G. Knepley     next = mesh->labels;
668bf5244c7SMatthew G. Knepley     for (l = 0; l < numLabels; ++l) {
669bf5244c7SMatthew G. Knepley       DMLabel   labelNew;
670bf5244c7SMatthew G. Knepley       PetscBool isdepth;
671bf5244c7SMatthew G. Knepley 
672bf5244c7SMatthew G. Knepley       /* Skip "depth" because it is recreated */
673bf5244c7SMatthew G. Knepley       if (!rank) {ierr = PetscStrcmp(next->name, "depth", &isdepth);CHKERRQ(ierr);}
674bf5244c7SMatthew G. Knepley       ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr);
675bf5244c7SMatthew G. Knepley       if (isdepth) {if (!rank) next = next->next; continue;}
676bf5244c7SMatthew G. Knepley       ierr = DMLabelDistribute(next, partSection, part, renumbering, &labelNew);CHKERRQ(ierr);
677bf5244c7SMatthew G. Knepley       /* Insert into list */
678bf5244c7SMatthew G. Knepley       if (newNext) newNext->next = labelNew;
679bf5244c7SMatthew G. Knepley       else         pmesh->labels = labelNew;
680bf5244c7SMatthew G. Knepley       newNext = labelNew;
681bf5244c7SMatthew G. Knepley       if (!rank) next = next->next;
682bf5244c7SMatthew G. Knepley     }
683bf5244c7SMatthew G. Knepley   }
684bf5244c7SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr);
685bf5244c7SMatthew G. Knepley   /* Setup hybrid structure */
686bf5244c7SMatthew G. Knepley   {
687bf5244c7SMatthew G. Knepley     const PetscInt *gpoints;
688bf5244c7SMatthew G. Knepley     PetscInt        depth, n, d;
689bf5244c7SMatthew G. Knepley 
690bf5244c7SMatthew G. Knepley     for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];}
691bf5244c7SMatthew G. Knepley     ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr);
692bf5244c7SMatthew G. Knepley     ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr);
693bf5244c7SMatthew G. Knepley     ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr);
694bf5244c7SMatthew G. Knepley     ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr);
695bf5244c7SMatthew G. Knepley     for (d = 0; d <= dim; ++d) {
696bf5244c7SMatthew G. Knepley       PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p;
697bf5244c7SMatthew G. Knepley 
698bf5244c7SMatthew G. Knepley       if (pmax < 0) continue;
699bf5244c7SMatthew G. Knepley       ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr);
700bf5244c7SMatthew G. Knepley       ierr = DMPlexGetDepthStratum(*dmParallel, d, NULL, &pEnd);CHKERRQ(ierr);
701bf5244c7SMatthew G. Knepley       ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr);
702bf5244c7SMatthew G. Knepley       for (p = 0; p < n; ++p) {
703bf5244c7SMatthew G. Knepley         const PetscInt point = gpoints[p];
704bf5244c7SMatthew G. Knepley 
705bf5244c7SMatthew G. Knepley         if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax;
706bf5244c7SMatthew G. Knepley       }
707bf5244c7SMatthew G. Knepley       if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax;
708bf5244c7SMatthew G. Knepley       else            pmesh->hybridPointMax[d] = -1;
709bf5244c7SMatthew G. Knepley     }
710bf5244c7SMatthew G. Knepley     ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr);
711bf5244c7SMatthew G. Knepley   }
712bf5244c7SMatthew G. Knepley   /* Cleanup Partition */
713bf5244c7SMatthew G. Knepley   ierr = ISLocalToGlobalMappingDestroy(&renumbering);CHKERRQ(ierr);
714bf5244c7SMatthew G. Knepley   ierr = PetscSFDestroy(&partSF);CHKERRQ(ierr);
715bf5244c7SMatthew G. Knepley   ierr = PetscSectionDestroy(&partSection);CHKERRQ(ierr);
716bf5244c7SMatthew G. Knepley   ierr = ISDestroy(&part);CHKERRQ(ierr);
71786719b60SMatthew G. Knepley   /* Copy BC */
71886719b60SMatthew G. Knepley   ierr = DMPlexCopyBoundary(dm, *dmParallel);CHKERRQ(ierr);
71970034214SMatthew G. Knepley   /* Cleanup */
72070034214SMatthew G. Knepley   if (sf) {*sf = pointSF;}
72170034214SMatthew G. Knepley   else    {ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr);}
72270034214SMatthew G. Knepley   ierr = DMSetFromOptions(*dmParallel);CHKERRQ(ierr);
72370034214SMatthew G. Knepley   ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr);
72470034214SMatthew G. Knepley   PetscFunctionReturn(0);
72570034214SMatthew G. Knepley }
726