1 #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2 3 #undef __FUNCT__ 4 #define __FUNCT__ "DMPlexSetAdjacencyUseCone" 5 /*@ 6 DMPlexSetAdjacencyUseCone - Define adjacency in the mesh using either the cone or the support first 7 8 Input Parameters: 9 + dm - The DM object 10 - useCone - Flag to use the cone first 11 12 Level: intermediate 13 14 Notes: 15 $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 16 $ FVM: Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 17 $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 18 19 .seealso: DMPlexGetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator() 20 @*/ 21 PetscErrorCode DMPlexSetAdjacencyUseCone(DM dm, PetscBool useCone) 22 { 23 DM_Plex *mesh = (DM_Plex *) dm->data; 24 25 PetscFunctionBegin; 26 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 27 mesh->useCone = useCone; 28 PetscFunctionReturn(0); 29 } 30 31 #undef __FUNCT__ 32 #define __FUNCT__ "DMPlexGetAdjacencyUseCone" 33 /*@ 34 DMPlexGetAdjacencyUseCone - Query whether adjacency in the mesh uses the cone or the support first 35 36 Input Parameter: 37 . dm - The DM object 38 39 Output Parameter: 40 . useCone - Flag to use the cone first 41 42 Level: intermediate 43 44 Notes: 45 $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 46 $ FVM: Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 47 $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 48 49 .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexGetAdjacencyUseClosure(), DMPlexDistribute(), DMPlexPreallocateOperator() 50 @*/ 51 PetscErrorCode DMPlexGetAdjacencyUseCone(DM dm, PetscBool *useCone) 52 { 53 DM_Plex *mesh = (DM_Plex *) dm->data; 54 55 PetscFunctionBegin; 56 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 57 PetscValidIntPointer(useCone, 2); 58 *useCone = mesh->useCone; 59 PetscFunctionReturn(0); 60 } 61 62 #undef __FUNCT__ 63 #define __FUNCT__ "DMPlexSetAdjacencyUseClosure" 64 /*@ 65 DMPlexSetAdjacencyUseClosure - Define adjacency in the mesh using the transitive closure 66 67 Input Parameters: 68 + dm - The DM object 69 - useClosure - Flag to use the closure 70 71 Level: intermediate 72 73 Notes: 74 $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 75 $ FVM: Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 76 $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 77 78 .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator() 79 @*/ 80 PetscErrorCode DMPlexSetAdjacencyUseClosure(DM dm, PetscBool useClosure) 81 { 82 DM_Plex *mesh = (DM_Plex *) dm->data; 83 84 PetscFunctionBegin; 85 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 86 mesh->useClosure = useClosure; 87 PetscFunctionReturn(0); 88 } 89 90 #undef __FUNCT__ 91 #define __FUNCT__ "DMPlexGetAdjacencyUseClosure" 92 /*@ 93 DMPlexGetAdjacencyUseClosure - Query whether adjacency in the mesh uses the transitive closure 94 95 Input Parameter: 96 . dm - The DM object 97 98 Output Parameter: 99 . useClosure - Flag to use the closure 100 101 Level: intermediate 102 103 Notes: 104 $ FEM: Two points p and q are adjacent if q \in closure(star(p)), useCone = PETSC_FALSE, useClosure = PETSC_TRUE 105 $ FVM: Two points p and q are adjacent if q \in support(p+cone(p)), useCone = PETSC_TRUE, useClosure = PETSC_FALSE 106 $ FVM++: Two points p and q are adjacent if q \in star(closure(p)), useCone = PETSC_TRUE, useClosure = PETSC_TRUE 107 108 .seealso: DMPlexSetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator() 109 @*/ 110 PetscErrorCode DMPlexGetAdjacencyUseClosure(DM dm, PetscBool *useClosure) 111 { 112 DM_Plex *mesh = (DM_Plex *) dm->data; 113 114 PetscFunctionBegin; 115 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 116 PetscValidIntPointer(useClosure, 2); 117 *useClosure = mesh->useClosure; 118 PetscFunctionReturn(0); 119 } 120 121 #undef __FUNCT__ 122 #define __FUNCT__ "DMPlexSetAdjacencyUseAnchors" 123 /*@ 124 DMPlexSetAdjacencyUseAnchors - Define adjacency in the mesh using the point-to-point constraints. 125 126 Input Parameters: 127 + dm - The DM object 128 - useAnchors - Flag to use the constraints. If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place. 129 130 Level: intermediate 131 132 .seealso: DMPlexGetAdjacencyUseClosure(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors() 133 @*/ 134 PetscErrorCode DMPlexSetAdjacencyUseAnchors(DM dm, PetscBool useAnchors) 135 { 136 DM_Plex *mesh = (DM_Plex *) dm->data; 137 138 PetscFunctionBegin; 139 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 140 mesh->useAnchors = useAnchors; 141 PetscFunctionReturn(0); 142 } 143 144 #undef __FUNCT__ 145 #define __FUNCT__ "DMPlexGetAdjacencyUseAnchors" 146 /*@ 147 DMPlexGetAdjacencyUseAnchors - Query whether adjacency in the mesh uses the point-to-point constraints. 148 149 Input Parameter: 150 . dm - The DM object 151 152 Output Parameter: 153 . useAnchors - Flag to use the closure. If PETSC_TRUE, then constrained points are omitted from DMPlexGetAdjacency(), and their anchor points appear in their place. 154 155 Level: intermediate 156 157 .seealso: DMPlexSetAdjacencyUseAnchors(), DMPlexSetAdjacencyUseCone(), DMPlexGetAdjacencyUseCone(), DMPlexDistribute(), DMPlexPreallocateOperator(), DMPlexSetAnchors() 158 @*/ 159 PetscErrorCode DMPlexGetAdjacencyUseAnchors(DM dm, PetscBool *useAnchors) 160 { 161 DM_Plex *mesh = (DM_Plex *) dm->data; 162 163 PetscFunctionBegin; 164 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 165 PetscValidIntPointer(useAnchors, 2); 166 *useAnchors = mesh->useAnchors; 167 PetscFunctionReturn(0); 168 } 169 170 #undef __FUNCT__ 171 #define __FUNCT__ "DMPlexGetAdjacency_Cone_Internal" 172 static PetscErrorCode DMPlexGetAdjacency_Cone_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[]) 173 { 174 const PetscInt *cone = NULL; 175 PetscInt numAdj = 0, maxAdjSize = *adjSize, coneSize, c; 176 PetscErrorCode ierr; 177 178 PetscFunctionBeginHot; 179 ierr = DMPlexGetConeSize(dm, p, &coneSize);CHKERRQ(ierr); 180 ierr = DMPlexGetCone(dm, p, &cone);CHKERRQ(ierr); 181 for (c = 0; c <= coneSize; ++c) { 182 const PetscInt point = !c ? p : cone[c-1]; 183 const PetscInt *support = NULL; 184 PetscInt supportSize, s, q; 185 186 ierr = DMPlexGetSupportSize(dm, point, &supportSize);CHKERRQ(ierr); 187 ierr = DMPlexGetSupport(dm, point, &support);CHKERRQ(ierr); 188 for (s = 0; s < supportSize; ++s) { 189 for (q = 0; q < numAdj || (adj[numAdj++] = support[s],0); ++q) { 190 if (support[s] == adj[q]) break; 191 } 192 if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize); 193 } 194 } 195 *adjSize = numAdj; 196 PetscFunctionReturn(0); 197 } 198 199 #undef __FUNCT__ 200 #define __FUNCT__ "DMPlexGetAdjacency_Support_Internal" 201 static PetscErrorCode DMPlexGetAdjacency_Support_Internal(DM dm, PetscInt p, PetscInt *adjSize, PetscInt adj[]) 202 { 203 const PetscInt *support = NULL; 204 PetscInt numAdj = 0, maxAdjSize = *adjSize, supportSize, s; 205 PetscErrorCode ierr; 206 207 PetscFunctionBeginHot; 208 ierr = DMPlexGetSupportSize(dm, p, &supportSize);CHKERRQ(ierr); 209 ierr = DMPlexGetSupport(dm, p, &support);CHKERRQ(ierr); 210 for (s = 0; s <= supportSize; ++s) { 211 const PetscInt point = !s ? p : support[s-1]; 212 const PetscInt *cone = NULL; 213 PetscInt coneSize, c, q; 214 215 ierr = DMPlexGetConeSize(dm, point, &coneSize);CHKERRQ(ierr); 216 ierr = DMPlexGetCone(dm, point, &cone);CHKERRQ(ierr); 217 for (c = 0; c < coneSize; ++c) { 218 for (q = 0; q < numAdj || (adj[numAdj++] = cone[c],0); ++q) { 219 if (cone[c] == adj[q]) break; 220 } 221 if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize); 222 } 223 } 224 *adjSize = numAdj; 225 PetscFunctionReturn(0); 226 } 227 228 #undef __FUNCT__ 229 #define __FUNCT__ "DMPlexGetAdjacency_Transitive_Internal" 230 static PetscErrorCode DMPlexGetAdjacency_Transitive_Internal(DM dm, PetscInt p, PetscBool useClosure, PetscInt *adjSize, PetscInt adj[]) 231 { 232 PetscInt *star = NULL; 233 PetscInt numAdj = 0, maxAdjSize = *adjSize, starSize, s; 234 PetscErrorCode ierr; 235 236 PetscFunctionBeginHot; 237 ierr = DMPlexGetTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr); 238 for (s = 0; s < starSize*2; s += 2) { 239 const PetscInt *closure = NULL; 240 PetscInt closureSize, c, q; 241 242 ierr = DMPlexGetTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr); 243 for (c = 0; c < closureSize*2; c += 2) { 244 for (q = 0; q < numAdj || (adj[numAdj++] = closure[c],0); ++q) { 245 if (closure[c] == adj[q]) break; 246 } 247 if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize); 248 } 249 ierr = DMPlexRestoreTransitiveClosure(dm, star[s], (PetscBool)!useClosure, &closureSize, (PetscInt**) &closure);CHKERRQ(ierr); 250 } 251 ierr = DMPlexRestoreTransitiveClosure(dm, p, useClosure, &starSize, &star);CHKERRQ(ierr); 252 *adjSize = numAdj; 253 PetscFunctionReturn(0); 254 } 255 256 #undef __FUNCT__ 257 #define __FUNCT__ "DMPlexGetAdjacency_Internal" 258 PetscErrorCode DMPlexGetAdjacency_Internal(DM dm, PetscInt p, PetscBool useCone, PetscBool useTransitiveClosure, PetscBool useAnchors, PetscInt *adjSize, PetscInt *adj[]) 259 { 260 static PetscInt asiz = 0; 261 PetscInt maxAnchors = 1; 262 PetscInt aStart = -1, aEnd = -1; 263 PetscInt maxAdjSize; 264 PetscSection aSec = NULL; 265 IS aIS = NULL; 266 const PetscInt *anchors; 267 PetscErrorCode ierr; 268 269 PetscFunctionBeginHot; 270 if (useAnchors) { 271 ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr); 272 if (aSec) { 273 ierr = PetscSectionGetMaxDof(aSec,&maxAnchors);CHKERRQ(ierr); 274 maxAnchors = PetscMax(1,maxAnchors); 275 ierr = PetscSectionGetChart(aSec,&aStart,&aEnd);CHKERRQ(ierr); 276 ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr); 277 } 278 } 279 if (!*adj) { 280 PetscInt depth, coneSeries, supportSeries, maxC, maxS, pStart, pEnd; 281 282 ierr = DMPlexGetChart(dm, &pStart,&pEnd);CHKERRQ(ierr); 283 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 284 ierr = DMPlexGetMaxSizes(dm, &maxC, &maxS);CHKERRQ(ierr); 285 coneSeries = (maxC > 1) ? ((PetscPowInt(maxC,depth+1)-1)/(maxC-1)) : depth+1; 286 supportSeries = (maxS > 1) ? ((PetscPowInt(maxS,depth+1)-1)/(maxS-1)) : depth+1; 287 asiz = PetscMax(PetscPowInt(maxS,depth)*coneSeries,PetscPowInt(maxC,depth)*supportSeries); 288 asiz *= maxAnchors; 289 asiz = PetscMin(asiz,pEnd-pStart); 290 ierr = PetscMalloc1(asiz,adj);CHKERRQ(ierr); 291 } 292 if (*adjSize < 0) *adjSize = asiz; 293 maxAdjSize = *adjSize; 294 if (useTransitiveClosure) { 295 ierr = DMPlexGetAdjacency_Transitive_Internal(dm, p, useCone, adjSize, *adj);CHKERRQ(ierr); 296 } else if (useCone) { 297 ierr = DMPlexGetAdjacency_Cone_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr); 298 } else { 299 ierr = DMPlexGetAdjacency_Support_Internal(dm, p, adjSize, *adj);CHKERRQ(ierr); 300 } 301 if (useAnchors && aSec) { 302 PetscInt origSize = *adjSize; 303 PetscInt numAdj = origSize; 304 PetscInt i = 0, j; 305 PetscInt *orig = *adj; 306 307 while (i < origSize) { 308 PetscInt p = orig[i]; 309 PetscInt aDof = 0; 310 311 if (p >= aStart && p < aEnd) { 312 ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr); 313 } 314 if (aDof) { 315 PetscInt aOff; 316 PetscInt s, q; 317 318 for (j = i + 1; j < numAdj; j++) { 319 orig[j - 1] = orig[j]; 320 } 321 origSize--; 322 numAdj--; 323 ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr); 324 for (s = 0; s < aDof; ++s) { 325 for (q = 0; q < numAdj || (orig[numAdj++] = anchors[aOff+s],0); ++q) { 326 if (anchors[aOff+s] == orig[q]) break; 327 } 328 if (numAdj > maxAdjSize) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Invalid mesh exceeded adjacency allocation (%D)", maxAdjSize); 329 } 330 } 331 else { 332 i++; 333 } 334 } 335 *adjSize = numAdj; 336 ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr); 337 } 338 PetscFunctionReturn(0); 339 } 340 341 #undef __FUNCT__ 342 #define __FUNCT__ "DMPlexGetAdjacency" 343 /*@ 344 DMPlexGetAdjacency - Return all points adjacent to the given point 345 346 Input Parameters: 347 + dm - The DM object 348 . p - The point 349 . adjSize - The maximum size of adj if it is non-NULL, or PETSC_DETERMINE 350 - adj - Either NULL so that the array is allocated, or an existing array with size adjSize 351 352 Output Parameters: 353 + adjSize - The number of adjacent points 354 - adj - The adjacent points 355 356 Level: advanced 357 358 Notes: The user must PetscFree the adj array if it was not passed in. 359 360 .seealso: DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure(), DMPlexDistribute(), DMCreateMatrix(), DMPlexPreallocateOperator() 361 @*/ 362 PetscErrorCode DMPlexGetAdjacency(DM dm, PetscInt p, PetscInt *adjSize, PetscInt *adj[]) 363 { 364 DM_Plex *mesh = (DM_Plex *) dm->data; 365 PetscErrorCode ierr; 366 367 PetscFunctionBeginHot; 368 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 369 PetscValidPointer(adjSize,3); 370 PetscValidPointer(adj,4); 371 ierr = DMPlexGetAdjacency_Internal(dm, p, mesh->useCone, mesh->useClosure, mesh->useAnchors, adjSize, adj);CHKERRQ(ierr); 372 PetscFunctionReturn(0); 373 } 374 #undef __FUNCT__ 375 #define __FUNCT__ "DMPlexCreateTwoSidedProcessSF" 376 /*@ 377 DMPlexCreateTwoSidedProcessSF - Create an SF which just has process connectivity 378 379 Collective on DM 380 381 Input Parameters: 382 + dm - The DM 383 - sfPoint - The PetscSF which encodes point connectivity 384 385 Output Parameters: 386 + processRanks - A list of process neighbors, or NULL 387 - sfProcess - An SF encoding the two-sided process connectivity, or NULL 388 389 Level: developer 390 391 .seealso: PetscSFCreate(), DMPlexCreateProcessSF() 392 @*/ 393 PetscErrorCode DMPlexCreateTwoSidedProcessSF(DM dm, PetscSF sfPoint, PetscSection rootRankSection, IS rootRanks, PetscSection leafRankSection, IS leafRanks, IS *processRanks, PetscSF *sfProcess) 394 { 395 const PetscSFNode *remotePoints; 396 PetscInt *localPointsNew; 397 PetscSFNode *remotePointsNew; 398 const PetscInt *nranks; 399 PetscInt *ranksNew; 400 PetscBT neighbors; 401 PetscInt pStart, pEnd, p, numLeaves, l, numNeighbors, n; 402 PetscMPIInt numProcs, proc, rank; 403 PetscErrorCode ierr; 404 405 PetscFunctionBegin; 406 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 407 PetscValidHeaderSpecific(sfPoint, PETSCSF_CLASSID, 2); 408 if (processRanks) {PetscValidPointer(processRanks, 3);} 409 if (sfProcess) {PetscValidPointer(sfProcess, 4);} 410 ierr = MPI_Comm_size(PetscObjectComm((PetscObject) dm), &numProcs);CHKERRQ(ierr); 411 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 412 ierr = PetscSFGetGraph(sfPoint, NULL, &numLeaves, NULL, &remotePoints);CHKERRQ(ierr); 413 ierr = PetscBTCreate(numProcs, &neighbors);CHKERRQ(ierr); 414 ierr = PetscBTMemzero(numProcs, neighbors);CHKERRQ(ierr); 415 /* Compute root-to-leaf process connectivity */ 416 ierr = PetscSectionGetChart(rootRankSection, &pStart, &pEnd);CHKERRQ(ierr); 417 ierr = ISGetIndices(rootRanks, &nranks);CHKERRQ(ierr); 418 for (p = pStart; p < pEnd; ++p) { 419 PetscInt ndof, noff, n; 420 421 ierr = PetscSectionGetDof(rootRankSection, p, &ndof);CHKERRQ(ierr); 422 ierr = PetscSectionGetOffset(rootRankSection, p, &noff);CHKERRQ(ierr); 423 for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);} 424 } 425 ierr = ISRestoreIndices(rootRanks, &nranks);CHKERRQ(ierr); 426 /* Compute leaf-to-neighbor process connectivity */ 427 ierr = PetscSectionGetChart(leafRankSection, &pStart, &pEnd);CHKERRQ(ierr); 428 ierr = ISGetIndices(leafRanks, &nranks);CHKERRQ(ierr); 429 for (p = pStart; p < pEnd; ++p) { 430 PetscInt ndof, noff, n; 431 432 ierr = PetscSectionGetDof(leafRankSection, p, &ndof);CHKERRQ(ierr); 433 ierr = PetscSectionGetOffset(leafRankSection, p, &noff);CHKERRQ(ierr); 434 for (n = 0; n < ndof; ++n) {ierr = PetscBTSet(neighbors, nranks[noff+n]);} 435 } 436 ierr = ISRestoreIndices(leafRanks, &nranks);CHKERRQ(ierr); 437 /* Compute leaf-to-root process connectivity */ 438 for (l = 0; l < numLeaves; ++l) {PetscBTSet(neighbors, remotePoints[l].rank);} 439 /* Calculate edges */ 440 PetscBTClear(neighbors, rank); 441 for(proc = 0, numNeighbors = 0; proc < numProcs; ++proc) {if (PetscBTLookup(neighbors, proc)) ++numNeighbors;} 442 ierr = PetscMalloc1(numNeighbors, &ranksNew);CHKERRQ(ierr); 443 ierr = PetscMalloc1(numNeighbors, &localPointsNew);CHKERRQ(ierr); 444 ierr = PetscMalloc1(numNeighbors, &remotePointsNew);CHKERRQ(ierr); 445 for(proc = 0, n = 0; proc < numProcs; ++proc) { 446 if (PetscBTLookup(neighbors, proc)) { 447 ranksNew[n] = proc; 448 localPointsNew[n] = proc; 449 remotePointsNew[n].index = rank; 450 remotePointsNew[n].rank = proc; 451 ++n; 452 } 453 } 454 ierr = PetscBTDestroy(&neighbors);CHKERRQ(ierr); 455 if (processRanks) {ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm), numNeighbors, ranksNew, PETSC_OWN_POINTER, processRanks);CHKERRQ(ierr);} 456 else {ierr = PetscFree(ranksNew);CHKERRQ(ierr);} 457 if (sfProcess) { 458 ierr = PetscSFCreate(PetscObjectComm((PetscObject)dm), sfProcess);CHKERRQ(ierr); 459 ierr = PetscObjectSetName((PetscObject) *sfProcess, "Two-Sided Process SF");CHKERRQ(ierr); 460 ierr = PetscSFSetFromOptions(*sfProcess);CHKERRQ(ierr); 461 ierr = PetscSFSetGraph(*sfProcess, numProcs, numNeighbors, localPointsNew, PETSC_OWN_POINTER, remotePointsNew, PETSC_OWN_POINTER);CHKERRQ(ierr); 462 } 463 PetscFunctionReturn(0); 464 } 465 466 #undef __FUNCT__ 467 #define __FUNCT__ "DMPlexDistributeOwnership" 468 /*@ 469 DMPlexDistributeOwnership - Compute owner information for shared points. This basically gets two-sided for an SF. 470 471 Collective on DM 472 473 Input Parameter: 474 . dm - The DM 475 476 Output Parameters: 477 + rootSection - The number of leaves for a given root point 478 . rootrank - The rank of each edge into the root point 479 . leafSection - The number of processes sharing a given leaf point 480 - leafrank - The rank of each process sharing a leaf point 481 482 Level: developer 483 484 .seealso: DMPlexCreateOverlap() 485 @*/ 486 PetscErrorCode DMPlexDistributeOwnership(DM dm, PetscSection rootSection, IS *rootrank, PetscSection leafSection, IS *leafrank) 487 { 488 MPI_Comm comm; 489 PetscSF sfPoint; 490 const PetscInt *rootdegree; 491 PetscInt *myrank, *remoterank; 492 PetscInt pStart, pEnd, p, nedges; 493 PetscMPIInt rank; 494 PetscErrorCode ierr; 495 496 PetscFunctionBegin; 497 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 498 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 499 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 500 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 501 /* Compute number of leaves for each root */ 502 ierr = PetscObjectSetName((PetscObject) rootSection, "Root Section");CHKERRQ(ierr); 503 ierr = PetscSectionSetChart(rootSection, pStart, pEnd);CHKERRQ(ierr); 504 ierr = PetscSFComputeDegreeBegin(sfPoint, &rootdegree);CHKERRQ(ierr); 505 ierr = PetscSFComputeDegreeEnd(sfPoint, &rootdegree);CHKERRQ(ierr); 506 for (p = pStart; p < pEnd; ++p) {ierr = PetscSectionSetDof(rootSection, p, rootdegree[p-pStart]);CHKERRQ(ierr);} 507 ierr = PetscSectionSetUp(rootSection);CHKERRQ(ierr); 508 /* Gather rank of each leaf to root */ 509 ierr = PetscSectionGetStorageSize(rootSection, &nedges);CHKERRQ(ierr); 510 ierr = PetscMalloc1(pEnd-pStart, &myrank);CHKERRQ(ierr); 511 ierr = PetscMalloc1(nedges, &remoterank);CHKERRQ(ierr); 512 for (p = 0; p < pEnd-pStart; ++p) myrank[p] = rank; 513 ierr = PetscSFGatherBegin(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr); 514 ierr = PetscSFGatherEnd(sfPoint, MPIU_INT, myrank, remoterank);CHKERRQ(ierr); 515 ierr = PetscFree(myrank);CHKERRQ(ierr); 516 ierr = ISCreateGeneral(comm, nedges, remoterank, PETSC_OWN_POINTER, rootrank);CHKERRQ(ierr); 517 /* Distribute remote ranks to leaves */ 518 ierr = PetscObjectSetName((PetscObject) leafSection, "Leaf Section");CHKERRQ(ierr); 519 ierr = DMPlexDistributeFieldIS(dm, sfPoint, rootSection, *rootrank, leafSection, leafrank);CHKERRQ(ierr); 520 PetscFunctionReturn(0); 521 } 522 523 #undef __FUNCT__ 524 #define __FUNCT__ "DMPlexCreateOverlap" 525 /*@C 526 DMPlexCreateOverlap - Compute owner information for shared points. This basically gets two-sided for an SF. 527 528 Collective on DM 529 530 Input Parameters: 531 + dm - The DM 532 . levels - Number of overlap levels 533 . rootSection - The number of leaves for a given root point 534 . rootrank - The rank of each edge into the root point 535 . leafSection - The number of processes sharing a given leaf point 536 - leafrank - The rank of each process sharing a leaf point 537 538 Output Parameters: 539 + ovLabel - DMLabel containing remote overlap contributions as point/rank pairings 540 541 Level: developer 542 543 .seealso: DMPlexDistributeOwnership(), DMPlexDistribute() 544 @*/ 545 PetscErrorCode DMPlexCreateOverlap(DM dm, PetscInt levels, PetscSection rootSection, IS rootrank, PetscSection leafSection, IS leafrank, DMLabel *ovLabel) 546 { 547 MPI_Comm comm; 548 DMLabel ovAdjByRank; /* A DMLabel containing all points adjacent to shared points, separated by rank (value in label) */ 549 PetscSF sfPoint, sfProc; 550 const PetscSFNode *remote; 551 const PetscInt *local; 552 const PetscInt *nrank, *rrank; 553 PetscInt *adj = NULL; 554 PetscInt pStart, pEnd, p, sStart, sEnd, nleaves, l; 555 PetscMPIInt rank, numProcs; 556 PetscBool useCone, useClosure, flg; 557 PetscErrorCode ierr; 558 559 PetscFunctionBegin; 560 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 561 ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 562 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 563 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 564 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 565 ierr = PetscSectionGetChart(leafSection, &sStart, &sEnd);CHKERRQ(ierr); 566 ierr = PetscSFGetGraph(sfPoint, NULL, &nleaves, &local, &remote);CHKERRQ(ierr); 567 ierr = DMLabelCreate("Overlap adjacency", &ovAdjByRank);CHKERRQ(ierr); 568 /* Handle leaves: shared with the root point */ 569 for (l = 0; l < nleaves; ++l) { 570 PetscInt adjSize = PETSC_DETERMINE, a; 571 572 ierr = DMPlexGetAdjacency(dm, local[l], &adjSize, &adj);CHKERRQ(ierr); 573 for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remote[l].rank);CHKERRQ(ierr);} 574 } 575 ierr = ISGetIndices(rootrank, &rrank);CHKERRQ(ierr); 576 ierr = ISGetIndices(leafrank, &nrank);CHKERRQ(ierr); 577 /* Handle roots */ 578 for (p = pStart; p < pEnd; ++p) { 579 PetscInt adjSize = PETSC_DETERMINE, neighbors = 0, noff, n, a; 580 581 if ((p >= sStart) && (p < sEnd)) { 582 /* Some leaves share a root with other leaves on different processes */ 583 ierr = PetscSectionGetDof(leafSection, p, &neighbors);CHKERRQ(ierr); 584 if (neighbors) { 585 ierr = PetscSectionGetOffset(leafSection, p, &noff);CHKERRQ(ierr); 586 ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr); 587 for (n = 0; n < neighbors; ++n) { 588 const PetscInt remoteRank = nrank[noff+n]; 589 590 if (remoteRank == rank) continue; 591 for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);} 592 } 593 } 594 } 595 /* Roots are shared with leaves */ 596 ierr = PetscSectionGetDof(rootSection, p, &neighbors);CHKERRQ(ierr); 597 if (!neighbors) continue; 598 ierr = PetscSectionGetOffset(rootSection, p, &noff);CHKERRQ(ierr); 599 ierr = DMPlexGetAdjacency(dm, p, &adjSize, &adj);CHKERRQ(ierr); 600 for (n = 0; n < neighbors; ++n) { 601 const PetscInt remoteRank = rrank[noff+n]; 602 603 if (remoteRank == rank) continue; 604 for (a = 0; a < adjSize; ++a) {ierr = DMLabelSetValue(ovAdjByRank, adj[a], remoteRank);CHKERRQ(ierr);} 605 } 606 } 607 ierr = PetscFree(adj);CHKERRQ(ierr); 608 ierr = ISRestoreIndices(rootrank, &rrank);CHKERRQ(ierr); 609 ierr = ISRestoreIndices(leafrank, &nrank);CHKERRQ(ierr); 610 /* Add additional overlap levels */ 611 for (l = 1; l < levels; l++) {ierr = DMPlexPartitionLabelAdjacency(dm, ovAdjByRank);CHKERRQ(ierr);} 612 /* We require the closure in the overlap */ 613 ierr = DMPlexGetAdjacencyUseCone(dm, &useCone);CHKERRQ(ierr); 614 ierr = DMPlexGetAdjacencyUseClosure(dm, &useClosure);CHKERRQ(ierr); 615 if (useCone || !useClosure) { 616 ierr = DMPlexPartitionLabelClosure(dm, ovAdjByRank);CHKERRQ(ierr); 617 } 618 ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-overlap_view", &flg);CHKERRQ(ierr); 619 if (flg) { 620 ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr); 621 ierr = DMLabelView(ovAdjByRank, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 622 } 623 /* Make process SF and invert sender to receiver label */ 624 ierr = DMPlexCreateTwoSidedProcessSF(dm, sfPoint, rootSection, rootrank, leafSection, leafrank, NULL, &sfProc);CHKERRQ(ierr); 625 ierr = DMLabelCreate("Overlap label", ovLabel);CHKERRQ(ierr); 626 ierr = DMPlexPartitionLabelInvert(dm, ovAdjByRank, sfProc, *ovLabel);CHKERRQ(ierr); 627 /* Add owned points, except for shared local points */ 628 for (p = pStart; p < pEnd; ++p) {ierr = DMLabelSetValue(*ovLabel, p, rank);CHKERRQ(ierr);} 629 for (l = 0; l < nleaves; ++l) { 630 ierr = DMLabelClearValue(*ovLabel, local[l], rank);CHKERRQ(ierr); 631 ierr = DMLabelSetValue(*ovLabel, remote[l].index, remote[l].rank);CHKERRQ(ierr); 632 } 633 /* Clean up */ 634 ierr = DMLabelDestroy(&ovAdjByRank);CHKERRQ(ierr); 635 ierr = PetscSFDestroy(&sfProc);CHKERRQ(ierr); 636 PetscFunctionReturn(0); 637 } 638 639 #undef __FUNCT__ 640 #define __FUNCT__ "DMPlexCreateOverlapMigrationSF" 641 PetscErrorCode DMPlexCreateOverlapMigrationSF(DM dm, PetscSF overlapSF, PetscSF *migrationSF) 642 { 643 MPI_Comm comm; 644 PetscMPIInt rank, numProcs; 645 PetscInt d, dim, p, pStart, pEnd, nroots, nleaves, newLeaves, point, numSharedPoints; 646 PetscInt *pointDepths, *remoteDepths, *ilocal; 647 PetscInt *depthRecv, *depthShift, *depthIdx; 648 PetscSFNode *iremote; 649 PetscSF pointSF; 650 const PetscInt *sharedLocal; 651 const PetscSFNode *overlapRemote, *sharedRemote; 652 PetscErrorCode ierr; 653 654 PetscFunctionBegin; 655 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 656 ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 657 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 658 ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 659 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 660 661 /* Before building the migration SF we need to know the new stratum offsets */ 662 ierr = PetscSFGetGraph(overlapSF, &nroots, &nleaves, NULL, &overlapRemote);CHKERRQ(ierr); 663 ierr = PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);CHKERRQ(ierr); 664 for (d=0; d<dim+1; d++) { 665 ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 666 for (p=pStart; p<pEnd; p++) pointDepths[p] = d; 667 } 668 for (p=0; p<nleaves; p++) remoteDepths[p] = -1; 669 ierr = PetscSFBcastBegin(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr); 670 ierr = PetscSFBcastEnd(overlapSF, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr); 671 672 /* Count recevied points in each stratum and compute the internal strata shift */ 673 ierr = PetscMalloc3(dim+1, &depthRecv, dim+1, &depthShift, dim+1, &depthIdx);CHKERRQ(ierr); 674 for (d=0; d<dim+1; d++) depthRecv[d]=0; 675 for (p=0; p<nleaves; p++) depthRecv[remoteDepths[p]]++; 676 depthShift[dim] = 0; 677 for (d=0; d<dim; d++) depthShift[d] = depthRecv[dim]; 678 for (d=1; d<dim; d++) depthShift[d] += depthRecv[0]; 679 for (d=dim-2; d>0; d--) depthShift[d] += depthRecv[d+1]; 680 for (d=0; d<dim+1; d++) { 681 ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 682 depthIdx[d] = pStart + depthShift[d]; 683 } 684 685 /* Form the overlap SF build an SF that describes the full overlap migration SF */ 686 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 687 newLeaves = pEnd - pStart + nleaves; 688 ierr = PetscMalloc1(newLeaves, &ilocal);CHKERRQ(ierr); 689 ierr = PetscMalloc1(newLeaves, &iremote);CHKERRQ(ierr); 690 /* First map local points to themselves */ 691 for (d=0; d<dim+1; d++) { 692 ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 693 for (p=pStart; p<pEnd; p++) { 694 point = p + depthShift[d]; 695 ilocal[point] = point; 696 iremote[point].index = p; 697 iremote[point].rank = rank; 698 depthIdx[d]++; 699 } 700 } 701 702 /* Add in the remote roots for currently shared points */ 703 ierr = DMGetPointSF(dm, &pointSF);CHKERRQ(ierr); 704 ierr = PetscSFGetGraph(pointSF, NULL, &numSharedPoints, &sharedLocal, &sharedRemote);CHKERRQ(ierr); 705 for (d=0; d<dim+1; d++) { 706 ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 707 for (p=0; p<numSharedPoints; p++) { 708 if (pStart <= sharedLocal[p] && sharedLocal[p] < pEnd) { 709 point = sharedLocal[p] + depthShift[d]; 710 iremote[point].index = sharedRemote[p].index; 711 iremote[point].rank = sharedRemote[p].rank; 712 } 713 } 714 } 715 716 /* Now add the incoming overlap points */ 717 for (p=0; p<nleaves; p++) { 718 point = depthIdx[remoteDepths[p]]; 719 ilocal[point] = point; 720 iremote[point].index = overlapRemote[p].index; 721 iremote[point].rank = overlapRemote[p].rank; 722 depthIdx[remoteDepths[p]]++; 723 } 724 ierr = PetscFree2(pointDepths,remoteDepths);CHKERRQ(ierr); 725 726 ierr = PetscSFCreate(comm, migrationSF);CHKERRQ(ierr); 727 ierr = PetscObjectSetName((PetscObject) *migrationSF, "Overlap Migration SF");CHKERRQ(ierr); 728 ierr = PetscSFSetFromOptions(*migrationSF);CHKERRQ(ierr); 729 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 730 ierr = PetscSFSetGraph(*migrationSF, pEnd-pStart, newLeaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_OWN_POINTER);CHKERRQ(ierr); 731 732 ierr = PetscFree3(depthRecv, depthShift, depthIdx);CHKERRQ(ierr); 733 PetscFunctionReturn(0); 734 } 735 736 #undef __FUNCT__ 737 #define __FUNCT__ "DMPlexStratifyMigrationSF" 738 /*@ 739 DMPlexStratifyMigrationSF - Add partition overlap to a distributed non-overlapping DM. 740 741 Input Parameter: 742 + dm - The DM 743 - sf - A star forest with non-ordered leaves, usually defining a DM point migration 744 745 Output Parameter: 746 . migrationSF - A star forest with added leaf indirection that ensures the resulting DM is stratified 747 748 Level: developer 749 750 .seealso: DMPlexPartitionLabelCreateSF(), DMPlexDistribute(), DMPlexDistributeOverlap() 751 @*/ 752 PetscErrorCode DMPlexStratifyMigrationSF(DM dm, PetscSF sf, PetscSF *migrationSF) 753 { 754 MPI_Comm comm; 755 PetscMPIInt rank, numProcs; 756 PetscInt d, dim, depth, p, pStart, pEnd, nroots, nleaves; 757 PetscInt *pointDepths, *remoteDepths, *ilocal; 758 PetscInt *depthRecv, *depthShift, *depthIdx; 759 const PetscSFNode *iremote; 760 PetscErrorCode ierr; 761 762 PetscFunctionBegin; 763 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 764 ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 765 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 766 ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 767 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 768 769 /* Before building the migration SF we need to know the new stratum offsets */ 770 ierr = PetscSFGetGraph(sf, &nroots, &nleaves, NULL, &iremote);CHKERRQ(ierr); 771 ierr = PetscMalloc2(nroots, &pointDepths, nleaves, &remoteDepths);CHKERRQ(ierr); 772 for (d=0; d<dim+1; d++) { 773 ierr = DMPlexGetDepthStratum(dm, d, &pStart, &pEnd);CHKERRQ(ierr); 774 for (p=pStart; p<pEnd; p++) pointDepths[p] = d; 775 } 776 for (p=0; p<nleaves; p++) remoteDepths[p] = -1; 777 ierr = PetscSFBcastBegin(sf, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr); 778 ierr = PetscSFBcastEnd(sf, MPIU_INT, pointDepths, remoteDepths);CHKERRQ(ierr); 779 /* Count recevied points in each stratum and compute the internal strata shift */ 780 ierr = PetscMalloc3(dim+1, &depthRecv, dim+1, &depthShift, dim+1, &depthIdx);CHKERRQ(ierr); 781 for (d=0; d<dim+1; d++) depthRecv[d]=0; 782 for (p=0; p<nleaves; p++) depthRecv[remoteDepths[p]]++; 783 depthShift[dim] = 0; 784 for (d=0; d<dim; d++) depthShift[d] = depthRecv[dim]; 785 for (d=1; d<dim; d++) depthShift[d] += depthRecv[0]; 786 for (d=dim-2; d>0; d--) depthShift[d] += depthRecv[d+1]; 787 for (d=0; d<dim+1; d++) {depthIdx[d] = 0;} 788 /* Derive a new local permutation based on stratified indices */ 789 ierr = PetscMalloc1(nleaves, &ilocal);CHKERRQ(ierr); 790 for (p=0; p<nleaves; p++) { 791 depth = remoteDepths[p]; 792 ilocal[p] = depthShift[depth] + depthIdx[depth]; 793 depthIdx[depth]++; 794 } 795 ierr = PetscSFCreate(comm, migrationSF);CHKERRQ(ierr); 796 ierr = PetscObjectSetName((PetscObject) *migrationSF, "Migration SF");CHKERRQ(ierr); 797 ierr = PetscSFSetGraph(*migrationSF, nroots, nleaves, ilocal, PETSC_OWN_POINTER, iremote, PETSC_COPY_VALUES);CHKERRQ(ierr); 798 ierr = PetscFree2(pointDepths,remoteDepths);CHKERRQ(ierr); 799 ierr = PetscFree3(depthRecv, depthShift, depthIdx);CHKERRQ(ierr); 800 PetscFunctionReturn(0); 801 } 802 803 #undef __FUNCT__ 804 #define __FUNCT__ "DMPlexDistributeField" 805 /*@ 806 DMPlexDistributeField - Distribute field data to match a given PetscSF, usually the SF from mesh distribution 807 808 Collective on DM 809 810 Input Parameters: 811 + dm - The DMPlex object 812 . pointSF - The PetscSF describing the communication pattern 813 . originalSection - The PetscSection for existing data layout 814 - originalVec - The existing data 815 816 Output Parameters: 817 + newSection - The PetscSF describing the new data layout 818 - newVec - The new data 819 820 Level: developer 821 822 .seealso: DMPlexDistribute(), DMPlexDistributeFieldIS(), DMPlexDistributeData() 823 @*/ 824 PetscErrorCode DMPlexDistributeField(DM dm, PetscSF pointSF, PetscSection originalSection, Vec originalVec, PetscSection newSection, Vec newVec) 825 { 826 PetscSF fieldSF; 827 PetscInt *remoteOffsets, fieldSize; 828 PetscScalar *originalValues, *newValues; 829 PetscErrorCode ierr; 830 831 PetscFunctionBegin; 832 ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 833 ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); 834 835 ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); 836 ierr = VecSetSizes(newVec, fieldSize, PETSC_DETERMINE);CHKERRQ(ierr); 837 ierr = VecSetType(newVec,dm->vectype);CHKERRQ(ierr); 838 839 ierr = VecGetArray(originalVec, &originalValues);CHKERRQ(ierr); 840 ierr = VecGetArray(newVec, &newValues);CHKERRQ(ierr); 841 ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); 842 ierr = PetscSFBcastBegin(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr); 843 ierr = PetscSFBcastEnd(fieldSF, MPIU_SCALAR, originalValues, newValues);CHKERRQ(ierr); 844 ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); 845 ierr = VecRestoreArray(newVec, &newValues);CHKERRQ(ierr); 846 ierr = VecRestoreArray(originalVec, &originalValues);CHKERRQ(ierr); 847 ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 848 PetscFunctionReturn(0); 849 } 850 851 #undef __FUNCT__ 852 #define __FUNCT__ "DMPlexDistributeFieldIS" 853 /*@ 854 DMPlexDistributeFieldIS - Distribute field data to match a given PetscSF, usually the SF from mesh distribution 855 856 Collective on DM 857 858 Input Parameters: 859 + dm - The DMPlex object 860 . pointSF - The PetscSF describing the communication pattern 861 . originalSection - The PetscSection for existing data layout 862 - originalIS - The existing data 863 864 Output Parameters: 865 + newSection - The PetscSF describing the new data layout 866 - newIS - The new data 867 868 Level: developer 869 870 .seealso: DMPlexDistribute(), DMPlexDistributeField(), DMPlexDistributeData() 871 @*/ 872 PetscErrorCode DMPlexDistributeFieldIS(DM dm, PetscSF pointSF, PetscSection originalSection, IS originalIS, PetscSection newSection, IS *newIS) 873 { 874 PetscSF fieldSF; 875 PetscInt *newValues, *remoteOffsets, fieldSize; 876 const PetscInt *originalValues; 877 PetscErrorCode ierr; 878 879 PetscFunctionBegin; 880 ierr = PetscLogEventBegin(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 881 ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); 882 883 ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); 884 ierr = PetscMalloc(fieldSize * sizeof(PetscInt), &newValues);CHKERRQ(ierr); 885 886 ierr = ISGetIndices(originalIS, &originalValues);CHKERRQ(ierr); 887 ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); 888 ierr = PetscSFBcastBegin(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr); 889 ierr = PetscSFBcastEnd(fieldSF, MPIU_INT, (PetscInt *) originalValues, newValues);CHKERRQ(ierr); 890 ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); 891 ierr = ISRestoreIndices(originalIS, &originalValues);CHKERRQ(ierr); 892 ierr = ISCreateGeneral(PetscObjectComm((PetscObject) pointSF), fieldSize, newValues, PETSC_OWN_POINTER, newIS);CHKERRQ(ierr); 893 ierr = PetscLogEventEnd(DMPLEX_DistributeField,dm,0,0,0);CHKERRQ(ierr); 894 PetscFunctionReturn(0); 895 } 896 897 #undef __FUNCT__ 898 #define __FUNCT__ "DMPlexDistributeData" 899 /*@ 900 DMPlexDistributeData - Distribute field data to match a given PetscSF, usually the SF from mesh distribution 901 902 Collective on DM 903 904 Input Parameters: 905 + dm - The DMPlex object 906 . pointSF - The PetscSF describing the communication pattern 907 . originalSection - The PetscSection for existing data layout 908 . datatype - The type of data 909 - originalData - The existing data 910 911 Output Parameters: 912 + newSection - The PetscSection describing the new data layout 913 - newData - The new data 914 915 Level: developer 916 917 .seealso: DMPlexDistribute(), DMPlexDistributeField() 918 @*/ 919 PetscErrorCode DMPlexDistributeData(DM dm, PetscSF pointSF, PetscSection originalSection, MPI_Datatype datatype, void *originalData, PetscSection newSection, void **newData) 920 { 921 PetscSF fieldSF; 922 PetscInt *remoteOffsets, fieldSize; 923 PetscMPIInt dataSize; 924 PetscErrorCode ierr; 925 926 PetscFunctionBegin; 927 ierr = PetscLogEventBegin(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr); 928 ierr = PetscSFDistributeSection(pointSF, originalSection, &remoteOffsets, newSection);CHKERRQ(ierr); 929 930 ierr = PetscSectionGetStorageSize(newSection, &fieldSize);CHKERRQ(ierr); 931 ierr = MPI_Type_size(datatype, &dataSize);CHKERRQ(ierr); 932 ierr = PetscMalloc(fieldSize * dataSize, newData);CHKERRQ(ierr); 933 934 ierr = PetscSFCreateSectionSF(pointSF, originalSection, remoteOffsets, newSection, &fieldSF);CHKERRQ(ierr); 935 ierr = PetscSFBcastBegin(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr); 936 ierr = PetscSFBcastEnd(fieldSF, datatype, originalData, *newData);CHKERRQ(ierr); 937 ierr = PetscSFDestroy(&fieldSF);CHKERRQ(ierr); 938 ierr = PetscLogEventEnd(DMPLEX_DistributeData,dm,0,0,0);CHKERRQ(ierr); 939 PetscFunctionReturn(0); 940 } 941 942 #undef __FUNCT__ 943 #define __FUNCT__ "DMPlexDistributeCones" 944 PetscErrorCode DMPlexDistributeCones(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel) 945 { 946 DM_Plex *mesh = (DM_Plex*) dm->data; 947 DM_Plex *pmesh = (DM_Plex*) (dmParallel)->data; 948 MPI_Comm comm; 949 PetscSF coneSF; 950 PetscSection originalConeSection, newConeSection; 951 PetscInt *remoteOffsets, *cones, *newCones, newConesSize; 952 PetscBool flg; 953 PetscErrorCode ierr; 954 955 PetscFunctionBegin; 956 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 957 PetscValidPointer(dmParallel,4); 958 ierr = PetscLogEventBegin(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr); 959 960 /* Distribute cone section */ 961 ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 962 ierr = DMPlexGetConeSection(dm, &originalConeSection);CHKERRQ(ierr); 963 ierr = DMPlexGetConeSection(dmParallel, &newConeSection);CHKERRQ(ierr); 964 ierr = PetscSFDistributeSection(migrationSF, originalConeSection, &remoteOffsets, newConeSection);CHKERRQ(ierr); 965 ierr = DMSetUp(dmParallel);CHKERRQ(ierr); 966 { 967 PetscInt pStart, pEnd, p; 968 969 ierr = PetscSectionGetChart(newConeSection, &pStart, &pEnd);CHKERRQ(ierr); 970 for (p = pStart; p < pEnd; ++p) { 971 PetscInt coneSize; 972 ierr = PetscSectionGetDof(newConeSection, p, &coneSize);CHKERRQ(ierr); 973 pmesh->maxConeSize = PetscMax(pmesh->maxConeSize, coneSize); 974 } 975 } 976 /* Communicate and renumber cones */ 977 ierr = PetscSFCreateSectionSF(migrationSF, originalConeSection, remoteOffsets, newConeSection, &coneSF);CHKERRQ(ierr); 978 ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr); 979 ierr = DMPlexGetCones(dmParallel, &newCones);CHKERRQ(ierr); 980 ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 981 ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 982 ierr = PetscSectionGetStorageSize(newConeSection, &newConesSize);CHKERRQ(ierr); 983 ierr = ISGlobalToLocalMappingApplyBlock(renumbering, IS_GTOLM_MASK, newConesSize, newCones, NULL, newCones);CHKERRQ(ierr); 984 #if PETSC_USE_DEBUG 985 { 986 PetscInt p; 987 PetscBool valid = PETSC_TRUE; 988 for (p = 0; p < newConesSize; ++p) { 989 if (newCones[p] < 0) {valid = PETSC_FALSE; ierr = PetscPrintf(PETSC_COMM_SELF, "Point %d not in overlap SF\n", p);CHKERRQ(ierr);} 990 } 991 if (!valid) SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid global to local map"); 992 } 993 #endif 994 ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-cones_view", &flg);CHKERRQ(ierr); 995 if (flg) { 996 ierr = PetscPrintf(comm, "Serial Cone Section:\n");CHKERRQ(ierr); 997 ierr = PetscSectionView(originalConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 998 ierr = PetscPrintf(comm, "Parallel Cone Section:\n");CHKERRQ(ierr); 999 ierr = PetscSectionView(newConeSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1000 ierr = PetscSFView(coneSF, NULL);CHKERRQ(ierr); 1001 } 1002 ierr = DMPlexGetConeOrientations(dm, &cones);CHKERRQ(ierr); 1003 ierr = DMPlexGetConeOrientations(dmParallel, &newCones);CHKERRQ(ierr); 1004 ierr = PetscSFBcastBegin(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 1005 ierr = PetscSFBcastEnd(coneSF, MPIU_INT, cones, newCones);CHKERRQ(ierr); 1006 ierr = PetscSFDestroy(&coneSF);CHKERRQ(ierr); 1007 ierr = PetscLogEventEnd(DMPLEX_DistributeCones,dm,0,0,0);CHKERRQ(ierr); 1008 /* Create supports and stratify sieve */ 1009 { 1010 PetscInt pStart, pEnd; 1011 1012 ierr = PetscSectionGetChart(pmesh->coneSection, &pStart, &pEnd);CHKERRQ(ierr); 1013 ierr = PetscSectionSetChart(pmesh->supportSection, pStart, pEnd);CHKERRQ(ierr); 1014 } 1015 ierr = DMPlexSymmetrize(dmParallel);CHKERRQ(ierr); 1016 ierr = DMPlexStratify(dmParallel);CHKERRQ(ierr); 1017 pmesh->useCone = mesh->useCone; 1018 pmesh->useClosure = mesh->useClosure; 1019 PetscFunctionReturn(0); 1020 } 1021 1022 #undef __FUNCT__ 1023 #define __FUNCT__ "DMPlexDistributeCoordinates" 1024 PetscErrorCode DMPlexDistributeCoordinates(DM dm, PetscSF migrationSF, DM dmParallel) 1025 { 1026 MPI_Comm comm; 1027 PetscSection originalCoordSection, newCoordSection; 1028 Vec originalCoordinates, newCoordinates; 1029 PetscInt bs; 1030 const char *name; 1031 const PetscReal *maxCell, *L; 1032 PetscErrorCode ierr; 1033 1034 PetscFunctionBegin; 1035 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1036 PetscValidPointer(dmParallel, 3); 1037 1038 ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 1039 ierr = DMGetCoordinateSection(dm, &originalCoordSection);CHKERRQ(ierr); 1040 ierr = DMGetCoordinateSection(dmParallel, &newCoordSection);CHKERRQ(ierr); 1041 ierr = DMGetCoordinatesLocal(dm, &originalCoordinates);CHKERRQ(ierr); 1042 if (originalCoordinates) { 1043 ierr = VecCreate(comm, &newCoordinates);CHKERRQ(ierr); 1044 ierr = PetscObjectGetName((PetscObject) originalCoordinates, &name);CHKERRQ(ierr); 1045 ierr = PetscObjectSetName((PetscObject) newCoordinates, name);CHKERRQ(ierr); 1046 1047 ierr = DMPlexDistributeField(dm, migrationSF, originalCoordSection, originalCoordinates, newCoordSection, newCoordinates);CHKERRQ(ierr); 1048 ierr = DMSetCoordinatesLocal(dmParallel, newCoordinates);CHKERRQ(ierr); 1049 ierr = VecGetBlockSize(originalCoordinates, &bs);CHKERRQ(ierr); 1050 ierr = VecSetBlockSize(newCoordinates, bs);CHKERRQ(ierr); 1051 ierr = VecDestroy(&newCoordinates);CHKERRQ(ierr); 1052 } 1053 ierr = DMGetPeriodicity(dm, &maxCell, &L);CHKERRQ(ierr); 1054 if (L) {ierr = DMSetPeriodicity(dmParallel, maxCell, L);CHKERRQ(ierr);} 1055 PetscFunctionReturn(0); 1056 } 1057 1058 #undef __FUNCT__ 1059 #define __FUNCT__ "DMPlexDistributeLabels" 1060 /* Here we are assuming that process 0 always has everything */ 1061 PetscErrorCode DMPlexDistributeLabels(DM dm, PetscSF migrationSF, DM dmParallel) 1062 { 1063 MPI_Comm comm; 1064 PetscMPIInt rank; 1065 PetscInt numLabels, numLocalLabels, l; 1066 PetscBool hasLabels = PETSC_FALSE; 1067 PetscErrorCode ierr; 1068 1069 PetscFunctionBegin; 1070 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1071 PetscValidHeaderSpecific(dm, DM_CLASSID, 3); 1072 ierr = PetscLogEventBegin(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr); 1073 ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 1074 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1075 1076 /* Everyone must have either the same number of labels, or none */ 1077 ierr = DMPlexGetNumLabels(dm, &numLocalLabels);CHKERRQ(ierr); 1078 numLabels = numLocalLabels; 1079 ierr = MPI_Bcast(&numLabels, 1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1080 if (numLabels == numLocalLabels) hasLabels = PETSC_TRUE; 1081 for (l = numLabels-1; l >= 0; --l) { 1082 DMLabel label = NULL, labelNew = NULL; 1083 PetscBool isdepth; 1084 1085 if (hasLabels) { 1086 ierr = DMPlexGetLabelByNum(dm, l, &label);CHKERRQ(ierr); 1087 /* Skip "depth" because it is recreated */ 1088 ierr = PetscStrcmp(label->name, "depth", &isdepth);CHKERRQ(ierr); 1089 } 1090 ierr = MPI_Bcast(&isdepth, 1, MPIU_BOOL, 0, comm);CHKERRQ(ierr); 1091 if (isdepth) continue; 1092 ierr = DMLabelDistribute(label, migrationSF, &labelNew);CHKERRQ(ierr); 1093 ierr = DMPlexAddLabel(dmParallel, labelNew);CHKERRQ(ierr); 1094 } 1095 ierr = PetscLogEventEnd(DMPLEX_DistributeLabels,dm,0,0,0);CHKERRQ(ierr); 1096 PetscFunctionReturn(0); 1097 } 1098 1099 #undef __FUNCT__ 1100 #define __FUNCT__ "DMPlexDistributeSetupHybrid" 1101 PetscErrorCode DMPlexDistributeSetupHybrid(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel) 1102 { 1103 DM_Plex *mesh = (DM_Plex*) dm->data; 1104 DM_Plex *pmesh = (DM_Plex*) (dmParallel)->data; 1105 MPI_Comm comm; 1106 const PetscInt *gpoints; 1107 PetscInt dim, depth, n, d; 1108 PetscErrorCode ierr; 1109 1110 PetscFunctionBegin; 1111 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1112 PetscValidPointer(dmParallel, 4); 1113 1114 ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 1115 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1116 1117 /* Setup hybrid structure */ 1118 for (d = 0; d <= dim; ++d) {pmesh->hybridPointMax[d] = mesh->hybridPointMax[d];} 1119 ierr = MPI_Bcast(pmesh->hybridPointMax, dim+1, MPIU_INT, 0, comm);CHKERRQ(ierr); 1120 ierr = ISLocalToGlobalMappingGetSize(renumbering, &n);CHKERRQ(ierr); 1121 ierr = ISLocalToGlobalMappingGetIndices(renumbering, &gpoints);CHKERRQ(ierr); 1122 ierr = DMPlexGetDepth(dm, &depth);CHKERRQ(ierr); 1123 for (d = 0; d <= dim; ++d) { 1124 PetscInt pmax = pmesh->hybridPointMax[d], newmax = 0, pEnd, stratum[2], p; 1125 1126 if (pmax < 0) continue; 1127 ierr = DMPlexGetDepthStratum(dm, d > depth ? depth : d, &stratum[0], &stratum[1]);CHKERRQ(ierr); 1128 ierr = DMPlexGetDepthStratum(dmParallel, d, NULL, &pEnd);CHKERRQ(ierr); 1129 ierr = MPI_Bcast(stratum, 2, MPIU_INT, 0, comm);CHKERRQ(ierr); 1130 for (p = 0; p < n; ++p) { 1131 const PetscInt point = gpoints[p]; 1132 1133 if ((point >= stratum[0]) && (point < stratum[1]) && (point >= pmax)) ++newmax; 1134 } 1135 if (newmax > 0) pmesh->hybridPointMax[d] = pEnd - newmax; 1136 else pmesh->hybridPointMax[d] = -1; 1137 } 1138 ierr = ISLocalToGlobalMappingRestoreIndices(renumbering, &gpoints);CHKERRQ(ierr); 1139 PetscFunctionReturn(0); 1140 } 1141 1142 #undef __FUNCT__ 1143 #define __FUNCT__ "DMPlexDistributeSetupTree" 1144 PetscErrorCode DMPlexDistributeSetupTree(DM dm, PetscSF migrationSF, ISLocalToGlobalMapping renumbering, DM dmParallel) 1145 { 1146 DM_Plex *mesh = (DM_Plex*) dm->data; 1147 DM_Plex *pmesh = (DM_Plex*) (dmParallel)->data; 1148 MPI_Comm comm; 1149 DM refTree; 1150 PetscSection origParentSection, newParentSection; 1151 PetscInt *origParents, *origChildIDs; 1152 PetscBool flg; 1153 PetscErrorCode ierr; 1154 1155 PetscFunctionBegin; 1156 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1157 PetscValidHeaderSpecific(dm, DM_CLASSID, 4); 1158 ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 1159 1160 /* Set up tree */ 1161 ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr); 1162 ierr = DMPlexSetReferenceTree(dmParallel,refTree);CHKERRQ(ierr); 1163 ierr = DMPlexGetTree(dm,&origParentSection,&origParents,&origChildIDs,NULL,NULL);CHKERRQ(ierr); 1164 if (origParentSection) { 1165 PetscInt pStart, pEnd; 1166 PetscInt *newParents, *newChildIDs; 1167 PetscInt *remoteOffsetsParents, newParentSize; 1168 PetscSF parentSF; 1169 1170 ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr); 1171 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dmParallel),&newParentSection);CHKERRQ(ierr); 1172 ierr = PetscSectionSetChart(newParentSection,pStart,pEnd);CHKERRQ(ierr); 1173 ierr = PetscSFDistributeSection(migrationSF, origParentSection, &remoteOffsetsParents, newParentSection);CHKERRQ(ierr); 1174 ierr = PetscSFCreateSectionSF(migrationSF, origParentSection, remoteOffsetsParents, newParentSection, &parentSF);CHKERRQ(ierr); 1175 ierr = PetscSectionGetStorageSize(newParentSection,&newParentSize);CHKERRQ(ierr); 1176 ierr = PetscMalloc2(newParentSize,&newParents,newParentSize,&newChildIDs);CHKERRQ(ierr); 1177 ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr); 1178 ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origParents, newParents);CHKERRQ(ierr); 1179 ierr = PetscSFBcastBegin(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr); 1180 ierr = PetscSFBcastEnd(parentSF, MPIU_INT, origChildIDs, newChildIDs);CHKERRQ(ierr); 1181 ierr = ISGlobalToLocalMappingApplyBlock(renumbering,IS_GTOLM_MASK, newParentSize, newParents, NULL, newParents);CHKERRQ(ierr); 1182 ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-parents_view", &flg);CHKERRQ(ierr); 1183 if (flg) { 1184 ierr = PetscPrintf(comm, "Serial Parent Section: \n");CHKERRQ(ierr); 1185 ierr = PetscSectionView(origParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1186 ierr = PetscPrintf(comm, "Parallel Parent Section: \n");CHKERRQ(ierr); 1187 ierr = PetscSectionView(newParentSection, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1188 ierr = PetscSFView(parentSF, NULL);CHKERRQ(ierr); 1189 } 1190 ierr = DMPlexSetTree(dmParallel,newParentSection,newParents,newChildIDs);CHKERRQ(ierr); 1191 ierr = PetscSectionDestroy(&newParentSection);CHKERRQ(ierr); 1192 ierr = PetscFree2(newParents,newChildIDs);CHKERRQ(ierr); 1193 ierr = PetscSFDestroy(&parentSF);CHKERRQ(ierr); 1194 } 1195 pmesh->useAnchors = mesh->useAnchors; 1196 PetscFunctionReturn(0); 1197 } 1198 1199 #undef __FUNCT__ 1200 #define __FUNCT__ "DMPlexDistributeSF" 1201 PetscErrorCode DMPlexDistributeSF(DM dm, PetscSF migrationSF, DM dmParallel) 1202 { 1203 DM_Plex *mesh = (DM_Plex*) dm->data; 1204 DM_Plex *pmesh = (DM_Plex*) (dmParallel)->data; 1205 PetscMPIInt rank, numProcs; 1206 MPI_Comm comm; 1207 PetscErrorCode ierr; 1208 1209 PetscFunctionBegin; 1210 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1211 PetscValidPointer(dmParallel,7); 1212 1213 /* Create point SF for parallel mesh */ 1214 ierr = PetscLogEventBegin(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr); 1215 ierr = PetscObjectGetComm((PetscObject)dm, &comm);CHKERRQ(ierr); 1216 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1217 ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 1218 { 1219 const PetscInt *leaves; 1220 PetscSFNode *remotePoints, *rowners, *lowners; 1221 PetscInt numRoots, numLeaves, numGhostPoints = 0, p, gp, *ghostPoints; 1222 PetscInt pStart, pEnd; 1223 1224 ierr = DMPlexGetChart(dmParallel, &pStart, &pEnd);CHKERRQ(ierr); 1225 ierr = PetscSFGetGraph(migrationSF, &numRoots, &numLeaves, &leaves, NULL);CHKERRQ(ierr); 1226 ierr = PetscMalloc2(numRoots,&rowners,numLeaves,&lowners);CHKERRQ(ierr); 1227 for (p=0; p<numRoots; p++) { 1228 rowners[p].rank = -1; 1229 rowners[p].index = -1; 1230 } 1231 ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 1232 ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 1233 for (p = 0; p < numLeaves; ++p) { 1234 if (lowners[p].rank < 0 || lowners[p].rank == rank) { /* Either put in a bid or we know we own it */ 1235 lowners[p].rank = rank; 1236 lowners[p].index = leaves ? leaves[p] : p; 1237 } else if (lowners[p].rank >= 0) { /* Point already claimed so flag so that MAXLOC does not listen to us */ 1238 lowners[p].rank = -2; 1239 lowners[p].index = -2; 1240 } 1241 } 1242 for (p=0; p<numRoots; p++) { /* Root must not participate in the rediction, flag so that MAXLOC does not use */ 1243 rowners[p].rank = -3; 1244 rowners[p].index = -3; 1245 } 1246 ierr = PetscSFReduceBegin(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr); 1247 ierr = PetscSFReduceEnd(migrationSF, MPIU_2INT, lowners, rowners, MPI_MAXLOC);CHKERRQ(ierr); 1248 ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 1249 ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rowners, lowners);CHKERRQ(ierr); 1250 for (p = 0; p < numLeaves; ++p) { 1251 if (lowners[p].rank < 0 || lowners[p].index < 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Cell partition corrupt: point not claimed"); 1252 if (lowners[p].rank != rank) ++numGhostPoints; 1253 } 1254 ierr = PetscMalloc1(numGhostPoints, &ghostPoints);CHKERRQ(ierr); 1255 ierr = PetscMalloc1(numGhostPoints, &remotePoints);CHKERRQ(ierr); 1256 for (p = 0, gp = 0; p < numLeaves; ++p) { 1257 if (lowners[p].rank != rank) { 1258 ghostPoints[gp] = leaves ? leaves[p] : p; 1259 remotePoints[gp].rank = lowners[p].rank; 1260 remotePoints[gp].index = lowners[p].index; 1261 ++gp; 1262 } 1263 } 1264 ierr = PetscFree2(rowners,lowners);CHKERRQ(ierr); 1265 ierr = PetscSFSetGraph((dmParallel)->sf, pEnd - pStart, numGhostPoints, ghostPoints, PETSC_OWN_POINTER, remotePoints, PETSC_OWN_POINTER);CHKERRQ(ierr); 1266 ierr = PetscSFSetFromOptions((dmParallel)->sf);CHKERRQ(ierr); 1267 } 1268 pmesh->useCone = mesh->useCone; 1269 pmesh->useClosure = mesh->useClosure; 1270 ierr = PetscLogEventEnd(DMPLEX_DistributeSF,dm,0,0,0);CHKERRQ(ierr); 1271 PetscFunctionReturn(0); 1272 } 1273 1274 #undef __FUNCT__ 1275 #define __FUNCT__ "DMPlexCreatePointSF" 1276 /*@C 1277 DMPlexDerivePointSF - Build a point SF from an SF describing a point migration 1278 1279 Input Parameter: 1280 + dm - The source DMPlex object 1281 . migrationSF - The star forest that describes the parallel point remapping 1282 1283 Output Parameter: 1284 - pointSF - The star forest describing the point overlap in the remapped DM 1285 1286 Level: developer 1287 1288 .seealso: DMPlexDistribute(), DMPlexDistributeOverlap() 1289 @*/ 1290 PetscErrorCode DMPlexCreatePointSF(DM dm, PetscSF migrationSF, PetscBool ownership, PetscSF *pointSF) 1291 { 1292 PetscMPIInt rank; 1293 PetscInt p, pStart, pEnd, nroots, nleaves, idx, npointLeaves; 1294 PetscInt *pointLocal; 1295 const PetscInt *leaves; 1296 const PetscSFNode *roots; 1297 PetscSFNode *rootNodes, *leafNodes, *pointRemote; 1298 PetscErrorCode ierr; 1299 1300 PetscFunctionBegin; 1301 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1302 ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRQ(ierr); 1303 1304 ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 1305 ierr = PetscSFGetGraph(migrationSF, &nroots, &nleaves, &leaves, &roots);CHKERRQ(ierr); 1306 ierr = PetscMalloc2(nroots, &rootNodes, nleaves, &leafNodes);CHKERRQ(ierr); 1307 if (ownership) { 1308 for (p = 0; p < nleaves; ++p) { 1309 /* Either put in a bid or we know we own it */ 1310 leafNodes[p].rank = rank; 1311 leafNodes[p].index = p; 1312 } 1313 for (p = 0; p < nroots; p++) { 1314 /* Root must not participate in the rediction, flag so that MAXLOC does not use */ 1315 rootNodes[p].rank = -3; 1316 rootNodes[p].index = -3; 1317 } 1318 ierr = PetscSFReduceBegin(migrationSF, MPIU_2INT, leafNodes, rootNodes, MPI_MAXLOC);CHKERRQ(ierr); 1319 ierr = PetscSFReduceEnd(migrationSF, MPIU_2INT, leafNodes, rootNodes, MPI_MAXLOC);CHKERRQ(ierr); 1320 } else { 1321 for (p = 0; p < nroots; p++) { 1322 rootNodes[p].index = -1; 1323 rootNodes[p].rank = rank; 1324 }; 1325 for (p = 0; p < nleaves; p++) { 1326 /* Write new local id into old location */ 1327 if (roots[p].rank == rank) { 1328 rootNodes[roots[p].index].index = leaves[p]; 1329 } 1330 } 1331 } 1332 ierr = PetscSFBcastBegin(migrationSF, MPIU_2INT, rootNodes, leafNodes);CHKERRQ(ierr); 1333 ierr = PetscSFBcastEnd(migrationSF, MPIU_2INT, rootNodes, leafNodes);CHKERRQ(ierr); 1334 1335 for (npointLeaves = 0, p = 0; p < nleaves; p++) {if (leafNodes[p].rank != rank) npointLeaves++;} 1336 ierr = PetscMalloc2(npointLeaves, &pointLocal, npointLeaves, &pointRemote);CHKERRQ(ierr); 1337 for (idx = 0, p = 0; p < nleaves; p++) { 1338 if (leafNodes[p].rank != rank) { 1339 pointLocal[idx] = p; 1340 pointRemote[idx] = leafNodes[p]; 1341 idx++; 1342 } 1343 } 1344 ierr = PetscSFCreate(PetscObjectComm((PetscObject) dm), pointSF);CHKERRQ(ierr); 1345 ierr = PetscObjectSetName((PetscObject) *pointSF, "Point SF");CHKERRQ(ierr); 1346 ierr = PetscSFSetGraph(*pointSF, nleaves, npointLeaves, pointLocal, PETSC_OWN_POINTER, pointRemote, PETSC_OWN_POINTER);CHKERRQ(ierr); 1347 ierr = PetscFree2(rootNodes, leafNodes);CHKERRQ(ierr); 1348 PetscFunctionReturn(0); 1349 } 1350 1351 #undef __FUNCT__ 1352 #define __FUNCT__ "DMPlexMigrate" 1353 /*@C 1354 DMPlexMigrate - Migrates internal DM data over the supplied star forest 1355 1356 Input Parameter: 1357 + dm - The source DMPlex object 1358 . sf - The start forest communication context describing the migration pattern 1359 1360 Output Parameter: 1361 - targetDM - The target DMPlex object 1362 1363 Level: intermediate 1364 1365 .seealso: DMPlexDistribute(), DMPlexDistributeOverlap() 1366 @*/ 1367 PetscErrorCode DMPlexMigrate(DM dm, PetscSF sf, DM targetDM) 1368 { 1369 MPI_Comm comm; 1370 PetscInt dim, nroots; 1371 PetscSF sfPoint; 1372 ISLocalToGlobalMapping ltogMigration; 1373 PetscBool flg; 1374 PetscErrorCode ierr; 1375 1376 PetscFunctionBegin; 1377 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1378 ierr = PetscObjectGetComm((PetscObject) dm, &comm);CHKERRQ(ierr); 1379 ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 1380 ierr = DMSetDimension(targetDM, dim);CHKERRQ(ierr); 1381 1382 ierr = DMGetPointSF(dm, &sfPoint);CHKERRQ(ierr); 1383 ierr = PetscSFGetGraph(sfPoint, &nroots, NULL, NULL, NULL);CHKERRQ(ierr); 1384 if (nroots > 0) { 1385 IS isOriginal; 1386 ISLocalToGlobalMapping ltogOriginal; 1387 PetscInt n, size, nleaves, conesSize; 1388 PetscInt *numbering_orig, *numbering_new, *cones; 1389 PetscSection coneSection; 1390 /* Get the original point numbering */ 1391 ierr = DMPlexCreatePointNumbering(dm, &isOriginal);CHKERRQ(ierr); 1392 ierr = ISLocalToGlobalMappingCreateIS(isOriginal, <ogOriginal);CHKERRQ(ierr); 1393 ierr = ISLocalToGlobalMappingGetSize(ltogOriginal, &size);CHKERRQ(ierr); 1394 ierr = ISLocalToGlobalMappingGetBlockIndices(ltogOriginal, (const PetscInt**)&numbering_orig);CHKERRQ(ierr); 1395 /* Convert to positive global numbers */ 1396 for (n=0; n<size; n++) {if (numbering_orig[n] < 0) numbering_orig[n] = -(numbering_orig[n]+1);} 1397 /* Derive the new local-to-global mapping from the old one */ 1398 ierr = PetscSFGetGraph(sf, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr); 1399 ierr = PetscMalloc1(nleaves, &numbering_new);CHKERRQ(ierr); 1400 ierr = PetscSFBcastBegin(sf, MPIU_INT, (PetscInt *) numbering_orig, numbering_new);CHKERRQ(ierr); 1401 ierr = PetscSFBcastEnd(sf, MPIU_INT, (PetscInt *) numbering_orig, numbering_new);CHKERRQ(ierr); 1402 ierr = ISLocalToGlobalMappingCreate(comm, 1, nleaves, (const PetscInt*) numbering_new, PETSC_OWN_POINTER, <ogMigration);CHKERRQ(ierr); 1403 ierr = ISLocalToGlobalMappingRestoreIndices(ltogOriginal, (const PetscInt**)&numbering_orig);CHKERRQ(ierr); 1404 /* Convert cones to global numbering before migrating them */ 1405 ierr = DMPlexGetConeSection(dm, &coneSection);CHKERRQ(ierr); 1406 ierr = PetscSectionGetStorageSize(coneSection, &conesSize);CHKERRQ(ierr); 1407 ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr); 1408 ierr = ISLocalToGlobalMappingApplyBlock(ltogOriginal, conesSize, cones, cones);CHKERRQ(ierr); 1409 ierr = ISDestroy(&isOriginal);CHKERRQ(ierr); 1410 ierr = ISLocalToGlobalMappingDestroy(<ogOriginal); 1411 } else { 1412 /* Assume zero-to-all-ranks pattern and derive LToG from SF */ 1413 ierr = ISLocalToGlobalMappingCreateSF(sf, 0, <ogMigration);CHKERRQ(ierr); 1414 } 1415 ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr); 1416 if (flg) { 1417 ierr = PetscPrintf(comm, "Point renumbering for DM migration:\n");CHKERRQ(ierr); 1418 ierr = ISLocalToGlobalMappingView(ltogMigration, NULL);CHKERRQ(ierr); 1419 } 1420 /* Migrate DM data to target DM */ 1421 ierr = DMPlexDistributeCones(dm, sf, ltogMigration, targetDM);CHKERRQ(ierr); 1422 ierr = DMPlexDistributeCoordinates(dm, sf, targetDM);CHKERRQ(ierr); 1423 ierr = DMPlexDistributeLabels(dm, sf, targetDM);CHKERRQ(ierr); 1424 ierr = DMPlexDistributeSetupHybrid(dm, sf, ltogMigration, targetDM);CHKERRQ(ierr); 1425 ierr = DMPlexDistributeSetupTree(dm, sf, ltogMigration, targetDM);CHKERRQ(ierr); 1426 ierr = DMPlexCopyBoundary(dm, targetDM);CHKERRQ(ierr); 1427 ierr = ISLocalToGlobalMappingDestroy(<ogMigration); 1428 PetscFunctionReturn(0); 1429 } 1430 1431 #undef __FUNCT__ 1432 #define __FUNCT__ "DMPlexDistribute" 1433 /*@C 1434 DMPlexDistribute - Distributes the mesh and any associated sections. 1435 1436 Not Collective 1437 1438 Input Parameter: 1439 + dm - The original DMPlex object 1440 - overlap - The overlap of partitions, 0 is the default 1441 1442 Output Parameter: 1443 + sf - The PetscSF used for point distribution 1444 - parallelMesh - The distributed DMPlex object, or NULL 1445 1446 Note: If the mesh was not distributed, the return value is NULL. 1447 1448 The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and 1449 DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function 1450 representation on the mesh. 1451 1452 Level: intermediate 1453 1454 .keywords: mesh, elements 1455 .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure() 1456 @*/ 1457 PetscErrorCode DMPlexDistribute(DM dm, PetscInt overlap, PetscSF *sf, DM *dmParallel) 1458 { 1459 MPI_Comm comm; 1460 PetscPartitioner partitioner; 1461 IS cellPart; 1462 PetscSection cellPartSection; 1463 DMLabel lblPartition, lblMigration; 1464 PetscSF sfProcess, sfMigration, sfStratified, sfPoint; 1465 PetscBool flg; 1466 PetscMPIInt rank, numProcs, p; 1467 PetscErrorCode ierr; 1468 1469 PetscFunctionBegin; 1470 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1471 if (sf) PetscValidPointer(sf,4); 1472 PetscValidPointer(dmParallel,5); 1473 1474 ierr = PetscLogEventBegin(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr); 1475 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1476 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1477 ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr); 1478 1479 *dmParallel = NULL; 1480 if (numProcs == 1) PetscFunctionReturn(0); 1481 1482 /* Create cell partition */ 1483 ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr); 1484 ierr = PetscSectionCreate(comm, &cellPartSection);CHKERRQ(ierr); 1485 ierr = DMPlexGetPartitioner(dm, &partitioner);CHKERRQ(ierr); 1486 ierr = PetscPartitionerPartition(partitioner, dm, cellPartSection, &cellPart);CHKERRQ(ierr); 1487 { 1488 /* Convert partition to DMLabel */ 1489 PetscInt proc, pStart, pEnd, npoints, poffset; 1490 const PetscInt *points; 1491 ierr = DMLabelCreate("Point Partition", &lblPartition);CHKERRQ(ierr); 1492 ierr = ISGetIndices(cellPart, &points);CHKERRQ(ierr); 1493 ierr = PetscSectionGetChart(cellPartSection, &pStart, &pEnd);CHKERRQ(ierr); 1494 for (proc = pStart; proc < pEnd; proc++) { 1495 ierr = PetscSectionGetDof(cellPartSection, proc, &npoints);CHKERRQ(ierr); 1496 ierr = PetscSectionGetOffset(cellPartSection, proc, &poffset);CHKERRQ(ierr); 1497 for (p = poffset; p < poffset+npoints; p++) { 1498 ierr = DMLabelSetValue(lblPartition, points[p], proc);CHKERRQ(ierr); 1499 } 1500 } 1501 ierr = ISRestoreIndices(cellPart, &points);CHKERRQ(ierr); 1502 } 1503 ierr = DMPlexPartitionLabelClosure(dm, lblPartition);CHKERRQ(ierr); 1504 { 1505 /* Build a global process SF */ 1506 PetscSFNode *remoteProc; 1507 ierr = PetscMalloc1(numProcs, &remoteProc);CHKERRQ(ierr); 1508 for (p = 0; p < numProcs; ++p) { 1509 remoteProc[p].rank = p; 1510 remoteProc[p].index = rank; 1511 } 1512 ierr = PetscSFCreate(comm, &sfProcess);CHKERRQ(ierr); 1513 ierr = PetscObjectSetName((PetscObject) sfProcess, "Process SF");CHKERRQ(ierr); 1514 ierr = PetscSFSetGraph(sfProcess, numProcs, numProcs, NULL, PETSC_OWN_POINTER, remoteProc, PETSC_OWN_POINTER);CHKERRQ(ierr); 1515 } 1516 ierr = DMLabelCreate("Point migration", &lblMigration);CHKERRQ(ierr); 1517 ierr = DMPlexPartitionLabelInvert(dm, lblPartition, sfProcess, lblMigration);CHKERRQ(ierr); 1518 ierr = DMPlexPartitionLabelCreateSF(dm, lblMigration, &sfMigration); 1519 /* Stratify the SF in case we are migrating an already parallel plex */ 1520 ierr = DMPlexStratifyMigrationSF(dm, sfMigration, &sfStratified);CHKERRQ(ierr); 1521 ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr); 1522 sfMigration = sfStratified; 1523 ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr); 1524 ierr = PetscOptionsHasName(((PetscObject) dm)->prefix, "-partition_view", &flg);CHKERRQ(ierr); 1525 if (flg) { 1526 ierr = PetscViewerASCIISynchronizedAllow(PETSC_VIEWER_STDOUT_WORLD, PETSC_TRUE);CHKERRQ(ierr); 1527 ierr = DMLabelView(lblPartition, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1528 ierr = PetscSFView(sfMigration, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 1529 } 1530 1531 /* Create non-overlapping parallel DM and migrate internal data */ 1532 ierr = DMPlexCreate(comm, dmParallel);CHKERRQ(ierr); 1533 ierr = PetscObjectSetName((PetscObject) *dmParallel, "Parallel Mesh");CHKERRQ(ierr); 1534 ierr = DMPlexMigrate(dm, sfMigration, *dmParallel);CHKERRQ(ierr); 1535 1536 /* Build the point SF without overlap */ 1537 ierr = DMPlexCreatePointSF(dm, sfMigration, PETSC_TRUE, &sfPoint);CHKERRQ(ierr); 1538 ierr = DMSetPointSF(*dmParallel, sfPoint);CHKERRQ(ierr); 1539 if (flg) {ierr = PetscSFView(sfPoint, PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);} 1540 1541 if (overlap > 0) { 1542 DM dmOverlap; 1543 PetscInt nroots, nleaves; 1544 PetscSFNode *newRemote; 1545 const PetscSFNode *oldRemote; 1546 PetscSF sfOverlap, sfOverlapPoint; 1547 ierr = PetscLogEventBegin(DMPLEX_DistributeOverlap,dm,0,0,0);CHKERRQ(ierr); 1548 /* Add the partition overlap to the distributed DM */ 1549 ierr = DMPlexDistributeOverlap(*dmParallel, overlap, &sfOverlap, &dmOverlap);CHKERRQ(ierr); 1550 ierr = DMDestroy(dmParallel);CHKERRQ(ierr); 1551 *dmParallel = dmOverlap; 1552 if (flg) { 1553 ierr = PetscPrintf(comm, "Overlap Migration SF:\n");CHKERRQ(ierr); 1554 ierr = PetscSFView(sfOverlap, NULL);CHKERRQ(ierr); 1555 } 1556 1557 /* Re-map the migration SF to establish the full migration pattern */ 1558 ierr = PetscSFGetGraph(sfMigration, &nroots, NULL, NULL, &oldRemote);CHKERRQ(ierr); 1559 ierr = PetscSFGetGraph(sfOverlap, NULL, &nleaves, NULL, NULL);CHKERRQ(ierr); 1560 ierr = PetscMalloc1(nleaves, &newRemote);CHKERRQ(ierr); 1561 ierr = PetscSFBcastBegin(sfOverlap, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr); 1562 ierr = PetscSFBcastEnd(sfOverlap, MPIU_2INT, oldRemote, newRemote);CHKERRQ(ierr); 1563 ierr = PetscSFCreate(comm, &sfOverlapPoint);CHKERRQ(ierr); 1564 ierr = PetscSFSetGraph(sfOverlapPoint, nroots, nleaves, NULL, PETSC_OWN_POINTER, newRemote, PETSC_OWN_POINTER);CHKERRQ(ierr); 1565 ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr); 1566 ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr); 1567 sfMigration = sfOverlapPoint; 1568 ierr = PetscLogEventEnd(DMPLEX_DistributeOverlap,dm,0,0,0);CHKERRQ(ierr); 1569 } 1570 /* Cleanup Partition */ 1571 ierr = PetscSFDestroy(&sfProcess);CHKERRQ(ierr); 1572 ierr = DMLabelDestroy(&lblPartition);CHKERRQ(ierr); 1573 ierr = DMLabelDestroy(&lblMigration);CHKERRQ(ierr); 1574 ierr = PetscSectionDestroy(&cellPartSection);CHKERRQ(ierr); 1575 ierr = ISDestroy(&cellPart);CHKERRQ(ierr); 1576 if (sf) {*sf = sfMigration;} 1577 else {ierr = PetscSFDestroy(&sfMigration);CHKERRQ(ierr);} 1578 ierr = PetscSFDestroy(&sfPoint);CHKERRQ(ierr); 1579 ierr = DMSetFromOptions(*dmParallel);CHKERRQ(ierr); 1580 ierr = PetscLogEventEnd(DMPLEX_Distribute,dm,0,0,0);CHKERRQ(ierr); 1581 PetscFunctionReturn(0); 1582 } 1583 1584 #undef __FUNCT__ 1585 #define __FUNCT__ "DMPlexDistributeOverlap" 1586 /*@C 1587 DMPlexDistribute - Add partition overlap to a distributed non-overlapping DM. 1588 1589 Not Collective 1590 1591 Input Parameter: 1592 + dm - The non-overlapping distrbuted DMPlex object 1593 - overlap - The overlap of partitions, 0 is the default 1594 1595 Output Parameter: 1596 + sf - The PetscSF used for point distribution 1597 - dmOverlap - The overlapping distributed DMPlex object, or NULL 1598 1599 Note: If the mesh was not distributed, the return value is NULL. 1600 1601 The user can control the definition of adjacency for the mesh using DMPlexGetAdjacencyUseCone() and 1602 DMPlexSetAdjacencyUseClosure(). They should choose the combination appropriate for the function 1603 representation on the mesh. 1604 1605 Level: intermediate 1606 1607 .keywords: mesh, elements 1608 .seealso: DMPlexCreate(), DMPlexDistributeByFace(), DMPlexSetAdjacencyUseCone(), DMPlexSetAdjacencyUseClosure() 1609 @*/ 1610 PetscErrorCode DMPlexDistributeOverlap(DM dm, PetscInt overlap, PetscSF *sf, DM *dmOverlap) 1611 { 1612 MPI_Comm comm; 1613 PetscMPIInt rank; 1614 PetscSection rootSection, leafSection; 1615 IS rootrank, leafrank; 1616 DMLabel lblOverlap; 1617 PetscSF sfOverlap, sfStratified, sfPoint; 1618 PetscErrorCode ierr; 1619 1620 PetscFunctionBegin; 1621 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1622 if (sf) PetscValidPointer(sf, 3); 1623 PetscValidPointer(dmOverlap, 4); 1624 1625 ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1626 ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 1627 1628 /* Compute point overlap with neighbouring processes on the distributed DM */ 1629 ierr = PetscLogEventBegin(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr); 1630 ierr = PetscSectionCreate(comm, &rootSection);CHKERRQ(ierr); 1631 ierr = PetscSectionCreate(comm, &leafSection);CHKERRQ(ierr); 1632 ierr = DMPlexDistributeOwnership(dm, rootSection, &rootrank, leafSection, &leafrank);CHKERRQ(ierr); 1633 ierr = DMPlexCreateOverlap(dm, overlap, rootSection, rootrank, leafSection, leafrank, &lblOverlap);CHKERRQ(ierr); 1634 /* Convert overlap label to stratified migration SF */ 1635 ierr = DMPlexPartitionLabelCreateSF(dm, lblOverlap, &sfOverlap);CHKERRQ(ierr); 1636 ierr = DMPlexStratifyMigrationSF(dm, sfOverlap, &sfStratified);CHKERRQ(ierr); 1637 ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr); 1638 sfOverlap = sfStratified; 1639 ierr = PetscObjectSetName((PetscObject) sfOverlap, "Overlap SF");CHKERRQ(ierr); 1640 ierr = PetscSFSetFromOptions(sfOverlap);CHKERRQ(ierr); 1641 1642 ierr = PetscSectionDestroy(&rootSection);CHKERRQ(ierr); 1643 ierr = PetscSectionDestroy(&leafSection);CHKERRQ(ierr); 1644 ierr = ISDestroy(&rootrank);CHKERRQ(ierr); 1645 ierr = ISDestroy(&leafrank);CHKERRQ(ierr); 1646 ierr = PetscLogEventEnd(PETSCPARTITIONER_Partition,dm,0,0,0);CHKERRQ(ierr); 1647 1648 /* Build the overlapping DM */ 1649 ierr = DMPlexCreate(comm, dmOverlap);CHKERRQ(ierr); 1650 ierr = PetscObjectSetName((PetscObject) *dmOverlap, "Parallel Mesh");CHKERRQ(ierr); 1651 ierr = DMPlexMigrate(dm, sfOverlap, *dmOverlap);CHKERRQ(ierr); 1652 /* Build the new point SF */ 1653 ierr = DMPlexCreatePointSF(dm, sfOverlap, PETSC_FALSE, &sfPoint);CHKERRQ(ierr); 1654 ierr = DMSetPointSF(*dmOverlap, sfPoint);CHKERRQ(ierr); 1655 ierr = PetscSFDestroy(&sfPoint);CHKERRQ(ierr); 1656 /* Cleanup overlap partition */ 1657 ierr = DMLabelDestroy(&lblOverlap);CHKERRQ(ierr); 1658 if (sf) *sf = sfOverlap; 1659 else {ierr = PetscSFDestroy(&sfOverlap);CHKERRQ(ierr);} 1660 PetscFunctionReturn(0); 1661 } 1662