1*70034214SMatthew G. Knepley #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2*70034214SMatthew G. Knepley #include <petscsf.h> 3*70034214SMatthew G. Knepley 4*70034214SMatthew G. Knepley #undef __FUNCT__ 5*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexSetAdjacencyUseCone" 6*70034214SMatthew G. Knepley /*@ 7*70034214SMatthew G. Knepley DMPlexSetAdjacencyUseCone - Define adjacency in the mesh using either the cone or the support first 8*70034214SMatthew G. Knepley 9*70034214SMatthew G. Knepley Input Parameters: 10*70034214SMatthew G. Knepley + dm - The DM object 11*70034214SMatthew G. Knepley - useCone - Flag to use the cone first 12*70034214SMatthew G. Knepley 13*70034214SMatthew G. Knepley Level: intermediate 14*70034214SMatthew G. Knepley 15*70034214SMatthew G. Knepley Notes: 16*70034214SMatthew G. Knepley $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 17*70034214SMatthew G. Knepley $ FVM: Two points p and q are adjacent if q \in star(cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 18*70034214SMatthew G. Knepley $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 19*70034214SMatthew G. Knepley 20*70034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator() 21*70034214SMatthew G. Knepley @*/ 22*70034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseCone(DM dm, PetscBool useCone) 23*70034214SMatthew G. Knepley { 24*70034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 25*70034214SMatthew G. Knepley 26*70034214SMatthew G. Knepley PetscFunctionBegin; 27*70034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 28*70034214SMatthew G. Knepley mesh->useCone = useCone; 29*70034214SMatthew G. Knepley PetscFunctionReturn(0); 30*70034214SMatthew G. Knepley } 31*70034214SMatthew G. Knepley 32*70034214SMatthew G. Knepley #undef __FUNCT__ 33*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacencyUseCone" 34*70034214SMatthew G. Knepley /*@ 35*70034214SMatthew G. Knepley DMPlexGetAdjacencyUseCone - Query whether adjacency in the mesh uses the cone or the support first 36*70034214SMatthew G. Knepley 37*70034214SMatthew G. Knepley Input Parameter: 38*70034214SMatthew G. Knepley . dm - The DM object 39*70034214SMatthew G. Knepley 40*70034214SMatthew G. Knepley Output Parameter: 41*70034214SMatthew G. Knepley . useCone - Flag to use the cone first 42*70034214SMatthew G. Knepley 43*70034214SMatthew G. Knepley Level: intermediate 44*70034214SMatthew G. Knepley 45*70034214SMatthew G. Knepley Notes: 46*70034214SMatthew G. Knepley $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 47*70034214SMatthew G. Knepley $ FVM: Two points p and q are adjacent if q \in star(cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 48*70034214SMatthew G. Knepley $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 49*70034214SMatthew G. Knepley 50*70034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator() 51*70034214SMatthew G. Knepley @*/ 52*70034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseCone(DM dm, PetscBool *useCone) 53*70034214SMatthew G. Knepley { 54*70034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 55*70034214SMatthew G. Knepley 56*70034214SMatthew G. Knepley PetscFunctionBegin; 57*70034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 58*70034214SMatthew G. Knepley PetscValidIntPointer(useCone, 2); 59*70034214SMatthew G. Knepley *useCone = mesh->useCone; 60*70034214SMatthew G. Knepley PetscFunctionReturn(0); 61*70034214SMatthew G. Knepley } 62*70034214SMatthew G. Knepley 63*70034214SMatthew G. Knepley #undef __FUNCT__ 64*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexSetAdjacencyUseClosure" 65*70034214SMatthew G. Knepley /*@ 66*70034214SMatthew G. Knepley DMPlexSetAdjacencyUseClosure - Define adjacency in the mesh using the transitive closure 67*70034214SMatthew G. Knepley 68*70034214SMatthew G. Knepley Input Parameters: 69*70034214SMatthew G. Knepley + dm - The DM object 70*70034214SMatthew G. Knepley - useClosure - Flag to use the closure 71*70034214SMatthew G. Knepley 72*70034214SMatthew G. Knepley Level: intermediate 73*70034214SMatthew G. Knepley 74*70034214SMatthew G. Knepley Notes: 75*70034214SMatthew G. Knepley $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 76*70034214SMatthew G. Knepley $ FVM: Two points p and q are adjacent if q \in star(cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 77*70034214SMatthew G. Knepley $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 78*70034214SMatthew G. Knepley 79*70034214SMatthew G. Knepley .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator() 80*70034214SMatthew G. Knepley @*/ 81*70034214SMatthew G. Knepley PetscErrorCode DMPlexSetAdjacencyUseClosure(DM dm, PetscBool useClosure) 82*70034214SMatthew G. Knepley { 83*70034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 84*70034214SMatthew G. Knepley 85*70034214SMatthew G. Knepley PetscFunctionBegin; 86*70034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 87*70034214SMatthew G. Knepley mesh->useClosure = useClosure; 88*70034214SMatthew G. Knepley PetscFunctionReturn(0); 89*70034214SMatthew G. Knepley } 90*70034214SMatthew G. Knepley 91*70034214SMatthew G. Knepley #undef __FUNCT__ 92*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacencyUseClosure" 93*70034214SMatthew G. Knepley /*@ 94*70034214SMatthew G. Knepley DMPlexGetAdjacencyUseClosure - Query whether adjacency in the mesh uses the transitive closure 95*70034214SMatthew G. Knepley 96*70034214SMatthew G. Knepley Input Parameter: 97*70034214SMatthew G. Knepley . dm - The DM object 98*70034214SMatthew G. Knepley 99*70034214SMatthew G. Knepley Output Parameter: 100*70034214SMatthew G. Knepley . useClosure - Flag to use the closure 101*70034214SMatthew G. Knepley 102*70034214SMatthew G. Knepley Level: intermediate 103*70034214SMatthew G. Knepley 104*70034214SMatthew G. Knepley Notes: 105*70034214SMatthew G. Knepley $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 106*70034214SMatthew G. Knepley $ FVM: Two points p and q are adjacent if q \in star(cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 107*70034214SMatthew G. Knepley $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 108*70034214SMatthew G. Knepley 109*70034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator() 110*70034214SMatthew G. Knepley @*/ 111*70034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacencyUseClosure(DM dm, PetscBool *useClosure) 112*70034214SMatthew G. Knepley { 113*70034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 114*70034214SMatthew G. Knepley 115*70034214SMatthew G. Knepley PetscFunctionBegin; 116*70034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 117*70034214SMatthew G. Knepley PetscValidIntPointer(useClosure, 2); 118*70034214SMatthew G. Knepley *useClosure = mesh->useClosure; 119*70034214SMatthew G. Knepley PetscFunctionReturn(0); 120*70034214SMatthew G. Knepley } 121*70034214SMatthew G. Knepley 122*70034214SMatthew G. Knepley #undef __FUNCT__ 123*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Cone_Internal" 124*70034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Cone_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[]) 125*70034214SMatthew G. Knepley { 126*70034214SMatthew G. Knepley const PetscInt *cone = NULL; 127*70034214SMatthew G. Knepley PetscInt numAdj = 0, maxAdjSize = *adjSize, coneSize, c; 128*70034214SMatthew G. Knepley PetscErrorCode ierr; 129*70034214SMatthew G. Knepley 130*70034214SMatthew G. Knepley PetscFunctionBeginHot; 131*70034214SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 132*70034214SMatthew G. Knepley ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 133*70034214SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 134*70034214SMatthew G. Knepley const PetscInt *support = NULL; 135*70034214SMatthew G. Knepley PetscInt supportSize, s, q; 136*70034214SMatthew G. Knepley 137*70034214SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, cone[c], &supportSize);CHKERRQ(ierr); 138*70034214SMatthew G. Knepley ierr = DMPlexGetSupport(dm, cone[c], &support);CHKERRQ(ierr); 139*70034214SMatthew G. Knepley for (s = 0; s < supportSize; ++s) { 140*70034214SMatthew G. Knepley for (q = 0; q < numAdj || (adj[numAdj++] = support[s],0); ++q) { 141*70034214SMatthew G. Knepley if (support[s] == adj[q]) break; 142*70034214SMatthew G. Knepley } 143*70034214SMatthew G. Knepley if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize); 144*70034214SMatthew G. Knepley } 145*70034214SMatthew G. Knepley } 146*70034214SMatthew G. Knepley *adjSize = numAdj; 147*70034214SMatthew G. Knepley PetscFunctionReturn(0); 148*70034214SMatthew G. Knepley } 149*70034214SMatthew G. Knepley 150*70034214SMatthew G. Knepley #undef __FUNCT__ 151*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Support_Internal" 152*70034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Support_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[]) 153*70034214SMatthew G. Knepley { 154*70034214SMatthew G. Knepley const PetscInt *support = NULL; 155*70034214SMatthew G. Knepley PetscInt numAdj = 0, maxAdjSize = *adjSize, supportSize, s; 156*70034214SMatthew G. Knepley PetscErrorCode ierr; 157*70034214SMatthew G. Knepley 158*70034214SMatthew G. Knepley PetscFunctionBeginHot; 159*70034214SMatthew G. Knepley ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr); 160*70034214SMatthew G. Knepley ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr); 161*70034214SMatthew G. Knepley for (s = 0; s < supportSize; ++s) { 162*70034214SMatthew G. Knepley const PetscInt *cone = NULL; 163*70034214SMatthew G. Knepley PetscInt coneSize, c, q; 164*70034214SMatthew G. Knepley 165*70034214SMatthew G. Knepley ierr = DMPlexGetConeSize(dm, support[s], &coneSize);CHKERRQ(ierr); 166*70034214SMatthew G. Knepley ierr = DMPlexGetCone(dm, support[s], &cone);CHKERRQ(ierr); 167*70034214SMatthew G. Knepley for (c = 0; c < coneSize; ++c) { 168*70034214SMatthew G. Knepley for (q = 0; q < numAdj || (adj[numAdj++] = cone[c],0); ++q) { 169*70034214SMatthew G. Knepley if (cone[c] == adj[q]) break; 170*70034214SMatthew G. Knepley } 171*70034214SMatthew G. Knepley if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize); 172*70034214SMatthew G. Knepley } 173*70034214SMatthew G. Knepley } 174*70034214SMatthew G. Knepley *adjSize = numAdj; 175*70034214SMatthew G. Knepley PetscFunctionReturn(0); 176*70034214SMatthew G. Knepley } 177*70034214SMatthew G. Knepley 178*70034214SMatthew G. Knepley #undef __FUNCT__ 179*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Transitive_Internal" 180*70034214SMatthew G. Knepley static PetscErrorCode DMPlexGetAdjacency_Transitive_Internal(DM dm, PetscInt p, PetscBool useClosure, PetscInt *adjSize, PetscInt adj[]) 181*70034214SMatthew G. Knepley { 182*70034214SMatthew G. Knepley PetscInt *star = NULL; 183*70034214SMatthew G. Knepley PetscInt numAdj = 0, maxAdjSize = *adjSize, starSize, s; 184*70034214SMatthew G. Knepley PetscErrorCode ierr; 185*70034214SMatthew G. Knepley 186*70034214SMatthew G. Knepley PetscFunctionBeginHot; 187*70034214SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr); 188*70034214SMatthew G. Knepley for (s = 0; s < starSize*2; s += 2) { 189*70034214SMatthew G. Knepley const PetscInt *closure = NULL; 190*70034214SMatthew G. Knepley PetscInt closureSize, c, q; 191*70034214SMatthew G. Knepley 192*70034214SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr); 193*70034214SMatthew G. Knepley for (c = 0; c < closureSize*2; c += 2) { 194*70034214SMatthew G. Knepley for (q = 0; q < numAdj || (adj[numAdj++] = closure[c],0); ++q) { 195*70034214SMatthew G. Knepley if (closure[c] == adj[q]) break; 196*70034214SMatthew G. Knepley } 197*70034214SMatthew G. Knepley if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize); 198*70034214SMatthew G. Knepley } 199*70034214SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr); 200*70034214SMatthew G. Knepley } 201*70034214SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr); 202*70034214SMatthew G. Knepley *adjSize = numAdj; 203*70034214SMatthew G. Knepley PetscFunctionReturn(0); 204*70034214SMatthew G. Knepley } 205*70034214SMatthew G. Knepley 206*70034214SMatthew G. Knepley #undef __FUNCT__ 207*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency_Internal" 208*70034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacency_Internal(DM dm, PetscInt p, PetscBool useCone, PetscBool useTransitiveClosure, PetscInt *adjSize, PetscInt adj[]) 209*70034214SMatthew G. Knepley { 210*70034214SMatthew G. Knepley PetscErrorCode ierr; 211*70034214SMatthew G. Knepley 212*70034214SMatthew G. Knepley PetscFunctionBeginHot; 213*70034214SMatthew G. Knepley if (useTransitiveClosure) { 214*70034214SMatthew G. Knepley ierr = DMPlexGetAdjacency_Transitive_Internal(dm, p, useCone, adjSize, adj);CHKERRQ(ierr); 215*70034214SMatthew G. Knepley } else if (useCone) { 216*70034214SMatthew G. Knepley ierr = DMPlexGetAdjacency_Cone_Internal(dm, p, adjSize, adj);CHKERRQ(ierr); 217*70034214SMatthew G. Knepley } else { 218*70034214SMatthew G. Knepley ierr = DMPlexGetAdjacency_Support_Internal(dm, p, adjSize, adj);CHKERRQ(ierr); 219*70034214SMatthew G. Knepley } 220*70034214SMatthew G. Knepley PetscFunctionReturn(0); 221*70034214SMatthew G. Knepley } 222*70034214SMatthew G. Knepley 223*70034214SMatthew G. Knepley #undef __FUNCT__ 224*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexGetAdjacency" 225*70034214SMatthew G. Knepley /*@ 226*70034214SMatthew G. Knepley DMPlexGetAdjacency - Return all points adjacent to the given point 227*70034214SMatthew G. Knepley 228*70034214SMatthew G. Knepley Input Parameters: 229*70034214SMatthew G. Knepley + dm - The DM object 230*70034214SMatthew G. Knepley . p - The point 231*70034214SMatthew G. Knepley . adjSize - The maximum size of adj if it is non-NULL, or PETSC_DETERMINE 232*70034214SMatthew G. Knepley - adj - Either NULL so that the array is allocated, or an existing array with size adjSize 233*70034214SMatthew G. Knepley 234*70034214SMatthew G. Knepley Output Parameters: 235*70034214SMatthew G. Knepley + adjSize - The number of adjacent points 236*70034214SMatthew G. Knepley - adj - The adjacent points 237*70034214SMatthew G. Knepley 238*70034214SMatthew G. Knepley Level: advanced 239*70034214SMatthew G. Knepley 240*70034214SMatthew G. Knepley Notes: The user must PetscFree the adj array if it was not passed in. 241*70034214SMatthew G. Knepley 242*70034214SMatthew G. Knepley .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexDistribute(), DMCreateMatrix(), DMPlexPreallocateOperator() 243*70034214SMatthew G. Knepley @*/ 244*70034214SMatthew G. Knepley PetscErrorCode DMPlexGetAdjacency(DM dm, PetscInt p, PetscInt *adjSize, PetscInt *adj[]) 245*70034214SMatthew G. Knepley { 246*70034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex *) dm->data; 247*70034214SMatthew G. Knepley static PetscInt asiz = 0; 248*70034214SMatthew G. Knepley PetscErrorCode ierr; 249*70034214SMatthew G. Knepley 250*70034214SMatthew G. Knepley PetscFunctionBeginHot; 251*70034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 252*70034214SMatthew G. Knepley PetscValidPointer(adjSize,3); 253*70034214SMatthew G. Knepley PetscValidPointer(adj,4); 254*70034214SMatthew G. Knepley if (!*adj) { 255*70034214SMatthew G. Knepley PetscInt depth, maxConeSize, maxSupportSize; 256*70034214SMatthew G. Knepley 257*70034214SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 258*70034214SMatthew G. Knepley ierr = DMPlexGetMaxSizes(dm, &maxConeSize, &maxSupportSize);CHKERRQ(ierr); 259*70034214SMatthew G. Knepley asiz = PetscPowInt(maxConeSize, depth+1) * PetscPowInt(maxSupportSize, depth+1) + 1; 260*70034214SMatthew G. Knepley ierr = PetscMalloc1(asiz,adj);CHKERRQ(ierr); 261*70034214SMatthew G. Knepley } 262*70034214SMatthew G. Knepley if (*adjSize < 0) *adjSize = asiz; 263*70034214SMatthew G. Knepley ierr = DMPlexGetAdjacency_Internal(dm, p, mesh->useCone, mesh->useClosure, adjSize, *adj);CHKERRQ(ierr); 264*70034214SMatthew G. Knepley PetscFunctionReturn(0); 265*70034214SMatthew G. Knepley } 266*70034214SMatthew G. Knepley 267*70034214SMatthew G. Knepley #undef __FUNCT__ 268*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeField" 269*70034214SMatthew G. Knepley /*@ 270*70034214SMatthew G. Knepley DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution 271*70034214SMatthew G. Knepley 272*70034214SMatthew G. Knepley Collective on DM 273*70034214SMatthew G. Knepley 274*70034214SMatthew G. Knepley Input Parameters: 275*70034214SMatthew G. Knepley + dm - The DMPlex object 276*70034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern 277*70034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout 278*70034214SMatthew G. Knepley - originalVec - The existing data 279*70034214SMatthew G. Knepley 280*70034214SMatthew G. Knepley Output Parameters: 281*70034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout 282*70034214SMatthew G. Knepley - newVec - The new data 283*70034214SMatthew G. Knepley 284*70034214SMatthew G. Knepley Level: developer 285*70034214SMatthew G. Knepley 286*70034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeData() 287*70034214SMatthew G. Knepley @*/ 288*70034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec) 289*70034214SMatthew G. Knepley { 290*70034214SMatthew G. Knepley PetscSF fieldSF; 291*70034214SMatthew G. Knepley PetscInt *remoteOffsets, fieldSize; 292*70034214SMatthew G. Knepley PetscScalar *originalValues, *newValues; 293*70034214SMatthew G. Knepley PetscErrorCode ierr; 294*70034214SMatthew G. Knepley 295*70034214SMatthew G. Knepley PetscFunctionBegin; 296*70034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 297*70034214SMatthew G. Knepley ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); 298*70034214SMatthew G. Knepley 299*70034214SMatthew G. Knepley ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); 300*70034214SMatthew G. Knepley ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr); 301*70034214SMatthew G. Knepley ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr); 302*70034214SMatthew G. Knepley 303*70034214SMatthew G. Knepley ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr); 304*70034214SMatthew G. Knepley ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr); 305*70034214SMatthew G. Knepley ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); 306*70034214SMatthew G. Knepley ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr); 307*70034214SMatthew G. Knepley ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr); 308*70034214SMatthew G. Knepley ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); 309*70034214SMatthew G. Knepley ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr); 310*70034214SMatthew G. Knepley ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr); 311*70034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 312*70034214SMatthew G. Knepley PetscFunctionReturn(0); 313*70034214SMatthew G. Knepley } 314*70034214SMatthew G. Knepley 315*70034214SMatthew G. Knepley #undef __FUNCT__ 316*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistributeData" 317*70034214SMatthew G. Knepley /*@ 318*70034214SMatthew G. Knepley DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution 319*70034214SMatthew G. Knepley 320*70034214SMatthew G. Knepley Collective on DM 321*70034214SMatthew G. Knepley 322*70034214SMatthew G. Knepley Input Parameters: 323*70034214SMatthew G. Knepley + dm - The DMPlex object 324*70034214SMatthew G. Knepley . pointSF - The PetscSF describing the communication pattern 325*70034214SMatthew G. Knepley . originalSection - The PetscSection for existing data layout 326*70034214SMatthew G. Knepley . datatype - The type of data 327*70034214SMatthew G. Knepley - originalData - The existing data 328*70034214SMatthew G. Knepley 329*70034214SMatthew G. Knepley Output Parameters: 330*70034214SMatthew G. Knepley + newSection - The PetscSF describing the new data layout 331*70034214SMatthew G. Knepley - newData - The new data 332*70034214SMatthew G. Knepley 333*70034214SMatthew G. Knepley Level: developer 334*70034214SMatthew G. Knepley 335*70034214SMatthew G. Knepley .seealso: DMPlexDistribute(), DMPlexDistributeField() 336*70034214SMatthew G. Knepley @*/ 337*70034214SMatthew G. Knepley PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData) 338*70034214SMatthew G. Knepley { 339*70034214SMatthew G. Knepley PetscSF fieldSF; 340*70034214SMatthew G. Knepley PetscInt *remoteOffsets, fieldSize; 341*70034214SMatthew G. Knepley PetscMPIInt dataSize; 342*70034214SMatthew G. Knepley PetscErrorCode ierr; 343*70034214SMatthew G. Knepley 344*70034214SMatthew G. Knepley PetscFunctionBegin; 345*70034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr); 346*70034214SMatthew G. Knepley ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); 347*70034214SMatthew G. Knepley 348*70034214SMatthew G. Knepley ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); 349*70034214SMatthew G. Knepley ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr); 350*70034214SMatthew G. Knepley ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr); 351*70034214SMatthew G. Knepley 352*70034214SMatthew G. Knepley ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); 353*70034214SMatthew G. Knepley ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr); 354*70034214SMatthew G. Knepley ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr); 355*70034214SMatthew G. Knepley ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); 356*70034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr); 357*70034214SMatthew G. Knepley PetscFunctionReturn(0); 358*70034214SMatthew G. Knepley } 359*70034214SMatthew G. Knepley 360*70034214SMatthew G. Knepley #undef __FUNCT__ 361*70034214SMatthew G. Knepley #define __FUNCT__ "DMPlexDistribute" 362*70034214SMatthew G. Knepley /*@C 363*70034214SMatthew G. Knepley DMPlexDistribute - Distributes the mesh and any associated sections. 364*70034214SMatthew G. Knepley 365*70034214SMatthew G. Knepley Not Collective 366*70034214SMatthew G. Knepley 367*70034214SMatthew G. Knepley Input Parameter: 368*70034214SMatthew G. Knepley + dm - The original DMPlex object 369*70034214SMatthew G. Knepley . partitioner - The partitioning package, or NULL for the default 370*70034214SMatthew G. Knepley - overlap - The overlap of partitions, 0 is the default 371*70034214SMatthew G. Knepley 372*70034214SMatthew G. Knepley Output Parameter: 373*70034214SMatthew G. Knepley + sf - The PetscSF used for point distribution 374*70034214SMatthew G. Knepley - parallelMesh - The distributed DMPlex object, or NULL 375*70034214SMatthew G. Knepley 376*70034214SMatthew G. Knepley Note: If the mesh was not distributed, the return value is NULL 377*70034214SMatthew G. Knepley 378*70034214SMatthew G. Knepley Level: intermediate 379*70034214SMatthew G. Knepley 380*70034214SMatthew G. Knepley .keywords: mesh, elements 381*70034214SMatthew G. Knepley .seealso: DMPlexCreate(), DMPlexDistributeByFace() 382*70034214SMatthew G. Knepley @*/ 383*70034214SMatthew G. Knepley PetscErrorCode DMPlexDistribute(DM dm, const char partitioner[], PetscInt overlap, PetscSF *sf, DM *dmParallel) 384*70034214SMatthew G. Knepley { 385*70034214SMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data, *pmesh; 386*70034214SMatthew G. Knepley MPI_Comm comm; 387*70034214SMatthew G. Knepley const PetscInt height = 0; 388*70034214SMatthew G. Knepley PetscInt dim, numRemoteRanks; 389*70034214SMatthew G. Knepley IS origCellPart, origPart, cellPart, part; 390*70034214SMatthew G. Knepley PetscSection origCellPartSection, origPartSection, cellPartSection, partSection; 391*70034214SMatthew G. Knepley PetscSFNode *remoteRanks; 392*70034214SMatthew G. Knepley PetscSF partSF, pointSF, coneSF; 393*70034214SMatthew G. Knepley ISLocalToGlobalMapping renumbering; 394*70034214SMatthew G. Knepley PetscSection originalConeSection, newConeSection; 395*70034214SMatthew G. Knepley PetscInt *remoteOffsets; 396*70034214SMatthew G. Knepley PetscInt *cones, *newCones, newConesSize; 397*70034214SMatthew G. Knepley PetscBool flg; 398*70034214SMatthew G. Knepley PetscMPIInt rank, numProcs, p; 399*70034214SMatthew G. Knepley PetscErrorCode ierr; 400*70034214SMatthew G. Knepley 401*70034214SMatthew G. Knepley PetscFunctionBegin; 402*70034214SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 403*70034214SMatthew G. Knepley if (sf) PetscValidPointer(sf,4); 404*70034214SMatthew G. Knepley PetscValidPointer(dmParallel,5); 405*70034214SMatthew G. Knepley 406*70034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr); 407*70034214SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 408*70034214SMatthew G. Knepley ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 409*70034214SMatthew G. Knepley ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 410*70034214SMatthew G. Knepley 411*70034214SMatthew G. Knepley *dmParallel = NULL; 412*70034214SMatthew G. Knepley if (numProcs == 1) PetscFunctionReturn(0); 413*70034214SMatthew G. Knepley 414*70034214SMatthew G. Knepley ierr = DMPlexGetDimension(dm, &dim);CHKERRQ(ierr); 415*70034214SMatthew G. Knepley /* Create cell partition - We need to rewrite to use IS, use the MatPartition stuff */ 416*70034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_Partition,dm,0,0,0);CHKERRQ(ierr); 417*70034214SMatthew G. Knepley if (overlap > 1) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "Overlap > 1 not yet implemented"); 418*70034214SMatthew G. Knepley ierr = DMPlexCreatePartition(dm, partitioner, height, overlap > 0 ? PETSC_TRUE : PETSC_FALSE, &cellPartSection, &cellPart, &origCellPartSection, &origCellPart);CHKERRQ(ierr); 419*70034214SMatthew G. Knepley /* Create SF assuming a serial partition for all processes: Could check for IS length here */ 420*70034214SMatthew G. Knepley if (!rank) numRemoteRanks = numProcs; 421*70034214SMatthew G. Knepley else numRemoteRanks = 0; 422*70034214SMatthew G. Knepley ierr = PetscMalloc1(numRemoteRanks, &remoteRanks);CHKERRQ(ierr); 423*70034214SMatthew G. Knepley for (p = 0; p < numRemoteRanks; ++p) { 424*70034214SMatthew G. Knepley remoteRanks[p].rank = p; 425*70034214SMatthew G. Knepley remoteRanks[p].index = 0; 426*70034214SMatthew G. Knepley } 427*70034214SMatthew G. Knepley ierr = PetscSFCreate(comm, &partSF);CHKERRQ(ierr); 428*70034214SMatthew G. Knepley ierr = PetscSFSetGraph(partSF, 1, numRemoteRanks, NULL, PETSC_OWN_POINTER, remoteRanks, PETSC_OWN_POINTER);CHKERRQ(ierr); 429*70034214SMatthew G. Knepley ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr); 430*70034214SMatthew G. Knepley if (flg) { 431*70034214SMatthew G. Knepley ierr = PetscPrintf(comm, "Cell Partition:\n");CHKERRQ(ierr); 432*70034214SMatthew G. Knepley ierr = PetscSectionView(cellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 433*70034214SMatthew G. Knepley ierr = ISView(cellPart, NULL);CHKERRQ(ierr); 434*70034214SMatthew G. Knepley if (origCellPart) { 435*70034214SMatthew G. Knepley ierr = PetscPrintf(comm, "Original Cell Partition:\n");CHKERRQ(ierr); 436*70034214SMatthew G. Knepley ierr = PetscSectionView(origCellPartSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 437*70034214SMatthew G. Knepley ierr = ISView(origCellPart, NULL);CHKERRQ(ierr); 438*70034214SMatthew G. Knepley } 439*70034214SMatthew G. Knepley ierr = PetscSFView(partSF, NULL);CHKERRQ(ierr); 440*70034214SMatthew G. Knepley } 441*70034214SMatthew G. Knepley /* Close the partition over the mesh */ 442*70034214SMatthew G. Knepley ierr = DMPlexCreatePartitionClosure(dm, cellPartSection, cellPart, &partSection, &part);CHKERRQ(ierr); 443*70034214SMatthew G. Knepley ierr = ISDestroy(&cellPart);CHKERRQ(ierr); 444*70034214SMatthew G. Knepley ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr); 445*70034214SMatthew G. Knepley /* Create new mesh */ 446*70034214SMatthew G. Knepley ierr = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr); 447*70034214SMatthew G. Knepley ierr = DMPlexSetDimension(*dmParallel, dim);CHKERRQ(ierr); 448*70034214SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr); 449*70034214SMatthew G. Knepley pmesh = (DM_Plex*) (*dmParallel)->data; 450*70034214SMatthew G. Knepley /* Distribute sieve points and the global point numbering (replaces creating remote bases) */ 451*70034214SMatthew G. Knepley ierr = PetscSFConvertPartition(partSF, partSection, part, &renumbering, &pointSF);CHKERRQ(ierr); 452*70034214SMatthew G. Knepley if (flg) { 453*70034214SMatthew G. Knepley ierr = PetscPrintf(comm, "Point Partition:\n");CHKERRQ(ierr); 454*70034214SMatthew G. Knepley ierr = PetscSectionView(partSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 455*70034214SMatthew G. Knepley ierr = ISView(part, NULL);CHKERRQ(ierr); 456*70034214SMatthew G. Knepley ierr = PetscSFView(pointSF, NULL);CHKERRQ(ierr); 457*70034214SMatthew G. Knepley ierr = PetscPrintf(comm, "Point Renumbering after partition:\n");CHKERRQ(ierr); 458*70034214SMatthew G. Knepley ierr = ISLocalToGlobalMappingView(renumbering, NULL);CHKERRQ(ierr); 459*70034214SMatthew G. Knepley } 460*70034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_Partition,dm,0,0,0);CHKERRQ(ierr); 461*70034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr); 462*70034214SMatthew G. Knepley /* Distribute cone section */ 463*70034214SMatthew G. Knepley ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr); 464*70034214SMatthew G. Knepley ierr = DMPlexGetConeSection(*dmParallel, &newConeSection);CHKERRQ(ierr); 465*70034214SMatthew G. Knepley ierr = PetscSFDistributeSection(pointSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr); 466*70034214SMatthew G. Knepley ierr = DMSetUp(*dmParallel);CHKERRQ(ierr); 467*70034214SMatthew G. Knepley { 468*70034214SMatthew G. Knepley PetscInt pStart, pEnd, p; 469*70034214SMatthew G. Knepley 470*70034214SMatthew G. Knepley ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr); 471*70034214SMatthew G. Knepley for (p = pStart; p < pEnd; ++p) { 472*70034214SMatthew G. Knepley PetscInt coneSize; 473*70034214SMatthew G. Knepley ierr = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr); 474*70034214SMatthew G. Knepley pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize); 475*70034214SMatthew G. Knepley } 476*70034214SMatthew G. Knepley } 477*70034214SMatthew G. Knepley /* Communicate and renumber cones */ 478*70034214SMatthew G. Knepley ierr = PetscSFCreateSectionSF(pointSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr); 479*70034214SMatthew G. Knepley ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr); 480*70034214SMatthew G. Knepley ierr = DMPlexGetCones(*dmParallel, &newCones);CHKERRQ(ierr); 481*70034214SMatthew G. Knepley ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 482*70034214SMatthew G. Knepley ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 483*70034214SMatthew G. Knepley ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr); 484*70034214SMatthew G. Knepley ierr = ISGlobalToLocalMappingApply(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr); 485*70034214SMatthew G. Knepley ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr); 486*70034214SMatthew G. Knepley if (flg) { 487*70034214SMatthew G. Knepley ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr); 488*70034214SMatthew G. Knepley ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 489*70034214SMatthew G. Knepley ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr); 490*70034214SMatthew G. Knepley ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 491*70034214SMatthew G. Knepley ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr); 492*70034214SMatthew G. Knepley } 493*70034214SMatthew G. Knepley ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr); 494*70034214SMatthew G. Knepley ierr = DMPlexGetConeOrientations(*dmParallel, &newCones);CHKERRQ(ierr); 495*70034214SMatthew G. Knepley ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 496*70034214SMatthew G. Knepley ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 497*70034214SMatthew G. Knepley ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr); 498*70034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr); 499*70034214SMatthew G. Knepley /* Create supports and stratify sieve */ 500*70034214SMatthew G. Knepley { 501*70034214SMatthew G. Knepley PetscInt pStart, pEnd; 502*70034214SMatthew G. Knepley 503*70034214SMatthew G. Knepley ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr); 504*70034214SMatthew G. Knepley ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr); 505*70034214SMatthew G. Knepley } 506*70034214SMatthew G. Knepley ierr = DMPlexSymmetrize(*dmParallel);CHKERRQ(ierr); 507*70034214SMatthew G. Knepley ierr = DMPlexStratify(*dmParallel);CHKERRQ(ierr); 508*70034214SMatthew G. Knepley /* Distribute Coordinates */ 509*70034214SMatthew G. Knepley { 510*70034214SMatthew G. Knepley PetscSection originalCoordSection, newCoordSection; 511*70034214SMatthew G. Knepley Vec originalCoordinates, newCoordinates; 512*70034214SMatthew G. Knepley const char *name; 513*70034214SMatthew G. Knepley 514*70034214SMatthew G. Knepley ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr); 515*70034214SMatthew G. Knepley ierr = DMGetCoordinateSection(*dmParallel, &newCoordSection);CHKERRQ(ierr); 516*70034214SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr); 517*70034214SMatthew G. Knepley ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr); 518*70034214SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr); 519*70034214SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr); 520*70034214SMatthew G. Knepley 521*70034214SMatthew G. Knepley ierr = DMPlexDistributeField(dm, pointSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr); 522*70034214SMatthew G. Knepley ierr = DMSetCoordinatesLocal(*dmParallel, newCoordinates);CHKERRQ(ierr); 523*70034214SMatthew G. Knepley ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr); 524*70034214SMatthew G. Knepley } 525*70034214SMatthew G. Knepley /* Distribute labels */ 526*70034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr); 527*70034214SMatthew G. Knepley { 528*70034214SMatthew G. Knepley DMLabel next = mesh->labels, newNext = pmesh->labels; 529*70034214SMatthew G. Knepley PetscInt numLabels = 0, l; 530*70034214SMatthew G. Knepley 531*70034214SMatthew G. Knepley /* Bcast number of labels */ 532*70034214SMatthew G. Knepley while (next) {++numLabels; next = next->next;} 533*70034214SMatthew G. Knepley ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 534*70034214SMatthew G. Knepley next = mesh->labels; 535*70034214SMatthew G. Knepley for (l = 0; l < numLabels; ++l) { 536*70034214SMatthew G. Knepley DMLabel labelNew; 537*70034214SMatthew G. Knepley PetscBool isdepth; 538*70034214SMatthew G. Knepley 539*70034214SMatthew G. Knepley /* Skip "depth" because it is recreated */ 540*70034214SMatthew G. Knepley if (!rank) {ierr = PetscStrcmp(next->name, "depth", &isdepth);CHKERRQ(ierr);} 541*70034214SMatthew G. Knepley ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr); 542*70034214SMatthew G. Knepley if (isdepth) {if (!rank) next = next->next; continue;} 543*70034214SMatthew G. Knepley ierr = DMLabelDistribute(next, partSection, part, renumbering, &labelNew);CHKERRQ(ierr); 544*70034214SMatthew G. Knepley /* Insert into list */ 545*70034214SMatthew G. Knepley if (newNext) newNext->next = labelNew; 546*70034214SMatthew G. Knepley else pmesh->labels = labelNew; 547*70034214SMatthew G. Knepley newNext = labelNew; 548*70034214SMatthew G. Knepley if (!rank) next = next->next; 549*70034214SMatthew G. Knepley } 550*70034214SMatthew G. Knepley } 551*70034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr); 552*70034214SMatthew G. Knepley /* Setup hybrid structure */ 553*70034214SMatthew G. Knepley { 554*70034214SMatthew G. Knepley const PetscInt *gpoints; 555*70034214SMatthew G. Knepley PetscInt depth, n, d; 556*70034214SMatthew G. Knepley 557*70034214SMatthew G. Knepley for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];} 558*70034214SMatthew G. Knepley ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr); 559*70034214SMatthew G. Knepley ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr); 560*70034214SMatthew G. Knepley ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr); 561*70034214SMatthew G. Knepley ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 562*70034214SMatthew G. Knepley for (d = 0; d <= dim; ++d) { 563*70034214SMatthew G. Knepley PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p; 564*70034214SMatthew G. Knepley 565*70034214SMatthew G. Knepley if (pmax < 0) continue; 566*70034214SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr); 567*70034214SMatthew G. Knepley ierr = DMPlexGetDepthStratum(*dmParallel, d, NULL, &pEnd);CHKERRQ(ierr); 568*70034214SMatthew G. Knepley ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr); 569*70034214SMatthew G. Knepley for (p = 0; p < n; ++p) { 570*70034214SMatthew G. Knepley const PetscInt point = gpoints[p]; 571*70034214SMatthew G. Knepley 572*70034214SMatthew G. Knepley if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax; 573*70034214SMatthew G. Knepley } 574*70034214SMatthew G. Knepley if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax; 575*70034214SMatthew G. Knepley else pmesh->hybridPointMax[d] = -1; 576*70034214SMatthew G. Knepley } 577*70034214SMatthew G. Knepley ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr); 578*70034214SMatthew G. Knepley } 579*70034214SMatthew G. Knepley /* Cleanup Partition */ 580*70034214SMatthew G. Knepley ierr = ISLocalToGlobalMappingDestroy(&renumbering);CHKERRQ(ierr); 581*70034214SMatthew G. Knepley ierr = PetscSFDestroy(&partSF);CHKERRQ(ierr); 582*70034214SMatthew G. Knepley ierr = PetscSectionDestroy(&partSection);CHKERRQ(ierr); 583*70034214SMatthew G. Knepley ierr = ISDestroy(&part);CHKERRQ(ierr); 584*70034214SMatthew G. Knepley /* Create point SF for parallel mesh */ 585*70034214SMatthew G. Knepley ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr); 586*70034214SMatthew G. Knepley { 587*70034214SMatthew G. Knepley const PetscInt *leaves; 588*70034214SMatthew G. Knepley PetscSFNode *remotePoints, *rowners, *lowners; 589*70034214SMatthew G. Knepley PetscInt numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints; 590*70034214SMatthew G. Knepley PetscInt pStart, pEnd; 591*70034214SMatthew G. Knepley 592*70034214SMatthew G. Knepley ierr = DMPlexGetChart(*dmParallel, &pStart, &pEnd);CHKERRQ(ierr); 593*70034214SMatthew G. Knepley ierr = PetscSFGetGraph(pointSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr); 594*70034214SMatthew G. Knepley ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr); 595*70034214SMatthew G. Knepley for (p=0; p<numRoots; p++) { 596*70034214SMatthew G. Knepley rowners[p].rank = -1; 597*70034214SMatthew G. Knepley rowners[p].index = -1; 598*70034214SMatthew G. Knepley } 599*70034214SMatthew G. Knepley if (origCellPart) { 600*70034214SMatthew G. Knepley /* Make sure points in the original partition are not assigned to other procs */ 601*70034214SMatthew G. Knepley const PetscInt *origPoints; 602*70034214SMatthew G. Knepley 603*70034214SMatthew G. Knepley ierr = DMPlexCreatePartitionClosure(dm, origCellPartSection, origCellPart, &origPartSection, &origPart);CHKERRQ(ierr); 604*70034214SMatthew G. Knepley ierr = ISGetIndices(origPart, &origPoints);CHKERRQ(ierr); 605*70034214SMatthew G. Knepley for (p = 0; p < numProcs; ++p) { 606*70034214SMatthew G. Knepley PetscInt dof, off, d; 607*70034214SMatthew G. Knepley 608*70034214SMatthew G. Knepley ierr = PetscSectionGetDof(origPartSection, p, &dof);CHKERRQ(ierr); 609*70034214SMatthew G. Knepley ierr = PetscSectionGetOffset(origPartSection, p, &off);CHKERRQ(ierr); 610*70034214SMatthew G. Knepley for (d = off; d < off+dof; ++d) { 611*70034214SMatthew G. Knepley rowners[origPoints[d]].rank = p; 612*70034214SMatthew G. Knepley } 613*70034214SMatthew G. Knepley } 614*70034214SMatthew G. Knepley ierr = ISRestoreIndices(origPart, &origPoints);CHKERRQ(ierr); 615*70034214SMatthew G. Knepley ierr = ISDestroy(&origPart);CHKERRQ(ierr); 616*70034214SMatthew G. Knepley ierr = PetscSectionDestroy(&origPartSection);CHKERRQ(ierr); 617*70034214SMatthew G. Knepley } 618*70034214SMatthew G. Knepley ierr = ISDestroy(&origCellPart);CHKERRQ(ierr); 619*70034214SMatthew G. Knepley ierr = PetscSectionDestroy(&origCellPartSection);CHKERRQ(ierr); 620*70034214SMatthew G. Knepley 621*70034214SMatthew G. Knepley ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 622*70034214SMatthew G. Knepley ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 623*70034214SMatthew G. Knepley for (p = 0; p < numLeaves; ++p) { 624*70034214SMatthew G. Knepley if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */ 625*70034214SMatthew G. Knepley lowners[p].rank = rank; 626*70034214SMatthew G. Knepley lowners[p].index = leaves ? leaves[p] : p; 627*70034214SMatthew G. Knepley } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */ 628*70034214SMatthew G. Knepley lowners[p].rank = -2; 629*70034214SMatthew G. Knepley lowners[p].index = -2; 630*70034214SMatthew G. Knepley } 631*70034214SMatthew G. Knepley } 632*70034214SMatthew G. Knepley for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */ 633*70034214SMatthew G. Knepley rowners[p].rank = -3; 634*70034214SMatthew G. Knepley rowners[p].index = -3; 635*70034214SMatthew G. Knepley } 636*70034214SMatthew G. Knepley ierr = PetscSFReduceBegin(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr); 637*70034214SMatthew G. Knepley ierr = PetscSFReduceEnd(pointSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr); 638*70034214SMatthew G. Knepley ierr = PetscSFBcastBegin(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 639*70034214SMatthew G. Knepley ierr = PetscSFBcastEnd(pointSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 640*70034214SMatthew G. Knepley for (p = 0; p < numLeaves; ++p) { 641*70034214SMatthew G. Knepley if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed"); 642*70034214SMatthew G. Knepley if (lowners[p].rank != rank) ++numGhostPoints; 643*70034214SMatthew G. Knepley } 644*70034214SMatthew G. Knepley ierr = PetscMalloc1(numGhostPoints, &ghostPoints);CHKERRQ(ierr); 645*70034214SMatthew G. Knepley ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr); 646*70034214SMatthew G. Knepley for (p = 0, gp = 0; p < numLeaves; ++p) { 647*70034214SMatthew G. Knepley if (lowners[p].rank != rank) { 648*70034214SMatthew G. Knepley ghostPoints[gp] = leaves ? leaves[p] : p; 649*70034214SMatthew G. Knepley remotePoints[gp].rank = lowners[p].rank; 650*70034214SMatthew G. Knepley remotePoints[gp].index = lowners[p].index; 651*70034214SMatthew G. Knepley ++gp; 652*70034214SMatthew G. Knepley } 653*70034214SMatthew G. Knepley } 654*70034214SMatthew G. Knepley ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr); 655*70034214SMatthew G. Knepley ierr = PetscSFSetGraph((*dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 656*70034214SMatthew G. Knepley ierr = PetscSFSetFromOptions((*dmParallel)->sf);CHKERRQ(ierr); 657*70034214SMatthew G. Knepley } 658*70034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr); 659*70034214SMatthew G. Knepley /* Cleanup */ 660*70034214SMatthew G. Knepley if (sf) {*sf = pointSF;} 661*70034214SMatthew G. Knepley else {ierr = PetscSFDestroy(&pointSF);CHKERRQ(ierr);} 662*70034214SMatthew G. Knepley ierr = DMSetFromOptions(*dmParallel);CHKERRQ(ierr); 663*70034214SMatthew G. Knepley ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr); 664*70034214SMatthew G. Knepley PetscFunctionReturn(0); 665*70034214SMatthew G. Knepley } 666