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 } 266b0a623aaSMatthew G. Knepley #undef __FUNCT__ 267b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexCreateTwoSidedProcessSF" 268b0a623aaSMatthew G. Knepley /*@ 269b0a623aaSMatthew G. Knepley DMPlexCreateTwoSidedProcessSF - Create an SF which just has process connectivity 270b0a623aaSMatthew G. Knepley 271b0a623aaSMatthew G. Knepley Collective on DM 272b0a623aaSMatthew G. Knepley 273b0a623aaSMatthew G. Knepley Input Parameters: 274b0a623aaSMatthew G. Knepley + dm - The DM 275b0a623aaSMatthew G. Knepley - sfPoint - The PetscSF which encodes point connectivity 276b0a623aaSMatthew G. Knepley 277b0a623aaSMatthew G. Knepley Output Parameters: 278b0a623aaSMatthew G. Knepley + processRanks - A list of process neighbors, or NULL 279b0a623aaSMatthew G. Knepley - sfProcess - An SF encoding the two-sided process connectivity, or NULL 280b0a623aaSMatthew G. Knepley 281b0a623aaSMatthew G. Knepley Level: developer 282b0a623aaSMatthew G. Knepley 283b0a623aaSMatthew G. Knepley .seealso: PetscSFCreate(), DMPlexCreateProcessSF() 284b0a623aaSMatthew G. Knepley @*/ 285b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM dm, PetscSF sfPoint, PetscSection rootRankSection, IS rootRanks, PetscSection leafRankSection, IS leafRanks, IS *processRanks, PetscSF *sfProcess) 286b0a623aaSMatthew G. Knepley { 287b0a623aaSMatthew G. Knepley const PetscSFNode *remotePoints; 288b0a623aaSMatthew G. Knepley PetscInt *localPointsNew; 289b0a623aaSMatthew G. Knepley PetscSFNode *remotePointsNew; 290b0a623aaSMatthew G. Knepley const PetscInt *nranks; 291b0a623aaSMatthew G. Knepley PetscInt *ranksNew; 292b0a623aaSMatthew G. Knepley PetscBT neighbors; 293b0a623aaSMatthew G. Knepley PetscInt pStart, pEnd, p, numLeaves, l, numNeighbors, n; 294b0a623aaSMatthew G. Knepley PetscMPIInt numProcs, proc, rank; 295b0a623aaSMatthew G. Knepley PetscErrorCode ierr; 296b0a623aaSMatthew G. Knepley 297b0a623aaSMatthew G. Knepley PetscFunctionBegin; 298b0a623aaSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 299b0a623aaSMatthew G. Knepley PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2); 300b0a623aaSMatthew G. Knepley if (processRanks) {PetscValidPointer(processRanks, 3);} 301b0a623aaSMatthew G. Knepley if (sfProcess) {PetscValidPointer(sfProcess, 4);} 302b0a623aaSMatthew G. Knepley ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr); 303b0a623aaSMatthew G. Knepley ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 304b0a623aaSMatthew G. Knepley ierr = PetscSFGetGraph(sfPoint, NULL, &numLeaves, NULL, &remotePoints);CHKERRQ(ierr); 305b0a623aaSMatthew G. Knepley ierr = PetscBTCreate(numProcs, &neighbors);CHKERRQ(ierr); 306b0a623aaSMatthew G. Knepley ierr = PetscBTMemzero(numProcs, neighbors);CHKERRQ(ierr); 307b0a623aaSMatthew G. Knepley /* Compute root-to-leaf process connectivity */ 308b0a623aaSMatthew G. Knepley ierr = PetscSectionGetChart(rootRankSection, &pStart, &pEnd);CHKERRQ(ierr); 309b0a623aaSMatthew G. Knepley ierr = ISGetIndices(rootRanks, &nranks);CHKERRQ(ierr); 310b0a623aaSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 311b0a623aaSMatthew G. Knepley PetscInt ndof, noff, n; 312b0a623aaSMatthew G. Knepley 313b0a623aaSMatthew G. Knepley ierr = PetscSectionGetDof(rootRankSection, p, &ndof);CHKERRQ(ierr); 314b0a623aaSMatthew G. Knepley ierr = PetscSectionGetOffset(rootRankSection, p, &noff);CHKERRQ(ierr); 315b0a623aaSMatthew G. Knepley for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);} 316b0a623aaSMatthew G. Knepley } 317b0a623aaSMatthew G. Knepley ierr = ISRestoreIndices(rootRanks, &nranks);CHKERRQ(ierr); 318b0a623aaSMatthew G. Knepley /* Compute leaf-to-neighbor process connectivity */ 319b0a623aaSMatthew G. Knepley ierr = PetscSectionGetChart(leafRankSection, &pStart, &pEnd);CHKERRQ(ierr); 320b0a623aaSMatthew G. Knepley ierr = ISGetIndices(leafRanks, &nranks);CHKERRQ(ierr); 321b0a623aaSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 322b0a623aaSMatthew G. Knepley PetscInt ndof, noff, n; 323b0a623aaSMatthew G. Knepley 324b0a623aaSMatthew G. Knepley ierr = PetscSectionGetDof(leafRankSection, p, &ndof);CHKERRQ(ierr); 325b0a623aaSMatthew G. Knepley ierr = PetscSectionGetOffset(leafRankSection, p, &noff);CHKERRQ(ierr); 326b0a623aaSMatthew G. Knepley for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);} 327b0a623aaSMatthew G. Knepley } 328b0a623aaSMatthew G. Knepley ierr = ISRestoreIndices(leafRanks, &nranks);CHKERRQ(ierr); 329b0a623aaSMatthew G. Knepley /* Compute leaf-to-root process connectivity */ 330b0a623aaSMatthew G. Knepley for (l = 0; l < numLeaves; ++l) {PetscBTSet(neighbors, remotePoints[l].rank);} 331b0a623aaSMatthew G. Knepley /* Calculate edges */ 332b0a623aaSMatthew G. Knepley PetscBTClear(neighbors, rank); 333b0a623aaSMatthew G. Knepley for(proc = 0, numNeighbors = 0; proc < numProcs; ++proc) {if (PetscBTLookup(neighbors, proc)) ++numNeighbors;} 334b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(numNeighbors, &ranksNew);CHKERRQ(ierr); 335b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(numNeighbors, &localPointsNew);CHKERRQ(ierr); 336b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(numNeighbors, &remotePointsNew);CHKERRQ(ierr); 337b0a623aaSMatthew G. Knepley for(proc = 0, n = 0; proc < numProcs; ++proc) { 338b0a623aaSMatthew G. Knepley if (PetscBTLookup(neighbors, proc)) { 339b0a623aaSMatthew G. Knepley ranksNew[n] = proc; 340b0a623aaSMatthew G. Knepley localPointsNew[n] = proc; 341b0a623aaSMatthew G. Knepley remotePointsNew[n].index = rank; 342b0a623aaSMatthew G. Knepley remotePointsNew[n].rank = proc; 343b0a623aaSMatthew G. Knepley ++n; 344b0a623aaSMatthew G. Knepley } 345b0a623aaSMatthew G. Knepley } 346b0a623aaSMatthew G. Knepley ierr = PetscBTDestroy(&neighbors);CHKERRQ(ierr); 347b0a623aaSMatthew G. Knepley if (processRanks) {ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), numNeighbors, ranksNew, PETSC_OWN_POINTER, processRanks);CHKERRQ(ierr);} 348b0a623aaSMatthew G. Knepley else {ierr = PetscFree(ranksNew);CHKERRQ(ierr);} 349b0a623aaSMatthew G. Knepley if (sfProcess) { 350b0a623aaSMatthew G. Knepley ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess);CHKERRQ(ierr); 351b0a623aaSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *sfProcess, "Two-Sided Process SF");CHKERRQ(ierr); 352b0a623aaSMatthew G. Knepley ierr = PetscSFSetFromOptions(*sfProcess);CHKERRQ(ierr); 353b0a623aaSMatthew G. Knepley ierr = PetscSFSetGraph(*sfProcess, numProcs, numNeighbors, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr); 354b0a623aaSMatthew G. Knepley } 355b0a623aaSMatthew G. Knepley PetscFunctionReturn(0); 356b0a623aaSMatthew G. Knepley } 357b0a623aaSMatthew G. Knepley 358b0a623aaSMatthew G. Knepley #undef __FUNCT__ 359b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeOwnership" 360b0a623aaSMatthew G. Knepley /*@ 361b0a623aaSMatthew G. Knepley DMPlexDistributeOwnership - Compute owner information for shared points. This basically gets two-sided for an SF. 362b0a623aaSMatthew G. Knepley 363b0a623aaSMatthew G. Knepley Collective on DM 364b0a623aaSMatthew G. Knepley 365b0a623aaSMatthew G. Knepley Input Parameter: 366b0a623aaSMatthew G. Knepley . dm - The DM 367b0a623aaSMatthew G. Knepley 368b0a623aaSMatthew G. Knepley Output Parameters: 369b0a623aaSMatthew G. Knepley + rootSection - The number of leaves for a given root point 370b0a623aaSMatthew G. Knepley . rootrank - The rank of each edge into the root point 371b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point 372b0a623aaSMatthew G. Knepley - leafrank - The rank of each process sharing a leaf point 373b0a623aaSMatthew G. Knepley 374b0a623aaSMatthew G. Knepley Level: developer 375b0a623aaSMatthew G. Knepley 376b0a623aaSMatthew G. Knepley .seealso: DMPlexCreateOverlap() 377b0a623aaSMatthew G. Knepley @*/ 378b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexDistributeOwnership(DM dm, PetscSection rootSection, IS *rootrank, PetscSection leafSection, IS *leafrank) 379b0a623aaSMatthew G. Knepley { 380b0a623aaSMatthew G. Knepley MPI_Comm comm; 381b0a623aaSMatthew G. Knepley PetscSF sfPoint; 382b0a623aaSMatthew G. Knepley const PetscInt *rootdegree; 383b0a623aaSMatthew G. Knepley PetscInt *myrank, *remoterank; 384b0a623aaSMatthew G. Knepley PetscInt pStart, pEnd, p, nedges; 385b0a623aaSMatthew G. Knepley PetscMPIInt rank; 386b0a623aaSMatthew G. Knepley PetscErrorCode ierr; 387b0a623aaSMatthew G. Knepley 388b0a623aaSMatthew G. Knepley PetscFunctionBegin; 389b0a623aaSMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 390b0a623aaSMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 391b0a623aaSMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 392b0a623aaSMatthew G. Knepley ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 393b0a623aaSMatthew G. Knepley /* Compute number of leaves for each root */ 394b0a623aaSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) rootSection, "Root Section");CHKERRQ(ierr); 395b0a623aaSMatthew G. Knepley ierr = PetscSectionSetChart(rootSection, pStart, pEnd);CHKERRQ(ierr); 396b0a623aaSMatthew G. Knepley ierr = PetscSFComputeDegreeBegin(sfPoint, &rootdegree);CHKERRQ(ierr); 397b0a623aaSMatthew G. Knepley ierr = PetscSFComputeDegreeEnd(sfPoint, &rootdegree);CHKERRQ(ierr); 398b0a623aaSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) {ierr = PetscSectionSetDof(rootSection, p, rootdegree[p-pStart]);CHKERRQ(ierr);} 399b0a623aaSMatthew G. Knepley ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 400b0a623aaSMatthew G. Knepley /* Gather rank of each leaf to root */ 401b0a623aaSMatthew G. Knepley ierr = PetscSectionGetStorageSize(rootSection, &nedges);CHKERRQ(ierr); 402b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(pEnd-pStart, &myrank);CHKERRQ(ierr); 403b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(nedges, &remoterank);CHKERRQ(ierr); 404b0a623aaSMatthew G. Knepley for (p = 0; p < pEnd-pStart; ++p) myrank[p] = rank; 405b0a623aaSMatthew G. Knepley ierr = PetscSFGatherBegin(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr); 406b0a623aaSMatthew G. Knepley ierr = PetscSFGatherEnd(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr); 407b0a623aaSMatthew G. Knepley ierr = PetscFree(myrank);CHKERRQ(ierr); 408b0a623aaSMatthew G. Knepley ierr = ISCreateGeneral(comm, nedges, remoterank, PETSC_OWN_POINTER, rootrank);CHKERRQ(ierr); 409b0a623aaSMatthew G. Knepley /* Distribute remote ranks to leaves */ 410b0a623aaSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) leafSection, "Leaf Section");CHKERRQ(ierr); 411b0a623aaSMatthew G. Knepley ierr = DMPlexDistributeFieldIS(dm, sfPoint, rootSection, *rootrank, leafSection, leafrank);CHKERRQ(ierr); 412b0a623aaSMatthew G. Knepley PetscFunctionReturn(0); 413b0a623aaSMatthew G. Knepley } 414b0a623aaSMatthew G. Knepley 415b0a623aaSMatthew G. Knepley #undef __FUNCT__ 416b0a623aaSMatthew G. Knepley #define __FUNCT__ "DMPlexCreateOverlap" 417*278397a0SMatthew G. Knepley /*@C 418b0a623aaSMatthew G. Knepley DMPlexCreateOverlap - Compute owner information for shared points. This basically gets two-sided for an SF. 419b0a623aaSMatthew G. Knepley 420b0a623aaSMatthew G. Knepley Collective on DM 421b0a623aaSMatthew G. Knepley 422b0a623aaSMatthew G. Knepley Input Parameters: 423b0a623aaSMatthew G. Knepley + dm - The DM 424b0a623aaSMatthew G. Knepley . rootSection - The number of leaves for a given root point 425b0a623aaSMatthew G. Knepley . rootrank - The rank of each edge into the root point 426b0a623aaSMatthew G. Knepley . leafSection - The number of processes sharing a given leaf point 427b0a623aaSMatthew G. Knepley - leafrank - The rank of each process sharing a leaf point 428b0a623aaSMatthew G. Knepley 429b0a623aaSMatthew G. Knepley Output Parameters: 430b0a623aaSMatthew G. Knepley + ovRootSection - The number of new overlap points for each neighboring process 431b0a623aaSMatthew G. Knepley . ovRootPoints - The new overlap points for each neighboring process 432b0a623aaSMatthew G. Knepley . ovLeafSection - The number of new overlap points from each neighboring process 433b0a623aaSMatthew G. Knepley - ovLeafPoints - The new overlap points from each neighboring process 434b0a623aaSMatthew G. Knepley 435b0a623aaSMatthew G. Knepley Level: developer 436b0a623aaSMatthew G. Knepley 437b0a623aaSMatthew G. Knepley .seealso: DMPlexDistributeOwnership() 438b0a623aaSMatthew G. Knepley @*/ 439b0a623aaSMatthew G. Knepley PetscErrorCode DMPlexCreateOverlap(DM dm, PetscSection rootSection, IS rootrank, PetscSection leafSection, IS leafrank, PetscSection ovRootSection, PetscSFNode **ovRootPoints, PetscSection ovLeafSection, PetscSFNode **ovLeafPoints) 440b0a623aaSMatthew G. Knepley { 441b0a623aaSMatthew G. Knepley DMLabel ovAdjByRank; /* A DMLabel containing all points adjacent to shared points, separated by rank (value in label) */ 442b0a623aaSMatthew G. Knepley PetscSF sfPoint, sfProc; 443b0a623aaSMatthew G. Knepley IS valueIS; 444b0a623aaSMatthew G. Knepley const PetscSFNode *remote; 445b0a623aaSMatthew G. Knepley const PetscInt *local; 446b0a623aaSMatthew G. Knepley const PetscInt *nrank, *rrank, *neighbors; 447b0a623aaSMatthew G. Knepley PetscInt *adj = NULL; 448b0a623aaSMatthew G. Knepley PetscInt pStart, pEnd, p, sStart, sEnd, nleaves, l, numNeighbors, n, ovSize; 449b0a623aaSMatthew G. Knepley PetscMPIInt rank, numProcs; 450b0a623aaSMatthew G. Knepley PetscErrorCode ierr; 451b0a623aaSMatthew G. Knepley 452b0a623aaSMatthew G. Knepley PetscFunctionBegin; 453b0a623aaSMatthew G. Knepley ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr); 454b0a623aaSMatthew G. Knepley ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 455b0a623aaSMatthew G. Knepley ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 456b0a623aaSMatthew G. Knepley ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 457b0a623aaSMatthew G. Knepley ierr = PetscSectionGetChart(leafSection, &sStart, &sEnd);CHKERRQ(ierr); 458b0a623aaSMatthew G. Knepley ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr); 459b0a623aaSMatthew G. Knepley ierr = DMLabelCreate("Overlap adjacency", &ovAdjByRank);CHKERRQ(ierr); 460b0a623aaSMatthew G. Knepley /* Handle leaves: shared with the root point */ 461b0a623aaSMatthew G. Knepley for (l = 0; l < nleaves; ++l) { 462b0a623aaSMatthew G. Knepley PetscInt adjSize = PETSC_DETERMINE, a; 463b0a623aaSMatthew G. Knepley 464b0a623aaSMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, local[l], &adjSize, &adj);CHKERRQ(ierr); 465b0a623aaSMatthew G. Knepley for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remote[l].rank);CHKERRQ(ierr);} 466b0a623aaSMatthew G. Knepley } 467b0a623aaSMatthew G. Knepley ierr = ISGetIndices(rootrank, &rrank);CHKERRQ(ierr); 468b0a623aaSMatthew G. Knepley ierr = ISGetIndices(leafrank, &nrank);CHKERRQ(ierr); 469b0a623aaSMatthew G. Knepley /* Handle roots */ 470b0a623aaSMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 471b0a623aaSMatthew G. Knepley PetscInt adjSize = PETSC_DETERMINE, neighbors = 0, noff, n, a; 472b0a623aaSMatthew G. Knepley 473b0a623aaSMatthew G. Knepley if ((p >= sStart) && (p < sEnd)) { 474b0a623aaSMatthew G. Knepley /* Some leaves share a root with other leaves on different processes */ 475b0a623aaSMatthew G. Knepley ierr = PetscSectionGetDof(leafSection, p, &neighbors);CHKERRQ(ierr); 476b0a623aaSMatthew G. Knepley if (neighbors) { 477b0a623aaSMatthew G. Knepley ierr = PetscSectionGetOffset(leafSection, p, &noff);CHKERRQ(ierr); 478b0a623aaSMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr); 479b0a623aaSMatthew G. Knepley for (n = 0; n < neighbors; ++n) { 480b0a623aaSMatthew G. Knepley const PetscInt remoteRank = nrank[noff+n]; 481b0a623aaSMatthew G. Knepley 482b0a623aaSMatthew G. Knepley if (remoteRank == rank) continue; 483b0a623aaSMatthew G. Knepley for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);} 484b0a623aaSMatthew G. Knepley } 485b0a623aaSMatthew G. Knepley } 486b0a623aaSMatthew G. Knepley } 487b0a623aaSMatthew G. Knepley /* Roots are shared with leaves */ 488b0a623aaSMatthew G. Knepley ierr = PetscSectionGetDof(rootSection, p, &neighbors);CHKERRQ(ierr); 489b0a623aaSMatthew G. Knepley if (!neighbors) continue; 490b0a623aaSMatthew G. Knepley ierr = PetscSectionGetOffset(rootSection, p, &noff);CHKERRQ(ierr); 491b0a623aaSMatthew G. Knepley ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr); 492b0a623aaSMatthew G. Knepley for (n = 0; n < neighbors; ++n) { 493b0a623aaSMatthew G. Knepley const PetscInt remoteRank = rrank[noff+n]; 494b0a623aaSMatthew G. Knepley 495b0a623aaSMatthew G. Knepley if (remoteRank == rank) continue; 496b0a623aaSMatthew G. Knepley for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);} 497b0a623aaSMatthew G. Knepley } 498b0a623aaSMatthew G. Knepley } 499b0a623aaSMatthew G. Knepley ierr = PetscFree(adj);CHKERRQ(ierr); 500b0a623aaSMatthew G. Knepley ierr = ISRestoreIndices(rootrank, &rrank);CHKERRQ(ierr); 501b0a623aaSMatthew G. Knepley ierr = ISRestoreIndices(leafrank, &nrank);CHKERRQ(ierr); 502b0a623aaSMatthew G. Knepley { 503b0a623aaSMatthew G. Knepley ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr); 504b0a623aaSMatthew G. Knepley ierr = DMLabelView(ovAdjByRank, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 505b0a623aaSMatthew G. Knepley } 506b0a623aaSMatthew G. Knepley /* Convert to (point, rank) and use actual owners */ 507b0a623aaSMatthew G. Knepley ierr = PetscSectionSetChart(ovRootSection, 0, numProcs);CHKERRQ(ierr); 508b0a623aaSMatthew G. Knepley ierr = DMLabelGetValueIS(ovAdjByRank, &valueIS);CHKERRQ(ierr); 509b0a623aaSMatthew G. Knepley ierr = ISGetLocalSize(valueIS, &numNeighbors);CHKERRQ(ierr); 510b0a623aaSMatthew G. Knepley ierr = ISGetIndices(valueIS, &neighbors);CHKERRQ(ierr); 511b0a623aaSMatthew G. Knepley for (n = 0; n < numNeighbors; ++n) { 512b0a623aaSMatthew G. Knepley PetscInt numPoints; 513b0a623aaSMatthew G. Knepley 514b0a623aaSMatthew G. Knepley ierr = DMLabelGetStratumSize(ovAdjByRank, neighbors[n], &numPoints);CHKERRQ(ierr); 515b0a623aaSMatthew G. Knepley ierr = PetscSectionAddDof(ovRootSection, neighbors[n], numPoints);CHKERRQ(ierr); 516b0a623aaSMatthew G. Knepley } 517b0a623aaSMatthew G. Knepley ierr = PetscSectionSetUp(ovRootSection);CHKERRQ(ierr); 518b0a623aaSMatthew G. Knepley ierr = PetscSectionGetStorageSize(ovRootSection, &ovSize);CHKERRQ(ierr); 519b0a623aaSMatthew G. Knepley ierr = PetscMalloc1(ovSize, ovRootPoints);CHKERRQ(ierr); 520b0a623aaSMatthew G. Knepley ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr); 521b0a623aaSMatthew G. Knepley for (n = 0; n < numNeighbors; ++n) { 522b0a623aaSMatthew G. Knepley IS pointIS; 523b0a623aaSMatthew G. Knepley const PetscInt *points; 524b0a623aaSMatthew G. Knepley PetscInt off, numPoints, p; 525b0a623aaSMatthew G. Knepley 526b0a623aaSMatthew G. Knepley ierr = PetscSectionGetOffset(ovRootSection, neighbors[n], &off);CHKERRQ(ierr); 527b0a623aaSMatthew G. Knepley ierr = DMLabelGetStratumIS(ovAdjByRank, neighbors[n], &pointIS);CHKERRQ(ierr); 528b0a623aaSMatthew G. Knepley ierr = ISGetLocalSize(pointIS, &numPoints);CHKERRQ(ierr); 529b0a623aaSMatthew G. Knepley ierr = ISGetIndices(pointIS, &points);CHKERRQ(ierr); 530b0a623aaSMatthew G. Knepley for (p = 0; p < numPoints; ++p) { 531b0a623aaSMatthew G. Knepley ierr = PetscFindInt(points[p], nleaves, local, &l);CHKERRQ(ierr); 532b0a623aaSMatthew G. Knepley if (l >= 0) {(*ovRootPoints)[off+p] = remote[l];} 533b0a623aaSMatthew G. Knepley else {(*ovRootPoints)[off+p].index = points[p]; (*ovRootPoints)[off+p].rank = rank;} 534b0a623aaSMatthew G. Knepley } 535b0a623aaSMatthew G. Knepley ierr = ISRestoreIndices(pointIS, &points);CHKERRQ(ierr); 536b0a623aaSMatthew G. Knepley ierr = ISDestroy(&pointIS);CHKERRQ(ierr); 537b0a623aaSMatthew G. Knepley } 538b0a623aaSMatthew G. Knepley ierr = ISRestoreIndices(valueIS, &neighbors);CHKERRQ(ierr); 539b0a623aaSMatthew G. Knepley ierr = ISDestroy(&valueIS);CHKERRQ(ierr); 540b0a623aaSMatthew G. Knepley ierr = DMLabelDestroy(&ovAdjByRank);CHKERRQ(ierr); 541b0a623aaSMatthew G. Knepley /* Make process SF */ 542b0a623aaSMatthew G. Knepley ierr = DMPlexCreateTwoSidedProcessSF(dm, sfPoint, rootSection, rootrank, leafSection, leafrank, NULL, &sfProc);CHKERRQ(ierr); 543b0a623aaSMatthew G. Knepley { 544b0a623aaSMatthew G. Knepley ierr = PetscSFView(sfProc, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 545b0a623aaSMatthew G. Knepley } 546b0a623aaSMatthew G. Knepley /* Communicate overlap */ 547b0a623aaSMatthew G. Knepley ierr = DMPlexDistributeData(dm, sfProc, ovRootSection, MPIU_2INT, (void *) *ovRootPoints, ovLeafSection, (void **) ovLeafPoints);CHKERRQ(ierr); 548b0a623aaSMatthew G. Knepley ierr = PetscSFDestroy(&sfProc);CHKERRQ(ierr); 549b0a623aaSMatthew G. Knepley PetscFunctionReturn(0); 550b0a623aaSMatthew G. Knepley } 55170034214SMatthew G. Knepley 55270034214SMatthew G. Knepley #undef __FUNCT__ 55370034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeField" 55470034214SMatthew G. Knepley /*@ 55570034214SMatthew G. Knepley DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution 55670034214SMatthew G. Knepley 55770034214SMatthew G. Knepley Collective on DM 55870034214SMatthew G. Knepley 55970034214SMatthew G. Knepley Input Parameters: 56070034214SMatthew G. Knepley + dm - The DMPlex object 56170034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern 56270034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout 56370034214SMatthew G. Knepley - originalVec - The existing data 56470034214SMatthew G. Knepley 56570034214SMatthew G. Knepley Output Parameters: 56670034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout 56770034214SMatthew G. Knepley - newVec - The new data 56870034214SMatthew G. Knepley 56970034214SMatthew G. Knepley Level: developer 57070034214SMatthew G. Knepley 5711e8fc25dSMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeFieldIS(), DMPlexDistributeData() 57270034214SMatthew G. Knepley @*/ 57370034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec) 57470034214SMatthew G. Knepley { 57570034214SMatthew G. Knepley PetscSF fieldSF; 57670034214SMatthew G. Knepley PetscInt *remoteOffsets, fieldSize; 57770034214SMatthew G. Knepley PetscScalar *originalValues, *newValues; 57870034214SMatthew G. Knepley PetscErrorCode ierr; 57970034214SMatthew G. Knepley 58070034214SMatthew G. Knepley PetscFunctionBegin; 58170034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 58270034214SMatthew G. Knepley ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); 58370034214SMatthew G. Knepley 58470034214SMatthew G. Knepley ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); 58570034214SMatthew G. Knepley ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr); 58670034214SMatthew G. Knepley ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr); 58770034214SMatthew G. Knepley 58870034214SMatthew G. Knepley ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr); 58970034214SMatthew G. Knepley ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr); 59070034214SMatthew G. Knepley ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); 59170034214SMatthew G. Knepley ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr); 59270034214SMatthew G. Knepley ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr); 59370034214SMatthew G. Knepley ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); 59470034214SMatthew G. Knepley ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr); 59570034214SMatthew G. Knepley ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr); 59670034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 59770034214SMatthew G. Knepley PetscFunctionReturn(0); 59870034214SMatthew G. Knepley } 59970034214SMatthew G. Knepley 60070034214SMatthew G. Knepley #undef __FUNCT__ 6011e8fc25dSMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeFieldIS" 6021e8fc25dSMatthew G. Knepley /*@ 6031e8fc25dSMatthew G. Knepley DMPlexDistributeFieldIS - Distribute field data to match a given PetscSF, usually the SF from mesh distribution 6041e8fc25dSMatthew G. Knepley 6051e8fc25dSMatthew G. Knepley Collective on DM 6061e8fc25dSMatthew G. Knepley 6071e8fc25dSMatthew G. Knepley Input Parameters: 6081e8fc25dSMatthew G. Knepley + dm - The DMPlex object 6091e8fc25dSMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern 6101e8fc25dSMatthew G. Knepley . originalSection - The PetscSection for existing data layout 6111e8fc25dSMatthew G. Knepley - originalIS - The existing data 6121e8fc25dSMatthew G. Knepley 6131e8fc25dSMatthew G. Knepley Output Parameters: 6141e8fc25dSMatthew G. Knepley + newSection - The PetscSF describing the new data layout 6151e8fc25dSMatthew G. Knepley - newIS - The new data 6161e8fc25dSMatthew G. Knepley 6171e8fc25dSMatthew G. Knepley Level: developer 6181e8fc25dSMatthew G. Knepley 6191e8fc25dSMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexDistributeData() 6201e8fc25dSMatthew G. Knepley @*/ 6211e8fc25dSMatthew G. Knepley PetscErrorCode DMPlexDistributeFieldIS(DM dm, PetscSF pointSF, PetscSection originalSection, IS originalIS, PetscSection newSection, IS *newIS) 6221e8fc25dSMatthew G. Knepley { 6231e8fc25dSMatthew G. Knepley PetscSF fieldSF; 6241e8fc25dSMatthew G. Knepley PetscInt *newValues, *remoteOffsets, fieldSize; 6251e8fc25dSMatthew G. Knepley const PetscInt *originalValues; 6261e8fc25dSMatthew G. Knepley PetscErrorCode ierr; 6271e8fc25dSMatthew G. Knepley 6281e8fc25dSMatthew G. Knepley PetscFunctionBegin; 6291e8fc25dSMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 6301e8fc25dSMatthew G. Knepley ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); 6311e8fc25dSMatthew G. Knepley 6321e8fc25dSMatthew G. Knepley ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); 6331e8fc25dSMatthew G. Knepley ierr = PetscMalloc(fieldSize * sizeof(PetscInt), &newValues);CHKERRQ(ierr); 6341e8fc25dSMatthew G. Knepley 6351e8fc25dSMatthew G. Knepley ierr = ISGetIndices(originalIS, &originalValues);CHKERRQ(ierr); 6361e8fc25dSMatthew G. Knepley ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); 6371e8fc25dSMatthew G. Knepley ierr = PetscSFBcastBegin(fieldSF, MPIU_INT, originalValues, newValues);CHKERRQ(ierr); 6381e8fc25dSMatthew G. Knepley ierr = PetscSFBcastEnd(fieldSF, MPIU_INT, originalValues, newValues);CHKERRQ(ierr); 6391e8fc25dSMatthew G. Knepley ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); 6401e8fc25dSMatthew G. Knepley ierr = ISRestoreIndices(originalIS, &originalValues);CHKERRQ(ierr); 6411e8fc25dSMatthew G. Knepley ierr = ISCreateGeneral(PetscObjectComm((PetscObject) pointSF), fieldSize, newValues, PETSC_OWN_POINTER, newIS);CHKERRQ(ierr); 6421e8fc25dSMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 6431e8fc25dSMatthew G. Knepley PetscFunctionReturn(0); 6441e8fc25dSMatthew G. Knepley } 6451e8fc25dSMatthew G. Knepley 6461e8fc25dSMatthew G. Knepley #undef __FUNCT__ 64770034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeData" 64870034214SMatthew G. Knepley /*@ 64970034214SMatthew G. Knepley DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution 65070034214SMatthew G. Knepley 65170034214SMatthew G. Knepley Collective on DM 65270034214SMatthew G. Knepley 65370034214SMatthew G. Knepley Input Parameters: 65470034214SMatthew G. Knepley + dm - The DMPlex object 65570034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern 65670034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout 65770034214SMatthew G. Knepley . datatype - The type of data 65870034214SMatthew G. Knepley - originalData - The existing data 65970034214SMatthew G. Knepley 66070034214SMatthew G. Knepley Output Parameters: 6611e8fc25dSMatthew G. Knepley + newSection - The PetscSection describing the new data layout 66270034214SMatthew G. Knepley - newData - The new data 66370034214SMatthew G. Knepley 66470034214SMatthew G. Knepley Level: developer 66570034214SMatthew G. Knepley 66670034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField() 66770034214SMatthew G. Knepley @*/ 66870034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData) 66970034214SMatthew G. Knepley { 67070034214SMatthew G. Knepley PetscSF fieldSF; 67170034214SMatthew G. Knepley PetscInt *remoteOffsets, fieldSize; 67270034214SMatthew G. Knepley PetscMPIInt dataSize; 67370034214SMatthew G. Knepley PetscErrorCode ierr; 67470034214SMatthew G. Knepley 67570034214SMatthew G. Knepley PetscFunctionBegin; 67670034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr); 67770034214SMatthew G. Knepley ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); 67870034214SMatthew G. Knepley 67970034214SMatthew G. Knepley ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); 68070034214SMatthew G. Knepley ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr); 68170034214SMatthew G. Knepley ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr); 68270034214SMatthew G. Knepley 68370034214SMatthew G. Knepley ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); 68470034214SMatthew G. Knepley ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr); 68570034214SMatthew G. Knepley ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr); 68670034214SMatthew G. Knepley ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); 68770034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr); 68870034214SMatthew G. Knepley PetscFunctionReturn(0); 68970034214SMatthew G. Knepley } 69070034214SMatthew G. Knepley 69170034214SMatthew G. Knepley #undef __FUNCT__ 69270034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistribute" 69370034214SMatthew G. Knepley /*@C 69470034214SMatthew G. Knepley DMPlexDistribute - Distributes the mesh and any associated sections. 69570034214SMatthew G. Knepley 69670034214SMatthew G. Knepley Not Collective 69770034214SMatthew G. Knepley 69870034214SMatthew G. Knepley Input Parameter: 69970034214SMatthew G. Knepley + dm - The original DMPlex object 70070034214SMatthew G. Knepley . partitioner - The partitioning package, or NULL for the default 70170034214SMatthew G. Knepley - overlap - The overlap of partitions, 0 is the default 70270034214SMatthew G. Knepley 70370034214SMatthew G. Knepley Output Parameter: 70470034214SMatthew G. Knepley + sf - The PetscSF used for point distribution 70570034214SMatthew G. Knepley - parallelMesh - The distributed DMPlex object, or NULL 70670034214SMatthew G. Knepley 707a9c22940SMatthew G. Knepley Note: If the mesh was not distributed, the return value is NULL. 708a9c22940SMatthew G. Knepley 709a9c22940SMatthew G. Knepley The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and 710a9c22940SMatthew G. Knepley DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function 711a9c22940SMatthew G. Knepley representation on the mesh. 71270034214SMatthew G. Knepley 71370034214SMatthew G. Knepley Level: intermediate 71470034214SMatthew G. Knepley 71570034214SMatthew G. Knepley .keywords: mesh, elements 716a9c22940SMatthew G. Knepley .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure() 71770034214SMatthew G. Knepley @*/ 71870034214SMatthew G. Knepley PetscErrorCode DMPlexDistribute(DM dm, const char partitioner[], PetscInt overlap, PetscSF *sf, DM *dmParallel) 71970034214SMatthew G. Knepley { 72070034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data, *pmesh; 72170034214SMatthew G. Knepley MPI_Comm comm; 72270034214SMatthew G. Knepley const PetscInt height = 0; 72370034214SMatthew G. Knepley PetscInt dim, numRemoteRanks; 72470034214SMatthew G. Knepley IS origCellPart, origPart, cellPart, part; 72570034214SMatthew G. Knepley PetscSection origCellPartSection, origPartSection, cellPartSection, partSection; 72670034214SMatthew G. Knepley PetscSFNode *remoteRanks; 72770034214SMatthew G. Knepley PetscSF partSF, pointSF, coneSF; 72870034214SMatthew G. Knepley ISLocalToGlobalMapping renumbering; 72970034214SMatthew G. Knepley PetscSection originalConeSection, newConeSection; 73070034214SMatthew G. Knepley PetscInt *remoteOffsets; 73170034214SMatthew G. Knepley PetscInt *cones, *newCones, newConesSize; 73270034214SMatthew G. Knepley PetscBool flg; 73370034214SMatthew G. Knepley PetscMPIInt rank, numProcs, p; 73470034214SMatthew G. Knepley PetscErrorCode ierr; 73570034214SMatthew G. Knepley 73670034214SMatthew G. Knepley PetscFunctionBegin; 73770034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 73870034214SMatthew G. Knepley if (sf) PetscValidPointer(sf,4); 73970034214SMatthew G. Knepley PetscValidPointer(dmParallel,5); 74070034214SMatthew G. Knepley 74170034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr); 74270034214SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 74370034214SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 74470034214SMatthew G. Knepley ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 74570034214SMatthew G. Knepley 74670034214SMatthew G. Knepley *dmParallel = NULL; 74770034214SMatthew G. Knepley if (numProcs == 1) PetscFunctionReturn(0); 74870034214SMatthew G. Knepley 749c73cfb54SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 75070034214SMatthew G. Knepley /* Create cell partition - We need to rewrite to use IS, use the MatPartition stuff */ 75170034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_Partition,dm,0,0,0);CHKERRQ(ierr); 75270034214SMatthew G. Knepley if (overlap > 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Overlap > 1 not yet implemented"); 75370034214SMatthew G. Knepley ierr = DMPlexCreatePartition(dm, partitioner, height, overlap > 0 ? PETSC_TRUE : PETSC_FALSE, &cellPartSection, &cellPart, &origCellPartSection, &origCellPart);CHKERRQ(ierr); 75470034214SMatthew G. Knepley /* Create SF assuming a serial partition for all processes: Could check for IS length here */ 75570034214SMatthew G. Knepley if (!rank) numRemoteRanks = numProcs; 75670034214SMatthew G. Knepley else numRemoteRanks = 0; 75770034214SMatthew G. Knepley ierr = PetscMalloc1(numRemoteRanks, &remoteRanks);CHKERRQ(ierr); 75870034214SMatthew G. Knepley for (p = 0; p < numRemoteRanks; ++p) { 75970034214SMatthew G. Knepley remoteRanks[p].rank = p; 76070034214SMatthew G. Knepley remoteRanks[p].index = 0; 76170034214SMatthew G. Knepley } 76270034214SMatthew G. Knepley ierr = PetscSFCreate(comm, &partSF);CHKERRQ(ierr); 76370034214SMatthew G. Knepley ierr = PetscSFSetGraph(partSF, 1, numRemoteRanks, NULL, PETSC_OWN_POINTER, remoteRanks, PETSC_OWN_POINTER);CHKERRQ(ierr); 76470034214SMatthew G. Knepley ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr); 76570034214SMatthew G. Knepley if (flg) { 76670034214SMatthew G. Knepley ierr = PetscPrintf(comm, "Cell Partition:\n");CHKERRQ(ierr); 76770034214SMatthew G. Knepley ierr = PetscSectionView(cellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 76870034214SMatthew G. Knepley ierr = ISView(cellPart, NULL);CHKERRQ(ierr); 76970034214SMatthew G. Knepley if (origCellPart) { 77070034214SMatthew G. Knepley ierr = PetscPrintf(comm, "Original Cell Partition:\n");CHKERRQ(ierr); 77170034214SMatthew G. Knepley ierr = PetscSectionView(origCellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 77270034214SMatthew G. Knepley ierr = ISView(origCellPart, NULL);CHKERRQ(ierr); 77370034214SMatthew G. Knepley } 77470034214SMatthew G. Knepley ierr = PetscSFView(partSF, NULL);CHKERRQ(ierr); 77570034214SMatthew G. Knepley } 77670034214SMatthew G. Knepley /* Close the partition over the mesh */ 77770034214SMatthew G. Knepley ierr = DMPlexCreatePartitionClosure(dm, cellPartSection, cellPart, &partSection, &part);CHKERRQ(ierr); 77870034214SMatthew G. Knepley ierr = ISDestroy(&cellPart);CHKERRQ(ierr); 77970034214SMatthew G. Knepley ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr); 78070034214SMatthew G. Knepley /* Create new mesh */ 78170034214SMatthew G. Knepley ierr = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr); 782c73cfb54SMatthew G. Knepley ierr = DMSetDimension(*dmParallel, dim);CHKERRQ(ierr); 78370034214SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr); 78470034214SMatthew G. Knepley pmesh = (DM_Plex*) (*dmParallel)->data; 78570034214SMatthew G. Knepley /* Distribute sieve points and the global point numbering (replaces creating remote bases) */ 78670034214SMatthew G. Knepley ierr = PetscSFConvertPartition(partSF, partSection, part, &renumbering, &pointSF);CHKERRQ(ierr); 78770034214SMatthew G. Knepley if (flg) { 78870034214SMatthew G. Knepley ierr = PetscPrintf(comm, "Point Partition:\n");CHKERRQ(ierr); 78970034214SMatthew G. Knepley ierr = PetscSectionView(partSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 79070034214SMatthew G. Knepley ierr = ISView(part, NULL);CHKERRQ(ierr); 79170034214SMatthew G. Knepley ierr = PetscSFView(pointSF, NULL);CHKERRQ(ierr); 79270034214SMatthew G. Knepley ierr = PetscPrintf(comm, "Point Renumbering after partition:\n");CHKERRQ(ierr); 79370034214SMatthew G. Knepley ierr = ISLocalToGlobalMappingView(renumbering, NULL);CHKERRQ(ierr); 79470034214SMatthew G. Knepley } 79570034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_Partition,dm,0,0,0);CHKERRQ(ierr); 79670034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr); 79770034214SMatthew G. Knepley /* Distribute cone section */ 79870034214SMatthew G. Knepley ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr); 79970034214SMatthew G. Knepley ierr = DMPlexGetConeSection(*dmParallel, &newConeSection);CHKERRQ(ierr); 80070034214SMatthew G. Knepley ierr = PetscSFDistributeSection(pointSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr); 80170034214SMatthew G. Knepley ierr = DMSetUp(*dmParallel);CHKERRQ(ierr); 80270034214SMatthew G. Knepley { 80370034214SMatthew G. Knepley PetscInt pStart, pEnd, p; 80470034214SMatthew G. Knepley 80570034214SMatthew G. Knepley ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr); 80670034214SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 80770034214SMatthew G. Knepley PetscInt coneSize; 80870034214SMatthew G. Knepley ierr = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr); 80970034214SMatthew G. Knepley pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize); 81070034214SMatthew G. Knepley } 81170034214SMatthew G. Knepley } 81270034214SMatthew G. Knepley /* Communicate and renumber cones */ 81370034214SMatthew G. Knepley ierr = PetscSFCreateSectionSF(pointSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr); 81470034214SMatthew G. Knepley ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr); 81570034214SMatthew G. Knepley ierr = DMPlexGetCones(*dmParallel, &newCones);CHKERRQ(ierr); 81670034214SMatthew G. Knepley ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 81770034214SMatthew G. Knepley ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 81870034214SMatthew G. Knepley ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr); 8199d90f715SBarry Smith ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr); 82070034214SMatthew G. Knepley ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr); 82170034214SMatthew G. Knepley if (flg) { 82270034214SMatthew G. Knepley ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr); 82370034214SMatthew G. Knepley ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 82470034214SMatthew G. Knepley ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr); 82570034214SMatthew G. Knepley ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 82670034214SMatthew G. Knepley ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr); 82770034214SMatthew G. Knepley } 82870034214SMatthew G. Knepley ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr); 82970034214SMatthew G. Knepley ierr = DMPlexGetConeOrientations(*dmParallel, &newCones);CHKERRQ(ierr); 83070034214SMatthew G. Knepley ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 83170034214SMatthew G. Knepley ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 83270034214SMatthew G. Knepley ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr); 83370034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr); 83470034214SMatthew G. Knepley /* Create supports and stratify sieve */ 83570034214SMatthew G. Knepley { 83670034214SMatthew G. Knepley PetscInt pStart, pEnd; 83770034214SMatthew G. Knepley 83870034214SMatthew G. Knepley ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr); 83970034214SMatthew G. Knepley ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr); 84070034214SMatthew G. Knepley } 84170034214SMatthew G. Knepley ierr = DMPlexSymmetrize(*dmParallel);CHKERRQ(ierr); 84270034214SMatthew G. Knepley ierr = DMPlexStratify(*dmParallel);CHKERRQ(ierr); 84370034214SMatthew G. Knepley /* Create point SF for parallel mesh */ 84470034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr); 84570034214SMatthew G. Knepley { 84670034214SMatthew G. Knepley const PetscInt *leaves; 84770034214SMatthew G. Knepley PetscSFNode *remotePoints, *rowners, *lowners; 84870034214SMatthew G. Knepley PetscInt numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints; 84970034214SMatthew G. Knepley PetscInt pStart, pEnd; 85070034214SMatthew G. Knepley 85170034214SMatthew G. Knepley ierr = DMPlexGetChart(*dmParallel, &pStart, &pEnd);CHKERRQ(ierr); 85270034214SMatthew G. Knepley ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr); 85370034214SMatthew G. Knepley ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr); 85470034214SMatthew G. Knepley for (p=0; p<numRoots; p++) { 85570034214SMatthew G. Knepley rowners[p].rank = -1; 85670034214SMatthew G. Knepley rowners[p].index = -1; 85770034214SMatthew G. Knepley } 85870034214SMatthew G. Knepley if (origCellPart) { 85970034214SMatthew G. Knepley /* Make sure points in the original partition are not assigned to other procs */ 86070034214SMatthew G. Knepley const PetscInt *origPoints; 86170034214SMatthew G. Knepley 86270034214SMatthew G. Knepley ierr = DMPlexCreatePartitionClosure(dm, origCellPartSection, origCellPart, &origPartSection, &origPart);CHKERRQ(ierr); 86370034214SMatthew G. Knepley ierr = ISGetIndices(origPart, &origPoints);CHKERRQ(ierr); 86470034214SMatthew G. Knepley for (p = 0; p < numProcs; ++p) { 86570034214SMatthew G. Knepley PetscInt dof, off, d; 86670034214SMatthew G. Knepley 86770034214SMatthew G. Knepley ierr = PetscSectionGetDof(origPartSection, p, &dof);CHKERRQ(ierr); 86870034214SMatthew G. Knepley ierr = PetscSectionGetOffset(origPartSection, p, &off);CHKERRQ(ierr); 86970034214SMatthew G. Knepley for (d = off; d < off+dof; ++d) { 87070034214SMatthew G. Knepley rowners[origPoints[d]].rank = p; 87170034214SMatthew G. Knepley } 87270034214SMatthew G. Knepley } 87370034214SMatthew G. Knepley ierr = ISRestoreIndices(origPart, &origPoints);CHKERRQ(ierr); 87470034214SMatthew G. Knepley ierr = ISDestroy(&origPart);CHKERRQ(ierr); 87570034214SMatthew G. Knepley ierr = PetscSectionDestroy(&origPartSection);CHKERRQ(ierr); 87670034214SMatthew G. Knepley } 87770034214SMatthew G. Knepley ierr = ISDestroy(&origCellPart);CHKERRQ(ierr); 87870034214SMatthew G. Knepley ierr = PetscSectionDestroy(&origCellPartSection);CHKERRQ(ierr); 87970034214SMatthew G. Knepley 88070034214SMatthew G. Knepley ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 88170034214SMatthew G. Knepley ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 88270034214SMatthew G. Knepley for (p = 0; p < numLeaves; ++p) { 88370034214SMatthew G. Knepley if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */ 88470034214SMatthew G. Knepley lowners[p].rank = rank; 88570034214SMatthew G. Knepley lowners[p].index = leaves ? leaves[p] : p; 88670034214SMatthew G. Knepley } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */ 88770034214SMatthew G. Knepley lowners[p].rank = -2; 88870034214SMatthew G. Knepley lowners[p].index = -2; 88970034214SMatthew G. Knepley } 89070034214SMatthew G. Knepley } 89170034214SMatthew G. Knepley for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */ 89270034214SMatthew G. Knepley rowners[p].rank = -3; 89370034214SMatthew G. Knepley rowners[p].index = -3; 89470034214SMatthew G. Knepley } 89570034214SMatthew G. Knepley ierr = PetscSFReduceBegin(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr); 89670034214SMatthew G. Knepley ierr = PetscSFReduceEnd(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr); 89770034214SMatthew G. Knepley ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 89870034214SMatthew G. Knepley ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 89970034214SMatthew G. Knepley for (p = 0; p < numLeaves; ++p) { 90070034214SMatthew G. Knepley if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed"); 90170034214SMatthew G. Knepley if (lowners[p].rank != rank) ++numGhostPoints; 90270034214SMatthew G. Knepley } 90370034214SMatthew G. Knepley ierr = PetscMalloc1(numGhostPoints, &ghostPoints);CHKERRQ(ierr); 90470034214SMatthew G. Knepley ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr); 90570034214SMatthew G. Knepley for (p = 0, gp = 0; p < numLeaves; ++p) { 90670034214SMatthew G. Knepley if (lowners[p].rank != rank) { 90770034214SMatthew G. Knepley ghostPoints[gp] = leaves ? leaves[p] : p; 90870034214SMatthew G. Knepley remotePoints[gp].rank = lowners[p].rank; 90970034214SMatthew G. Knepley remotePoints[gp].index = lowners[p].index; 91070034214SMatthew G. Knepley ++gp; 91170034214SMatthew G. Knepley } 91270034214SMatthew G. Knepley } 91370034214SMatthew G. Knepley ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr); 91470034214SMatthew G. Knepley ierr = PetscSFSetGraph((*dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 91570034214SMatthew G. Knepley ierr = PetscSFSetFromOptions((*dmParallel)->sf);CHKERRQ(ierr); 91670034214SMatthew G. Knepley } 917a42b08eeSMatthew G. Knepley pmesh->useCone = mesh->useCone; 918a42b08eeSMatthew G. Knepley pmesh->useClosure = mesh->useClosure; 91970034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr); 920bf5244c7SMatthew G. Knepley /* Distribute Coordinates */ 921bf5244c7SMatthew G. Knepley { 922bf5244c7SMatthew G. Knepley PetscSection originalCoordSection, newCoordSection; 923bf5244c7SMatthew G. Knepley Vec originalCoordinates, newCoordinates; 924bf5244c7SMatthew G. Knepley PetscInt bs; 925bf5244c7SMatthew G. Knepley const char *name; 926bf5244c7SMatthew G. Knepley const PetscReal *maxCell, *L; 927bf5244c7SMatthew G. Knepley 928bf5244c7SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr); 929bf5244c7SMatthew G. Knepley ierr = DMGetCoordinateSection(*dmParallel, &newCoordSection);CHKERRQ(ierr); 930bf5244c7SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr); 9311c57f3caSMatthew G. Knepley if (originalCoordinates) { 932bf5244c7SMatthew G. Knepley ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr); 933bf5244c7SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr); 934bf5244c7SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr); 935bf5244c7SMatthew G. Knepley 936bf5244c7SMatthew G. Knepley ierr = DMPlexDistributeField(dm, pointSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr); 937bf5244c7SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dmParallel, newCoordinates);CHKERRQ(ierr); 938bf5244c7SMatthew G. Knepley ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr); 939bf5244c7SMatthew G. Knepley ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr); 940bf5244c7SMatthew G. Knepley ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr); 9411c57f3caSMatthew G. Knepley } 942bf5244c7SMatthew G. Knepley ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr); 943bf5244c7SMatthew G. Knepley if (L) {ierr = DMSetPeriodicity(*dmParallel, maxCell, L);CHKERRQ(ierr);} 944bf5244c7SMatthew G. Knepley } 945bf5244c7SMatthew G. Knepley /* Distribute labels */ 946bf5244c7SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr); 947bf5244c7SMatthew G. Knepley { 948bf5244c7SMatthew G. Knepley DMLabel next = mesh->labels, newNext = pmesh->labels; 949bf5244c7SMatthew G. Knepley PetscInt numLabels = 0, l; 950bf5244c7SMatthew G. Knepley 951bf5244c7SMatthew G. Knepley /* Bcast number of labels */ 952bf5244c7SMatthew G. Knepley while (next) {++numLabels; next = next->next;} 953bf5244c7SMatthew G. Knepley ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 954bf5244c7SMatthew G. Knepley next = mesh->labels; 955bf5244c7SMatthew G. Knepley for (l = 0; l < numLabels; ++l) { 956bf5244c7SMatthew G. Knepley DMLabel labelNew; 957bf5244c7SMatthew G. Knepley PetscBool isdepth; 958bf5244c7SMatthew G. Knepley 959bf5244c7SMatthew G. Knepley /* Skip "depth" because it is recreated */ 960bf5244c7SMatthew G. Knepley if (!rank) {ierr = PetscStrcmp(next->name, "depth", &isdepth);CHKERRQ(ierr);} 961bf5244c7SMatthew G. Knepley ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr); 962bf5244c7SMatthew G. Knepley if (isdepth) {if (!rank) next = next->next; continue;} 963bf5244c7SMatthew G. Knepley ierr = DMLabelDistribute(next, partSection, part, renumbering, &labelNew);CHKERRQ(ierr); 964bf5244c7SMatthew G. Knepley /* Insert into list */ 965bf5244c7SMatthew G. Knepley if (newNext) newNext->next = labelNew; 966bf5244c7SMatthew G. Knepley else pmesh->labels = labelNew; 967bf5244c7SMatthew G. Knepley newNext = labelNew; 968bf5244c7SMatthew G. Knepley if (!rank) next = next->next; 969bf5244c7SMatthew G. Knepley } 970bf5244c7SMatthew G. Knepley } 971bf5244c7SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr); 972bf5244c7SMatthew G. Knepley /* Setup hybrid structure */ 973bf5244c7SMatthew G. Knepley { 974bf5244c7SMatthew G. Knepley const PetscInt *gpoints; 975bf5244c7SMatthew G. Knepley PetscInt depth, n, d; 976bf5244c7SMatthew G. Knepley 977bf5244c7SMatthew G. Knepley for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];} 978bf5244c7SMatthew G. Knepley ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr); 979bf5244c7SMatthew G. Knepley ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr); 980bf5244c7SMatthew G. Knepley ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr); 981bf5244c7SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 982bf5244c7SMatthew G. Knepley for (d = 0; d <= dim; ++d) { 983bf5244c7SMatthew G. Knepley PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p; 984bf5244c7SMatthew G. Knepley 985bf5244c7SMatthew G. Knepley if (pmax < 0) continue; 986bf5244c7SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr); 987bf5244c7SMatthew G. Knepley ierr = DMPlexGetDepthStratum(*dmParallel, d, NULL, &pEnd);CHKERRQ(ierr); 988bf5244c7SMatthew G. Knepley ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr); 989bf5244c7SMatthew G. Knepley for (p = 0; p < n; ++p) { 990bf5244c7SMatthew G. Knepley const PetscInt point = gpoints[p]; 991bf5244c7SMatthew G. Knepley 992bf5244c7SMatthew G. Knepley if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax; 993bf5244c7SMatthew G. Knepley } 994bf5244c7SMatthew G. Knepley if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax; 995bf5244c7SMatthew G. Knepley else pmesh->hybridPointMax[d] = -1; 996bf5244c7SMatthew G. Knepley } 997bf5244c7SMatthew G. Knepley ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr); 998bf5244c7SMatthew G. Knepley } 999bf5244c7SMatthew G. Knepley /* Cleanup Partition */ 1000bf5244c7SMatthew G. Knepley ierr = ISLocalToGlobalMappingDestroy(&renumbering);CHKERRQ(ierr); 1001bf5244c7SMatthew G. Knepley ierr = PetscSFDestroy(&partSF);CHKERRQ(ierr); 1002bf5244c7SMatthew G. Knepley ierr = PetscSectionDestroy(&partSection);CHKERRQ(ierr); 1003bf5244c7SMatthew G. Knepley ierr = ISDestroy(&part);CHKERRQ(ierr); 100486719b60SMatthew G. Knepley /* Copy BC */ 100586719b60SMatthew G. Knepley ierr = DMPlexCopyBoundary(dm, *dmParallel);CHKERRQ(ierr); 100670034214SMatthew G. Knepley /* Cleanup */ 100770034214SMatthew G. Knepley if (sf) {*sf = pointSF;} 100870034214SMatthew G. Knepley else {ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr);} 100970034214SMatthew G. Knepley ierr = DMSetFromOptions(*dmParallel);CHKERRQ(ierr); 101070034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr); 101170034214SMatthew G. Knepley PetscFunctionReturn(0); 101270034214SMatthew G. Knepley } 1013