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