1 #include <petsc/private/dmforestimpl.h> /*I "petscdmforest.h" I*/ 2 #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 3 #include <petsc/private/dmlabelimpl.h> /*I "petscdmlabel.h" I*/ 4 #include <petscsf.h> 5 6 PetscBool DMForestPackageInitialized = PETSC_FALSE; 7 8 typedef struct _DMForestTypeLink *DMForestTypeLink; 9 10 struct _DMForestTypeLink { 11 char *name; 12 DMForestTypeLink next; 13 }; 14 15 DMForestTypeLink DMForestTypeList; 16 17 static PetscErrorCode DMForestPackageFinalize(void) 18 { 19 DMForestTypeLink oldLink, link = DMForestTypeList; 20 21 PetscFunctionBegin; 22 while (link) { 23 oldLink = link; 24 PetscCall(PetscFree(oldLink->name)); 25 link = oldLink->next; 26 PetscCall(PetscFree(oldLink)); 27 } 28 PetscFunctionReturn(PETSC_SUCCESS); 29 } 30 31 static PetscErrorCode DMForestPackageInitialize(void) 32 { 33 PetscFunctionBegin; 34 if (DMForestPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 35 DMForestPackageInitialized = PETSC_TRUE; 36 37 PetscCall(DMForestRegisterType(DMFOREST)); 38 PetscCall(PetscRegisterFinalize(DMForestPackageFinalize)); 39 PetscFunctionReturn(PETSC_SUCCESS); 40 } 41 42 /*@C 43 DMForestRegisterType - Registers a `DMType` as a subtype of `DMFOREST` (so that `DMIsForest()` will be correct) 44 45 Not Collective 46 47 Input Parameter: 48 . name - the name of the type 49 50 Level: advanced 51 52 .seealso: `DMFOREST`, `DMIsForest()` 53 @*/ 54 PetscErrorCode DMForestRegisterType(DMType name) 55 { 56 DMForestTypeLink link; 57 58 PetscFunctionBegin; 59 PetscCall(DMForestPackageInitialize()); 60 PetscCall(PetscNew(&link)); 61 PetscCall(PetscStrallocpy(name, &link->name)); 62 link->next = DMForestTypeList; 63 DMForestTypeList = link; 64 PetscFunctionReturn(PETSC_SUCCESS); 65 } 66 67 /*@ 68 DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes 69 70 Not Collective 71 72 Input Parameter: 73 . dm - the DM object 74 75 Output Parameter: 76 . isForest - whether dm is a subtype of DMFOREST 77 78 Level: intermediate 79 80 .seealso: `DMFOREST`, `DMForestRegisterType()` 81 @*/ 82 PetscErrorCode DMIsForest(DM dm, PetscBool *isForest) 83 { 84 DMForestTypeLink link = DMForestTypeList; 85 86 PetscFunctionBegin; 87 while (link) { 88 PetscBool sameType; 89 PetscCall(PetscObjectTypeCompare((PetscObject)dm, link->name, &sameType)); 90 if (sameType) { 91 *isForest = PETSC_TRUE; 92 PetscFunctionReturn(PETSC_SUCCESS); 93 } 94 link = link->next; 95 } 96 *isForest = PETSC_FALSE; 97 PetscFunctionReturn(PETSC_SUCCESS); 98 } 99 100 /*@ 101 DMForestTemplate - Create a new `DM` that will be adapted from a source `DM`. The new `DM` reproduces the configuration 102 of the source, but is not yet setup, so that the user can then define only the ways that the new `DM` should differ 103 (by, e.g., refinement or repartitioning). The source `DM` is also set as the adaptivity source `DM` of the new `DM` (see 104 `DMForestSetAdaptivityForest()`). 105 106 Collective 107 108 Input Parameters: 109 + dm - the source `DM` object 110 - comm - the communicator for the new `DM` (this communicator is currently ignored, but is present so that `DMForestTemplate()` can be used within `DMCoarsen()`) 111 112 Output Parameter: 113 . tdm - the new `DM` object 114 115 Level: intermediate 116 117 .seealso: `DM`, `DMFOREST`, `DMForestSetAdaptivityForest()` 118 @*/ 119 PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm) 120 { 121 DM_Forest *forest = (DM_Forest *)dm->data; 122 DMType type; 123 DM base; 124 DMForestTopology topology; 125 MatType mtype; 126 PetscInt dim, overlap, ref, factor; 127 DMForestAdaptivityStrategy strat; 128 void *ctx; 129 PetscErrorCode (*map)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void *); 130 void *mapCtx; 131 132 PetscFunctionBegin; 133 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 134 PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), tdm)); 135 PetscCall(DMGetType(dm, &type)); 136 PetscCall(DMSetType(*tdm, type)); 137 PetscCall(DMForestGetBaseDM(dm, &base)); 138 PetscCall(DMForestSetBaseDM(*tdm, base)); 139 PetscCall(DMForestGetTopology(dm, &topology)); 140 PetscCall(DMForestSetTopology(*tdm, topology)); 141 PetscCall(DMForestGetAdjacencyDimension(dm, &dim)); 142 PetscCall(DMForestSetAdjacencyDimension(*tdm, dim)); 143 PetscCall(DMForestGetPartitionOverlap(dm, &overlap)); 144 PetscCall(DMForestSetPartitionOverlap(*tdm, overlap)); 145 PetscCall(DMForestGetMinimumRefinement(dm, &ref)); 146 PetscCall(DMForestSetMinimumRefinement(*tdm, ref)); 147 PetscCall(DMForestGetMaximumRefinement(dm, &ref)); 148 PetscCall(DMForestSetMaximumRefinement(*tdm, ref)); 149 PetscCall(DMForestGetAdaptivityStrategy(dm, &strat)); 150 PetscCall(DMForestSetAdaptivityStrategy(*tdm, strat)); 151 PetscCall(DMForestGetGradeFactor(dm, &factor)); 152 PetscCall(DMForestSetGradeFactor(*tdm, factor)); 153 PetscCall(DMForestGetBaseCoordinateMapping(dm, &map, &mapCtx)); 154 PetscCall(DMForestSetBaseCoordinateMapping(*tdm, map, mapCtx)); 155 if (forest->ftemplate) PetscCall((*forest->ftemplate)(dm, *tdm)); 156 PetscCall(DMForestSetAdaptivityForest(*tdm, dm)); 157 PetscCall(DMCopyDisc(dm, *tdm)); 158 PetscCall(DMGetApplicationContext(dm, &ctx)); 159 PetscCall(DMSetApplicationContext(*tdm, &ctx)); 160 { 161 const PetscReal *maxCell, *L, *Lstart; 162 163 PetscCall(DMGetPeriodicity(dm, &maxCell, &Lstart, &L)); 164 PetscCall(DMSetPeriodicity(*tdm, maxCell, Lstart, L)); 165 } 166 PetscCall(DMGetMatType(dm, &mtype)); 167 PetscCall(DMSetMatType(*tdm, mtype)); 168 PetscFunctionReturn(PETSC_SUCCESS); 169 } 170 171 static PetscErrorCode DMInitialize_Forest(DM dm); 172 173 PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm) 174 { 175 DM_Forest *forest = (DM_Forest *)dm->data; 176 const char *type; 177 178 PetscFunctionBegin; 179 forest->refct++; 180 (*newdm)->data = forest; 181 PetscCall(PetscObjectGetType((PetscObject)dm, &type)); 182 PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, type)); 183 PetscCall(DMInitialize_Forest(*newdm)); 184 PetscFunctionReturn(PETSC_SUCCESS); 185 } 186 187 static PetscErrorCode DMDestroy_Forest(DM dm) 188 { 189 DM_Forest *forest = (DM_Forest *)dm->data; 190 191 PetscFunctionBegin; 192 if (--forest->refct > 0) PetscFunctionReturn(PETSC_SUCCESS); 193 if (forest->destroy) PetscCall((*forest->destroy)(dm)); 194 PetscCall(PetscSFDestroy(&forest->cellSF)); 195 PetscCall(PetscSFDestroy(&forest->preCoarseToFine)); 196 PetscCall(PetscSFDestroy(&forest->coarseToPreFine)); 197 PetscCall(DMLabelDestroy(&forest->adaptLabel)); 198 PetscCall(PetscFree(forest->adaptStrategy)); 199 PetscCall(DMDestroy(&forest->base)); 200 PetscCall(DMDestroy(&forest->adapt)); 201 PetscCall(PetscFree(forest->topology)); 202 PetscCall(PetscFree(forest)); 203 PetscFunctionReturn(PETSC_SUCCESS); 204 } 205 206 /*@C 207 DMForestSetTopology - Set the topology of a `DMFOREST` during the pre-setup phase. The topology is a string (e.g. 208 "cube", "shell") and can be interpreted by subtypes of `DMFOREST`) to construct the base DM of a forest during 209 `DMSetUp()`. 210 211 Logically collectiv 212 213 Input Parameters: 214 + dm - the forest 215 - topology - the topology of the forest 216 217 Level: intermediate 218 219 .seealso: `DM`, `DMFOREST`, `DMForestGetTopology()`, `DMForestSetBaseDM()` 220 @*/ 221 PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology) 222 { 223 DM_Forest *forest = (DM_Forest *)dm->data; 224 225 PetscFunctionBegin; 226 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 227 PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the topology after setup"); 228 PetscCall(PetscFree(forest->topology)); 229 PetscCall(PetscStrallocpy((const char *)topology, (char **)&forest->topology)); 230 PetscFunctionReturn(PETSC_SUCCESS); 231 } 232 233 /*@C 234 DMForestGetTopology - Get a string describing the topology of a `DMFOREST`. 235 236 Not Collective 237 238 Input Parameter: 239 . dm - the forest 240 241 Output Parameter: 242 . topology - the topology of the forest (e.g., 'cube', 'shell') 243 244 Level: intermediate 245 246 .seealso: `DM`, `DMFOREST`, `DMForestSetTopology()` 247 @*/ 248 PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology) 249 { 250 DM_Forest *forest = (DM_Forest *)dm->data; 251 252 PetscFunctionBegin; 253 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 254 PetscValidPointer(topology, 2); 255 *topology = forest->topology; 256 PetscFunctionReturn(PETSC_SUCCESS); 257 } 258 259 /*@ 260 DMForestSetBaseDM - During the pre-setup phase, set the `DM` that defines the base mesh of a `DMFOREST` forest. The 261 forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its 262 base. In general, two forest must share a base to be comparable, to do things like construct interpolators. 263 264 Logically Collective 265 266 Input Parameters: 267 + dm - the forest 268 - base - the base `DM` of the forest 269 270 Level: intermediate 271 272 Note: 273 Currently the base `DM` must be a `DMPLEX` 274 275 .seealso: `DM`, `DMFOREST`, `DMForestGetBaseDM()` 276 @*/ 277 PetscErrorCode DMForestSetBaseDM(DM dm, DM base) 278 { 279 DM_Forest *forest = (DM_Forest *)dm->data; 280 PetscInt dim, dimEmbed; 281 282 PetscFunctionBegin; 283 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 284 PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the base after setup"); 285 PetscCall(PetscObjectReference((PetscObject)base)); 286 PetscCall(DMDestroy(&forest->base)); 287 forest->base = base; 288 if (base) { 289 const PetscReal *maxCell, *Lstart, *L; 290 291 PetscValidHeaderSpecific(base, DM_CLASSID, 2); 292 PetscCall(DMGetDimension(base, &dim)); 293 PetscCall(DMSetDimension(dm, dim)); 294 PetscCall(DMGetCoordinateDim(base, &dimEmbed)); 295 PetscCall(DMSetCoordinateDim(dm, dimEmbed)); 296 PetscCall(DMGetPeriodicity(base, &maxCell, &Lstart, &L)); 297 PetscCall(DMSetPeriodicity(dm, maxCell, Lstart, L)); 298 } else PetscCall(DMSetPeriodicity(dm, NULL, NULL, NULL)); 299 PetscFunctionReturn(PETSC_SUCCESS); 300 } 301 302 /*@ 303 DMForestGetBaseDM - Get the base DM of a DMForest forest. The forest will be hierarchically refined from the base, 304 and all refinements/coarsenings of the forest will share its base. In general, two forest must share a base to be 305 comparable, to do things like construct interpolators. 306 307 Not Collective 308 309 Input Parameter: 310 . dm - the forest 311 312 Output Parameter: 313 . base - the base DM of the forest 314 315 Notes: 316 After DMSetUp(), the base DM will be redundantly distributed across MPI processes 317 318 Level: intermediate 319 320 .seealso: `DM`, `DMFOREST`, `DMForestSetBaseDM()` 321 @*/ 322 PetscErrorCode DMForestGetBaseDM(DM dm, DM *base) 323 { 324 DM_Forest *forest = (DM_Forest *)dm->data; 325 326 PetscFunctionBegin; 327 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 328 PetscValidPointer(base, 2); 329 *base = forest->base; 330 PetscFunctionReturn(PETSC_SUCCESS); 331 } 332 333 PetscErrorCode DMForestSetBaseCoordinateMapping(DM dm, PetscErrorCode (*func)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void *), void *ctx) 334 { 335 DM_Forest *forest = (DM_Forest *)dm->data; 336 337 PetscFunctionBegin; 338 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 339 forest->mapcoordinates = func; 340 forest->mapcoordinatesctx = ctx; 341 PetscFunctionReturn(PETSC_SUCCESS); 342 } 343 344 PetscErrorCode DMForestGetBaseCoordinateMapping(DM dm, PetscErrorCode (**func)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void *), void *ctx) 345 { 346 DM_Forest *forest = (DM_Forest *)dm->data; 347 348 PetscFunctionBegin; 349 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 350 if (func) *func = forest->mapcoordinates; 351 if (ctx) *((void **)ctx) = forest->mapcoordinatesctx; 352 PetscFunctionReturn(PETSC_SUCCESS); 353 } 354 355 /*@ 356 DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be 357 adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) in `DMSetUp()`. Usually not needed 358 by users directly: `DMForestTemplate()` constructs a new forest to be adapted from an old forest and calls this 359 routine. 360 361 Logically Collective 362 363 Input Parameters: 364 + dm - the new forest, which will be constructed from adapt 365 - adapt - the old forest 366 367 Level: intermediate 368 369 Note: 370 This can be called after setup with `adapt` = `NULL`, which will clear all internal data related to the 371 adaptivity forest from `dm`. This way, repeatedly adapting does not leave stale `DM` objects in memory. 372 373 .seealso: `DM`, `DMFOREST`, `DMForestGetAdaptivityForest()`, `DMForestSetAdaptivityPurpose()` 374 @*/ 375 PetscErrorCode DMForestSetAdaptivityForest(DM dm, DM adapt) 376 { 377 DM_Forest *forest, *adaptForest, *oldAdaptForest; 378 DM oldAdapt; 379 PetscBool isForest; 380 381 PetscFunctionBegin; 382 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 383 if (adapt) PetscValidHeaderSpecific(adapt, DM_CLASSID, 2); 384 PetscCall(DMIsForest(dm, &isForest)); 385 if (!isForest) PetscFunctionReturn(PETSC_SUCCESS); 386 PetscCheck(adapt == NULL || !dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the adaptation forest after setup"); 387 forest = (DM_Forest *)dm->data; 388 PetscCall(DMForestGetAdaptivityForest(dm, &oldAdapt)); 389 adaptForest = (DM_Forest *)(adapt ? adapt->data : NULL); 390 oldAdaptForest = (DM_Forest *)(oldAdapt ? oldAdapt->data : NULL); 391 if (adaptForest != oldAdaptForest) { 392 PetscCall(PetscSFDestroy(&forest->preCoarseToFine)); 393 PetscCall(PetscSFDestroy(&forest->coarseToPreFine)); 394 if (forest->clearadaptivityforest) PetscCall((*forest->clearadaptivityforest)(dm)); 395 } 396 switch (forest->adaptPurpose) { 397 case DM_ADAPT_DETERMINE: 398 PetscCall(PetscObjectReference((PetscObject)adapt)); 399 PetscCall(DMDestroy(&(forest->adapt))); 400 forest->adapt = adapt; 401 break; 402 case DM_ADAPT_REFINE: 403 PetscCall(DMSetCoarseDM(dm, adapt)); 404 break; 405 case DM_ADAPT_COARSEN: 406 case DM_ADAPT_COARSEN_LAST: 407 PetscCall(DMSetFineDM(dm, adapt)); 408 break; 409 default: 410 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid adaptivity purpose"); 411 } 412 PetscFunctionReturn(PETSC_SUCCESS); 413 } 414 415 /*@ 416 DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted. 417 418 Not Collective 419 420 Input Parameter: 421 . dm - the forest 422 423 Output Parameter: 424 . adapt - the forest from which `dm` is/was adapted 425 426 Level: intermediate 427 428 .seealso: `DM`, `DMFOREST`, `DMForestSetAdaptivityForest()`, `DMForestSetAdaptivityPurpose()` 429 @*/ 430 PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt) 431 { 432 DM_Forest *forest; 433 434 PetscFunctionBegin; 435 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 436 forest = (DM_Forest *)dm->data; 437 switch (forest->adaptPurpose) { 438 case DM_ADAPT_DETERMINE: 439 *adapt = forest->adapt; 440 break; 441 case DM_ADAPT_REFINE: 442 PetscCall(DMGetCoarseDM(dm, adapt)); 443 break; 444 case DM_ADAPT_COARSEN: 445 case DM_ADAPT_COARSEN_LAST: 446 PetscCall(DMGetFineDM(dm, adapt)); 447 break; 448 default: 449 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid adaptivity purpose"); 450 } 451 PetscFunctionReturn(PETSC_SUCCESS); 452 } 453 454 /*@ 455 DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current `DM` is being adapted from its 456 source (set with `DMForestSetAdaptivityForest()`) for the purpose of refinement (`DM_ADAPT_REFINE`), coarsening 457 (`DM_ADAPT_COARSEN`), or undefined (`DM_ADAPT_DETERMINE`). This only matters for the purposes of reference counting: 458 during `DMDestroy()`, cyclic references can be found between `DM`s only if the cyclic reference is due to a fine/coarse 459 relationship (see `DMSetFineDM()`/`DMSetCoarseDM()`). If the purpose is not refinement or coarsening, and the user does 460 not maintain a reference to the post-adaptation forest (i.e., the one created by `DMForestTemplate()`), then this can 461 cause a memory leak. This method is used by subtypes of `DMFOREST` when automatically constructing mesh hierarchies. 462 463 Logically Collective 464 465 Input Parameters: 466 + dm - the forest 467 - purpose - the adaptivity purpose 468 469 Level: advanced 470 471 .seealso: `DM`, `DMFOREST`, `DMForestTemplate()`, `DMForestSetAdaptivityForest()`, `DMForestGetAdaptivityForest()`, `DMAdaptFlag` 472 @*/ 473 PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMAdaptFlag purpose) 474 { 475 DM_Forest *forest; 476 477 PetscFunctionBegin; 478 forest = (DM_Forest *)dm->data; 479 PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the adaptation forest after setup"); 480 if (purpose != forest->adaptPurpose) { 481 DM adapt; 482 483 PetscCall(DMForestGetAdaptivityForest(dm, &adapt)); 484 PetscCall(PetscObjectReference((PetscObject)adapt)); 485 PetscCall(DMForestSetAdaptivityForest(dm, NULL)); 486 487 forest->adaptPurpose = purpose; 488 489 PetscCall(DMForestSetAdaptivityForest(dm, adapt)); 490 PetscCall(DMDestroy(&adapt)); 491 } 492 PetscFunctionReturn(PETSC_SUCCESS); 493 } 494 495 /*@ 496 DMForestGetAdaptivityPurpose - Get whether the current `DM` is being adapted from its source (set with 497 `DMForestSetAdaptivityForest()`) for the purpose of refinement (`DM_ADAPT_REFINE`), coarsening (`DM_ADAPT_COARSEN`), 498 coarsening only the last level (`DM_ADAPT_COARSEN_LAST`) or undefined (`DM_ADAPT_DETERMINE`). 499 This only matters for the purposes of reference counting: during `DMDestroy()`, cyclic 500 501 References Can Be Found Between `Dm`S Only If The Cyclic Reference Is Due To A Fine/Coarse Relationship (See 502 `DMSetFineDM()`/`DMSetCoarseDM()`). If the purpose is not refinement or coarsening, and the user does not maintain a 503 reference to the post-adaptation forest (i.e., the one created by `DMForestTemplate()`), then this can cause a memory 504 leak. This method is used by subtypes of `DMFOREST` when automatically constructing mesh hierarchies. 505 506 Not Collective 507 508 Input Parameter: 509 . dm - the forest 510 511 Output Parameter: 512 . purpose - the adaptivity purpose 513 514 Level: advanced 515 516 .seealso: `DM`, `DMFOREST`, `DMForestTemplate()`, `DMForestSetAdaptivityForest()`, `DMForestGetAdaptivityForest()`, `DMAdaptFlag` 517 @*/ 518 PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMAdaptFlag *purpose) 519 { 520 DM_Forest *forest; 521 522 PetscFunctionBegin; 523 forest = (DM_Forest *)dm->data; 524 *purpose = forest->adaptPurpose; 525 PetscFunctionReturn(PETSC_SUCCESS); 526 } 527 528 /*@ 529 DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine 530 cell adjacency (for the purposes of partitioning and overlap). 531 532 Logically Collective 533 534 Input Parameters: 535 + dm - the forest 536 - adjDim - default 0 (i.e., vertices determine adjacency) 537 538 Level: intermediate 539 540 .seealso: `DM`, `DMFOREST`, `DMForestGetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()`, `DMForestSetPartitionOverlap()` 541 @*/ 542 PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim) 543 { 544 PetscInt dim; 545 DM_Forest *forest = (DM_Forest *)dm->data; 546 547 PetscFunctionBegin; 548 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 549 PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the adjacency dimension after setup"); 550 PetscCheck(adjDim >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "adjacency dim cannot be < 0: %" PetscInt_FMT, adjDim); 551 PetscCall(DMGetDimension(dm, &dim)); 552 PetscCheck(adjDim <= dim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "adjacency dim cannot be > %" PetscInt_FMT ": %" PetscInt_FMT, dim, adjDim); 553 forest->adjDim = adjDim; 554 PetscFunctionReturn(PETSC_SUCCESS); 555 } 556 557 /*@ 558 DMForestSetAdjacencyCodimension - Like `DMForestSetAdjacencyDimension()`, but specified as a co-dimension (so that, 559 e.g., adjacency based on facets can be specified by codimension 1 in all cases) 560 561 Logically Collective 562 563 Input Parameters: 564 + dm - the forest 565 - adjCodim - default is the dimension of the forest (see `DMGetDimension()`), since this is the codimension of vertices 566 567 Level: intermediate 568 569 .seealso: `DM`, `DMFOREST`, `DMForestGetAdjacencyCodimension()`, `DMForestSetAdjacencyDimension()` 570 @*/ 571 PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim) 572 { 573 PetscInt dim; 574 575 PetscFunctionBegin; 576 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 577 PetscCall(DMGetDimension(dm, &dim)); 578 PetscCall(DMForestSetAdjacencyDimension(dm, dim - adjCodim)); 579 PetscFunctionReturn(PETSC_SUCCESS); 580 } 581 582 /*@ 583 DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the 584 purposes of partitioning and overlap). 585 586 Not Collective 587 588 Input Parameter: 589 . dm - the forest 590 591 Output Parameter: 592 . adjDim - default 0 (i.e., vertices determine adjacency) 593 594 Level: intermediate 595 596 .seealso: `DM`, `DMFOREST`, `DMForestSetAdjacencyDimension()`, `DMForestGetAdjacencyCodimension()`, `DMForestSetPartitionOverlap()` 597 @*/ 598 PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim) 599 { 600 DM_Forest *forest = (DM_Forest *)dm->data; 601 602 PetscFunctionBegin; 603 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 604 PetscValidPointer(adjDim, 2); 605 *adjDim = forest->adjDim; 606 PetscFunctionReturn(PETSC_SUCCESS); 607 } 608 609 /*@ 610 DMForestGetAdjacencyCodimension - Like `DMForestGetAdjacencyDimension()`, but specified as a co-dimension (so that, 611 e.g., adjacency based on facets can be specified by codimension 1 in all cases) 612 613 Not Collective 614 615 Input Parameter: 616 . dm - the forest 617 618 Output Parameter: 619 . adjCodim - default isthe dimension of the forest (see `DMGetDimension()`), since this is the codimension of vertices 620 621 Level: intermediate 622 623 .seealso: `DM`, `DMFOREST`, `DMForestSetAdjacencyCodimension()`, `DMForestGetAdjacencyDimension()` 624 @*/ 625 PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim) 626 { 627 DM_Forest *forest = (DM_Forest *)dm->data; 628 PetscInt dim; 629 630 PetscFunctionBegin; 631 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 632 PetscValidPointer(adjCodim, 2); 633 PetscCall(DMGetDimension(dm, &dim)); 634 *adjCodim = dim - forest->adjDim; 635 PetscFunctionReturn(PETSC_SUCCESS); 636 } 637 638 /*@ 639 DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel 640 partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding 641 adjacent cells 642 643 Logically Collective 644 645 Input Parameters: 646 + dm - the forest 647 - overlap - default 0 648 649 Level: intermediate 650 651 .seealso: `DM`, `DMFOREST`, `DMForestGetPartitionOverlap()`, `DMForestSetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()` 652 @*/ 653 PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap) 654 { 655 DM_Forest *forest = (DM_Forest *)dm->data; 656 657 PetscFunctionBegin; 658 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 659 PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the overlap after setup"); 660 PetscCheck(overlap >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "overlap cannot be < 0: %" PetscInt_FMT, overlap); 661 forest->overlap = overlap; 662 PetscFunctionReturn(PETSC_SUCCESS); 663 } 664 665 /*@ 666 DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values 667 > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells 668 669 Not Collective 670 671 Input Parameter: 672 . dm - the forest 673 674 Output Parameter: 675 . overlap - default 0 676 677 Level: intermediate 678 679 .seealso: `DM`, `DMFOREST`, `DMForestGetPartitionOverlap()`, `DMForestSetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()` 680 @*/ 681 PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap) 682 { 683 DM_Forest *forest = (DM_Forest *)dm->data; 684 685 PetscFunctionBegin; 686 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 687 PetscValidPointer(overlap, 2); 688 *overlap = forest->overlap; 689 PetscFunctionReturn(PETSC_SUCCESS); 690 } 691 692 /*@ 693 DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base 694 `DM`, see `DMForestGetBaseDM()`) allowed in the forest. If the forest is being created by coarsening a previous forest 695 (see `DMForestGetAdaptivityForest()`) this limits the amount of coarsening. 696 697 Logically Collective 698 699 Input Parameters: 700 + dm - the forest 701 - minRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`) 702 703 Level: intermediate 704 705 .seealso: `DM`, `DMFOREST`, `DMForestGetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestSetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()` 706 @*/ 707 PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement) 708 { 709 DM_Forest *forest = (DM_Forest *)dm->data; 710 711 PetscFunctionBegin; 712 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 713 PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the minimum refinement after setup"); 714 forest->minRefinement = minRefinement; 715 PetscFunctionReturn(PETSC_SUCCESS); 716 } 717 718 /*@ 719 DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base `DM`, see 720 `DMForestGetBaseDM()`) allowed in the forest. If the forest is being created by coarsening a previous forest (see 721 `DMForestGetAdaptivityForest()`), this limits the amount of coarsening. 722 723 Not Collective 724 725 Input Parameter: 726 . dm - the forest 727 728 Output Parameter: 729 . minRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`) 730 731 Level: intermediate 732 733 .seealso: `DM`, `DMFOREST`, `DMForestSetMinimumRefinement()`, `DMForestGetMaximumRefinement()`, `DMForestGetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()` 734 @*/ 735 PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement) 736 { 737 DM_Forest *forest = (DM_Forest *)dm->data; 738 739 PetscFunctionBegin; 740 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 741 PetscValidPointer(minRefinement, 2); 742 *minRefinement = forest->minRefinement; 743 PetscFunctionReturn(PETSC_SUCCESS); 744 } 745 746 /*@ 747 DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base 748 `DM`, see `DMForestGetBaseDM()`) allowed in the forest. 749 750 Logically Collective 751 752 Input Parameters: 753 + dm - the forest 754 - initRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`) 755 756 Level: intermediate 757 758 .seealso: `DM`, `DMFOREST`, `DMForestSetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestGetBaseDM()` 759 @*/ 760 PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement) 761 { 762 DM_Forest *forest = (DM_Forest *)dm->data; 763 764 PetscFunctionBegin; 765 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 766 PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the initial refinement after setup"); 767 forest->initRefinement = initRefinement; 768 PetscFunctionReturn(PETSC_SUCCESS); 769 } 770 771 /*@ 772 DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base `DM`, see 773 `DMForestGetBaseDM()`) allowed in the forest. 774 775 Not Collective 776 777 Input Parameter: 778 . dm - the forest 779 780 Output Parameter: 781 . initRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`) 782 783 Level: intermediate 784 785 .seealso: `DM`, `DMFOREST`, `DMForestSetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestGetBaseDM()` 786 @*/ 787 PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement) 788 { 789 DM_Forest *forest = (DM_Forest *)dm->data; 790 791 PetscFunctionBegin; 792 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 793 PetscValidPointer(initRefinement, 2); 794 *initRefinement = forest->initRefinement; 795 PetscFunctionReturn(PETSC_SUCCESS); 796 } 797 798 /*@ 799 DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base 800 `DM`, see `DMForestGetBaseDM()`) allowed in the forest. If the forest is being created by refining a previous forest 801 (see `DMForestGetAdaptivityForest()`), this limits the amount of refinement. 802 803 Logically Collective 804 805 Input Parameters: 806 + dm - the forest 807 - maxRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`) 808 809 Level: intermediate 810 811 .seealso: `DM`, `DMFOREST`, `DMForestGetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestSetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityDM()` 812 @*/ 813 PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement) 814 { 815 DM_Forest *forest = (DM_Forest *)dm->data; 816 817 PetscFunctionBegin; 818 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 819 PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the maximum refinement after setup"); 820 forest->maxRefinement = maxRefinement; 821 PetscFunctionReturn(PETSC_SUCCESS); 822 } 823 824 /*@ 825 DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base `DM`, see 826 `DMForestGetBaseDM()`) allowed in the forest. If the forest is being created by refining a previous forest (see 827 `DMForestGetAdaptivityForest()`), this limits the amount of refinement. 828 829 Not Collective 830 831 Input Parameter: 832 . dm - the forest 833 834 Output Parameter: 835 . maxRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`) 836 837 Level: intermediate 838 839 .seealso: `DM`, `DMFOREST`, `DMForestSetMaximumRefinement()`, `DMForestGetMinimumRefinement()`, `DMForestGetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()` 840 @*/ 841 PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement) 842 { 843 DM_Forest *forest = (DM_Forest *)dm->data; 844 845 PetscFunctionBegin; 846 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 847 PetscValidPointer(maxRefinement, 2); 848 *maxRefinement = forest->maxRefinement; 849 PetscFunctionReturn(PETSC_SUCCESS); 850 } 851 852 /*@C 853 DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes. 854 855 Logically Collective 856 857 Input Parameters: 858 + dm - the forest 859 - adaptStrategy - default `DMFORESTADAPTALL` 860 861 Level: advanced 862 863 Notes: 864 Subtypes of `DMFOREST` may define their own strategies. Two default strategies are `DMFORESTADAPTALL`, which indicates that all processes must agree 865 for a refinement/coarsening flag to be valid, and `DMFORESTADAPTANY`, which indicates that only one process needs to 866 specify refinement/coarsening. 867 868 .seealso: `DM`, `DMFOREST`, `DMForestGetAdaptivityStrategy()`, `DMFORESTADAPTALL`, `DMFORESTADAPTANY` 869 @*/ 870 PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy) 871 { 872 DM_Forest *forest = (DM_Forest *)dm->data; 873 874 PetscFunctionBegin; 875 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 876 PetscCall(PetscFree(forest->adaptStrategy)); 877 PetscCall(PetscStrallocpy((const char *)adaptStrategy, (char **)&forest->adaptStrategy)); 878 PetscFunctionReturn(PETSC_SUCCESS); 879 } 880 881 /*@C 882 DMForestGetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes. 883 884 Not Collective 885 886 Input Parameter: 887 . dm - the forest 888 889 Output Parameter: 890 . adaptStrategy - the adaptivity strategy (default `DMFORESTADAPTALL`) 891 892 Level: advanced 893 894 Note: 895 Subtypes 896 of `DMFOREST` may define their own strategies. Two default strategies are `DMFORESTADAPTALL`, which indicates that all 897 processes must agree for a refinement/coarsening flag to be valid, and `DMFORESTADAPTANY`, which indicates that only 898 one process needs to specify refinement/coarsening. 899 900 .seealso: `DM`, `DMFOREST`, `DMFORESTADAPTALL`, `DMFORESTADAPTANY`, `DMForestSetAdaptivityStrategy()` 901 @*/ 902 PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy) 903 { 904 DM_Forest *forest = (DM_Forest *)dm->data; 905 906 PetscFunctionBegin; 907 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 908 PetscValidPointer(adaptStrategy, 2); 909 *adaptStrategy = forest->adaptStrategy; 910 PetscFunctionReturn(PETSC_SUCCESS); 911 } 912 913 /*@ 914 DMForestGetAdaptivitySuccess - Return whether the requested adaptation (refinement, coarsening, repartitioning, 915 etc.) was successful. 916 917 Collective 918 919 Input Parameter: 920 . dm - the post-adaptation forest 921 922 Output Parameter: 923 . success - `PETSC_TRUE` if the post-adaptation forest is different from the pre-adaptation forest. 924 925 Level: intermediate 926 927 Notes: 928 `PETSC_FALSE` indicates that the post-adaptation forest is the same as the pre-adpatation 929 forest. A requested adaptation may have been unsuccessful if, for example, the requested refinement would have 930 exceeded the maximum refinement level. 931 932 .seealso: `DM`, `DMFOREST` 933 @*/ 934 PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success) 935 { 936 DM_Forest *forest; 937 938 PetscFunctionBegin; 939 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 940 PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() has not been called yet."); 941 forest = (DM_Forest *)dm->data; 942 PetscCall((forest->getadaptivitysuccess)(dm, success)); 943 PetscFunctionReturn(PETSC_SUCCESS); 944 } 945 946 /*@ 947 DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer `PetscSF`s should be computed 948 relating the cells of the pre-adaptation forest to the post-adaptiation forest. 949 950 Logically Collective 951 952 Input Parameters: 953 + dm - the post-adaptation forest 954 - computeSF - default `PETSC_TRUE` 955 956 Level: advanced 957 958 Note: 959 After `DMSetUp()` is called, the transfer `PetscSF`s can be accessed with `DMForestGetAdaptivitySF()`. 960 961 .seealso: `DM`, `DMFOREST`, `DMForestGetComputeAdaptivitySF()`, `DMForestGetAdaptivitySF()` 962 @*/ 963 PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF) 964 { 965 DM_Forest *forest; 966 967 PetscFunctionBegin; 968 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 969 PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot compute adaptivity PetscSFs after setup is called"); 970 forest = (DM_Forest *)dm->data; 971 forest->computeAdaptSF = computeSF; 972 PetscFunctionReturn(PETSC_SUCCESS); 973 } 974 975 PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time) 976 { 977 DM_Forest *forest; 978 979 PetscFunctionBegin; 980 PetscValidHeaderSpecific(dmIn, DM_CLASSID, 1); 981 PetscValidHeaderSpecific(vecIn, VEC_CLASSID, 2); 982 PetscValidHeaderSpecific(dmOut, DM_CLASSID, 3); 983 PetscValidHeaderSpecific(vecOut, VEC_CLASSID, 4); 984 forest = (DM_Forest *)dmIn->data; 985 PetscCheck(forest->transfervec, PetscObjectComm((PetscObject)dmIn), PETSC_ERR_SUP, "DMForestTransferVec() not implemented"); 986 PetscCall((forest->transfervec)(dmIn, vecIn, dmOut, vecOut, useBCs, time)); 987 PetscFunctionReturn(PETSC_SUCCESS); 988 } 989 990 PetscErrorCode DMForestTransferVecFromBase(DM dm, Vec vecIn, Vec vecOut) 991 { 992 DM_Forest *forest; 993 994 PetscFunctionBegin; 995 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 996 PetscValidHeaderSpecific(vecIn, VEC_CLASSID, 2); 997 PetscValidHeaderSpecific(vecOut, VEC_CLASSID, 3); 998 forest = (DM_Forest *)dm->data; 999 PetscCheck(forest->transfervecfrombase, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMForestTransferVecFromBase() not implemented"); 1000 PetscCall((forest->transfervecfrombase)(dm, vecIn, vecOut)); 1001 PetscFunctionReturn(PETSC_SUCCESS); 1002 } 1003 1004 /*@ 1005 DMForestGetComputeAdaptivitySF - Get whether transfer `PetscSF`s should be computed relating the cells of the 1006 pre-adaptation forest to the post-adaptiation forest. After `DMSetUp()` is called, these transfer PetscSFs can be 1007 accessed with `DMForestGetAdaptivitySF()`. 1008 1009 Not Collective 1010 1011 Input Parameter: 1012 . dm - the post-adaptation forest 1013 1014 Output Parameter: 1015 . computeSF - default `PETSC_TRUE` 1016 1017 Level: advanced 1018 1019 .seealso: `DM`, `DMFOREST`, `DMForestSetComputeAdaptivitySF()`, `DMForestGetAdaptivitySF()` 1020 @*/ 1021 PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF) 1022 { 1023 DM_Forest *forest; 1024 1025 PetscFunctionBegin; 1026 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1027 forest = (DM_Forest *)dm->data; 1028 *computeSF = forest->computeAdaptSF; 1029 PetscFunctionReturn(PETSC_SUCCESS); 1030 } 1031 1032 /*@ 1033 DMForestGetAdaptivitySF - Get `PetscSF`s that relate the pre-adaptation forest to the post-adaptation forest. 1034 Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be 1035 some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa. Therefore there are two 1036 `PetscSF`s: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates 1037 pre-adaptation fine cells to post-adaptation coarse cells. 1038 1039 Not Collective 1040 1041 Input Parameter: 1042 . dm - the post-adaptation forest 1043 1044 Output Parameters: 1045 + preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post- 1046 - coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre- 1047 1048 Level: advanced 1049 1050 .seealso: `DM`, `DMFOREST`, `DMForestGetComputeAdaptivitySF()`, `DMForestSetComputeAdaptivitySF()` 1051 @*/ 1052 PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine) 1053 { 1054 DM_Forest *forest; 1055 1056 PetscFunctionBegin; 1057 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1058 PetscCall(DMSetUp(dm)); 1059 forest = (DM_Forest *)dm->data; 1060 if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine; 1061 if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine; 1062 PetscFunctionReturn(PETSC_SUCCESS); 1063 } 1064 1065 /*@ 1066 DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to 1067 indicate that the diameter of neighboring cells should differ by at most a factor of 2. Subtypes of `DMFOREST` may 1068 only support one particular choice of grading factor. 1069 1070 Logically Collective 1071 1072 Input Parameters: 1073 + dm - the forest 1074 - grade - the grading factor 1075 1076 Level: advanced 1077 1078 .seealso: `DM`, `DMFOREST`, `DMForestGetGradeFactor()` 1079 @*/ 1080 PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade) 1081 { 1082 DM_Forest *forest = (DM_Forest *)dm->data; 1083 1084 PetscFunctionBegin; 1085 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1086 PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the grade factor after setup"); 1087 forest->gradeFactor = grade; 1088 PetscFunctionReturn(PETSC_SUCCESS); 1089 } 1090 1091 /*@ 1092 DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of 1093 neighboring cells should differ by at most a factor of 2. Subtypes of `DMFOREST` may only support one particular 1094 choice of grading factor. 1095 1096 Not Collective 1097 1098 Input Parameter: 1099 . dm - the forest 1100 1101 Output Parameter: 1102 . grade - the grading factor 1103 1104 Level: advanced 1105 1106 .seealso: `DM`, `DMFOREST`, `DMForestSetGradeFactor()` 1107 @*/ 1108 PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade) 1109 { 1110 DM_Forest *forest = (DM_Forest *)dm->data; 1111 1112 PetscFunctionBegin; 1113 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1114 PetscValidPointer(grade, 2); 1115 *grade = forest->gradeFactor; 1116 PetscFunctionReturn(PETSC_SUCCESS); 1117 } 1118 1119 /*@ 1120 DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes 1121 the cell weight (see `DMForestSetCellWeights()`) when calculating partitions. The final weight of a cell will be 1122 (cellWeight) * (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on 1123 its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the 1124 computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement. 1125 1126 Logically Collective 1127 1128 Input Parameters: 1129 + dm - the forest 1130 - weightsFactor - default 1. 1131 1132 Level: advanced 1133 1134 .seealso: `DM`, `DMFOREST`, `DMForestGetCellWeightFactor()`, `DMForestSetCellWeights()` 1135 @*/ 1136 PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor) 1137 { 1138 DM_Forest *forest = (DM_Forest *)dm->data; 1139 1140 PetscFunctionBegin; 1141 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1142 PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the weights factor after setup"); 1143 forest->weightsFactor = weightsFactor; 1144 PetscFunctionReturn(PETSC_SUCCESS); 1145 } 1146 1147 /*@ 1148 DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see 1149 `DMForestSetCellWeights()`) when calculating partitions. The final weight of a cell will be (cellWeight) * 1150 (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on its level; a 1151 factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation 1152 associated with a cell is multiplied by a factor of 2 for each additional level of refinement. 1153 1154 Not Collective 1155 1156 Input Parameter: 1157 . dm - the forest 1158 1159 Output Parameter: 1160 . weightsFactor - default 1. 1161 1162 Level: advanced 1163 1164 .seealso: `DM`, `DMFOREST`, `DMForestSetCellWeightFactor()`, `DMForestSetCellWeights()` 1165 @*/ 1166 PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor) 1167 { 1168 DM_Forest *forest = (DM_Forest *)dm->data; 1169 1170 PetscFunctionBegin; 1171 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1172 PetscValidPointer(weightsFactor, 2); 1173 *weightsFactor = forest->weightsFactor; 1174 PetscFunctionReturn(PETSC_SUCCESS); 1175 } 1176 1177 /*@ 1178 DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process 1179 1180 Not Collective 1181 1182 Input Parameter: 1183 . dm - the forest 1184 1185 Output Parameters: 1186 + cStart - the first cell on this process 1187 - cEnd - one after the final cell on this process 1188 1189 Level: intermediate 1190 1191 .seealso: `DM`, `DMFOREST`, `DMForestGetCellSF()` 1192 @*/ 1193 PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd) 1194 { 1195 DM_Forest *forest = (DM_Forest *)dm->data; 1196 1197 PetscFunctionBegin; 1198 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1199 PetscValidPointer(cStart, 2); 1200 PetscValidPointer(cEnd, 3); 1201 if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) PetscCall(forest->createcellchart(dm, &forest->cStart, &forest->cEnd)); 1202 *cStart = forest->cStart; 1203 *cEnd = forest->cEnd; 1204 PetscFunctionReturn(PETSC_SUCCESS); 1205 } 1206 1207 /*@ 1208 DMForestGetCellSF - After the setup phase, get the `PetscSF` for overlapping cells between processes 1209 1210 Not Collective 1211 1212 Input Parameter: 1213 . dm - the forest 1214 1215 Output Parameter: 1216 . cellSF - the `PetscSF` 1217 1218 Level: intermediate 1219 1220 .seealso: `DM`, `DMFOREST`, `DMForestGetCellChart()` 1221 @*/ 1222 PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF) 1223 { 1224 DM_Forest *forest = (DM_Forest *)dm->data; 1225 1226 PetscFunctionBegin; 1227 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1228 PetscValidPointer(cellSF, 2); 1229 if ((!forest->cellSF) && forest->createcellsf) PetscCall(forest->createcellsf(dm, &forest->cellSF)); 1230 *cellSF = forest->cellSF; 1231 PetscFunctionReturn(PETSC_SUCCESS); 1232 } 1233 1234 /*@C 1235 DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see 1236 `DMForestGetAdaptivityForest()`) that holds the adaptation flags (refinement, coarsening, or some combination). The 1237 interpretation of the label values is up to the subtype of `DMFOREST`, but `DM_ADAPT_DETERMINE`, `DM_ADAPT_KEEP`, 1238 `DM_ADAPT_REFINE`, and `DM_ADAPT_COARSEN` have been reserved as choices that should be accepted by all subtypes. 1239 1240 Logically Collective 1241 1242 Input Parameters: 1243 + dm - the forest 1244 - adaptLabel - the label in the pre-adaptation forest 1245 1246 Level: intermediate 1247 1248 .seealso: `DM`, `DMFOREST`, `DMForestGetAdaptivityLabel()` 1249 @*/ 1250 PetscErrorCode DMForestSetAdaptivityLabel(DM dm, DMLabel adaptLabel) 1251 { 1252 DM_Forest *forest = (DM_Forest *)dm->data; 1253 1254 PetscFunctionBegin; 1255 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1256 if (adaptLabel) PetscValidHeaderSpecific(adaptLabel, DMLABEL_CLASSID, 2); 1257 PetscCall(PetscObjectReference((PetscObject)adaptLabel)); 1258 PetscCall(DMLabelDestroy(&forest->adaptLabel)); 1259 forest->adaptLabel = adaptLabel; 1260 PetscFunctionReturn(PETSC_SUCCESS); 1261 } 1262 1263 /*@C 1264 DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see `DMForestGetAdaptivityForest()`) that 1265 holds the adaptation flags (refinement, coarsening, or some combination). The interpretation of the label values is 1266 up to the subtype of `DMFOREST`, but `DM_ADAPT_DETERMINE`, `DM_ADAPT_KEEP`, `DM_ADAPT_REFINE`, and `DM_ADAPT_COARSEN` have 1267 been reserved as choices that should be accepted by all subtypes. 1268 1269 Not Collective 1270 1271 Input Parameter: 1272 . dm - the forest 1273 1274 Output Parameter: 1275 . adaptLabel - the name of the label in the pre-adaptation forest 1276 1277 Level: intermediate 1278 1279 .seealso: `DM`, `DMFOREST`, `DMForestSetAdaptivityLabel()` 1280 @*/ 1281 PetscErrorCode DMForestGetAdaptivityLabel(DM dm, DMLabel *adaptLabel) 1282 { 1283 DM_Forest *forest = (DM_Forest *)dm->data; 1284 1285 PetscFunctionBegin; 1286 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1287 *adaptLabel = forest->adaptLabel; 1288 PetscFunctionReturn(PETSC_SUCCESS); 1289 } 1290 1291 /*@ 1292 DMForestSetCellWeights - Set the weights assigned to each of the cells (see `DMForestGetCellChart()`) of the current 1293 process: weights are used to determine parallel partitioning. 1294 1295 Logically Collective 1296 1297 Input Parameters: 1298 + dm - the forest 1299 . weights - the array of weights (see `DMForestSetWeightCapacity()`) for all cells, or `NULL` to indicate each cell has weight 1. 1300 - copyMode - how weights should reference weights 1301 1302 Level: advanced 1303 1304 .seealso: `DM`, `DMFOREST`, `DMForestGetCellWeights()`, `DMForestSetWeightCapacity()` 1305 @*/ 1306 PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode) 1307 { 1308 DM_Forest *forest = (DM_Forest *)dm->data; 1309 PetscInt cStart, cEnd; 1310 1311 PetscFunctionBegin; 1312 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1313 PetscCall(DMForestGetCellChart(dm, &cStart, &cEnd)); 1314 PetscCheck(cEnd >= cStart, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "cell chart [%" PetscInt_FMT ",%" PetscInt_FMT ") is not valid", cStart, cEnd); 1315 if (copyMode == PETSC_COPY_VALUES) { 1316 if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) PetscCall(PetscMalloc1(cEnd - cStart, &forest->cellWeights)); 1317 PetscCall(PetscArraycpy(forest->cellWeights, weights, cEnd - cStart)); 1318 forest->cellWeightsCopyMode = PETSC_OWN_POINTER; 1319 PetscFunctionReturn(PETSC_SUCCESS); 1320 } 1321 if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) PetscCall(PetscFree(forest->cellWeights)); 1322 forest->cellWeights = weights; 1323 forest->cellWeightsCopyMode = copyMode; 1324 PetscFunctionReturn(PETSC_SUCCESS); 1325 } 1326 1327 /*@ 1328 DMForestGetCellWeights - Get the weights assigned to each of the cells (see `DMForestGetCellChart()`) of the current 1329 process: weights are used to determine parallel partitioning. 1330 1331 Not Collective 1332 1333 Input Parameter: 1334 . dm - the forest 1335 1336 Output Parameter: 1337 . weights - the array of weights for all cells, or `NULL` to indicate each cell has weight 1. 1338 1339 Level: advanced 1340 1341 .seealso: `DM`, `DMFOREST`, `DMForestSetCellWeights()`, `DMForestSetWeightCapacity()` 1342 @*/ 1343 PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights) 1344 { 1345 DM_Forest *forest = (DM_Forest *)dm->data; 1346 1347 PetscFunctionBegin; 1348 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1349 PetscValidPointer(weights, 2); 1350 *weights = forest->cellWeights; 1351 PetscFunctionReturn(PETSC_SUCCESS); 1352 } 1353 1354 /*@ 1355 DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning 1356 a pre-adaptation forest (see `DMForestGetAdaptivityForest()`). After partitioning, the ratio of the weight of each 1357 process's cells to the process's capacity will be roughly equal for all processes. A capacity of 0 indicates that 1358 the current process should not have any cells after repartitioning. 1359 1360 Logically Collective 1361 1362 Input Parameters: 1363 + dm - the forest 1364 - capacity - this process's capacity 1365 1366 Level: advanced 1367 1368 .seealso `DM`, `DMFOREST`, `DMForestGetWeightCapacity()`, `DMForestSetCellWeights()`, `DMForestSetCellWeightFactor()` 1369 @*/ 1370 PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity) 1371 { 1372 DM_Forest *forest = (DM_Forest *)dm->data; 1373 1374 PetscFunctionBegin; 1375 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1376 PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the weight capacity after setup"); 1377 PetscCheck(capacity >= 0., PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot have negative weight capacity; %g", (double)capacity); 1378 forest->weightCapacity = capacity; 1379 PetscFunctionReturn(PETSC_SUCCESS); 1380 } 1381 1382 /*@ 1383 DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see 1384 `DMForestGetAdaptivityForest()`). After partitioning, the ratio of the weight of each process's cells to the 1385 process's capacity will be roughly equal for all processes. A capacity of 0 indicates that the current process 1386 should not have any cells after repartitioning. 1387 1388 Not Collective 1389 1390 Input Parameter: 1391 . dm - the forest 1392 1393 Output Parameter: 1394 . capacity - this process's capacity 1395 1396 Level: advanced 1397 1398 .seealso `DM`, `DMFOREST`, `DMForestSetWeightCapacity()`, `DMForestSetCellWeights()`, `DMForestSetCellWeightFactor()` 1399 @*/ 1400 PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity) 1401 { 1402 DM_Forest *forest = (DM_Forest *)dm->data; 1403 1404 PetscFunctionBegin; 1405 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1406 PetscValidPointer(capacity, 2); 1407 *capacity = forest->weightCapacity; 1408 PetscFunctionReturn(PETSC_SUCCESS); 1409 } 1410 1411 PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(DM dm, PetscOptionItems *PetscOptionsObject) 1412 { 1413 PetscBool flg, flg1, flg2, flg3, flg4; 1414 DMForestTopology oldTopo; 1415 char stringBuffer[256]; 1416 PetscViewer viewer; 1417 PetscViewerFormat format; 1418 PetscInt adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade; 1419 PetscReal weightsFactor; 1420 DMForestAdaptivityStrategy adaptStrategy; 1421 1422 PetscFunctionBegin; 1423 PetscCall(DMForestGetTopology(dm, &oldTopo)); 1424 PetscOptionsHeadBegin(PetscOptionsObject, "DMForest Options"); 1425 PetscCall(PetscOptionsString("-dm_forest_topology", "the topology of the forest's base mesh", "DMForestSetTopology", oldTopo, stringBuffer, sizeof(stringBuffer), &flg1)); 1426 PetscCall(PetscOptionsViewer("-dm_forest_base_dm", "load the base DM from a viewer specification", "DMForestSetBaseDM", &viewer, &format, &flg2)); 1427 PetscCall(PetscOptionsViewer("-dm_forest_coarse_forest", "load the coarse forest from a viewer specification", "DMForestSetCoarseForest", &viewer, &format, &flg3)); 1428 PetscCall(PetscOptionsViewer("-dm_forest_fine_forest", "load the fine forest from a viewer specification", "DMForestSetFineForest", &viewer, &format, &flg4)); 1429 PetscCheck((PetscInt)flg1 + (PetscInt)flg2 + (PetscInt)flg3 + (PetscInt)flg4 <= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_INCOMP, "Specify only one of -dm_forest_{topology,base_dm,coarse_forest,fine_forest}"); 1430 if (flg1) { 1431 PetscCall(DMForestSetTopology(dm, (DMForestTopology)stringBuffer)); 1432 PetscCall(DMForestSetBaseDM(dm, NULL)); 1433 PetscCall(DMForestSetAdaptivityForest(dm, NULL)); 1434 } 1435 if (flg2) { 1436 DM base; 1437 1438 PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &base)); 1439 PetscCall(PetscViewerPushFormat(viewer, format)); 1440 PetscCall(DMLoad(base, viewer)); 1441 PetscCall(PetscViewerDestroy(&viewer)); 1442 PetscCall(DMForestSetBaseDM(dm, base)); 1443 PetscCall(DMDestroy(&base)); 1444 PetscCall(DMForestSetTopology(dm, NULL)); 1445 PetscCall(DMForestSetAdaptivityForest(dm, NULL)); 1446 } 1447 if (flg3) { 1448 DM coarse; 1449 1450 PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &coarse)); 1451 PetscCall(PetscViewerPushFormat(viewer, format)); 1452 PetscCall(DMLoad(coarse, viewer)); 1453 PetscCall(PetscViewerDestroy(&viewer)); 1454 PetscCall(DMForestSetAdaptivityForest(dm, coarse)); 1455 PetscCall(DMDestroy(&coarse)); 1456 PetscCall(DMForestSetTopology(dm, NULL)); 1457 PetscCall(DMForestSetBaseDM(dm, NULL)); 1458 } 1459 if (flg4) { 1460 DM fine; 1461 1462 PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &fine)); 1463 PetscCall(PetscViewerPushFormat(viewer, format)); 1464 PetscCall(DMLoad(fine, viewer)); 1465 PetscCall(PetscViewerDestroy(&viewer)); 1466 PetscCall(DMForestSetAdaptivityForest(dm, fine)); 1467 PetscCall(DMDestroy(&fine)); 1468 PetscCall(DMForestSetTopology(dm, NULL)); 1469 PetscCall(DMForestSetBaseDM(dm, NULL)); 1470 } 1471 PetscCall(DMForestGetAdjacencyDimension(dm, &adjDim)); 1472 PetscCall(PetscOptionsBoundedInt("-dm_forest_adjacency_dimension", "set the dimension of points that define adjacency in the forest", "DMForestSetAdjacencyDimension", adjDim, &adjDim, &flg, 0)); 1473 if (flg) { 1474 PetscCall(DMForestSetAdjacencyDimension(dm, adjDim)); 1475 } else { 1476 PetscCall(DMForestGetAdjacencyCodimension(dm, &adjCodim)); 1477 PetscCall(PetscOptionsBoundedInt("-dm_forest_adjacency_codimension", "set the codimension of points that define adjacency in the forest", "DMForestSetAdjacencyCodimension", adjCodim, &adjCodim, &flg, 1)); 1478 if (flg) PetscCall(DMForestSetAdjacencyCodimension(dm, adjCodim)); 1479 } 1480 PetscCall(DMForestGetPartitionOverlap(dm, &overlap)); 1481 PetscCall(PetscOptionsBoundedInt("-dm_forest_partition_overlap", "set the degree of partition overlap", "DMForestSetPartitionOverlap", overlap, &overlap, &flg, 0)); 1482 if (flg) PetscCall(DMForestSetPartitionOverlap(dm, overlap)); 1483 #if 0 1484 PetscCall(PetscOptionsBoundedInt("-dm_refine","equivalent to -dm_forest_set_minimum_refinement and -dm_forest_set_initial_refinement with the same value",NULL,minRefinement,&minRefinement,&flg,0)); 1485 if (flg) { 1486 PetscCall(DMForestSetMinimumRefinement(dm,minRefinement)); 1487 PetscCall(DMForestSetInitialRefinement(dm,minRefinement)); 1488 } 1489 PetscCall(PetscOptionsBoundedInt("-dm_refine_hierarchy","equivalent to -dm_forest_set_minimum_refinement 0 and -dm_forest_set_initial_refinement",NULL,initRefinement,&initRefinement,&flg,0)); 1490 if (flg) { 1491 PetscCall(DMForestSetMinimumRefinement(dm,0)); 1492 PetscCall(DMForestSetInitialRefinement(dm,initRefinement)); 1493 } 1494 #endif 1495 PetscCall(DMForestGetMinimumRefinement(dm, &minRefinement)); 1496 PetscCall(PetscOptionsBoundedInt("-dm_forest_minimum_refinement", "set the minimum level of refinement in the forest", "DMForestSetMinimumRefinement", minRefinement, &minRefinement, &flg, 0)); 1497 if (flg) PetscCall(DMForestSetMinimumRefinement(dm, minRefinement)); 1498 PetscCall(DMForestGetInitialRefinement(dm, &initRefinement)); 1499 PetscCall(PetscOptionsBoundedInt("-dm_forest_initial_refinement", "set the initial level of refinement in the forest", "DMForestSetInitialRefinement", initRefinement, &initRefinement, &flg, 0)); 1500 if (flg) PetscCall(DMForestSetInitialRefinement(dm, initRefinement)); 1501 PetscCall(DMForestGetMaximumRefinement(dm, &maxRefinement)); 1502 PetscCall(PetscOptionsBoundedInt("-dm_forest_maximum_refinement", "set the maximum level of refinement in the forest", "DMForestSetMaximumRefinement", maxRefinement, &maxRefinement, &flg, 0)); 1503 if (flg) PetscCall(DMForestSetMaximumRefinement(dm, maxRefinement)); 1504 PetscCall(DMForestGetAdaptivityStrategy(dm, &adaptStrategy)); 1505 PetscCall(PetscOptionsString("-dm_forest_adaptivity_strategy", "the forest's adaptivity-flag resolution strategy", "DMForestSetAdaptivityStrategy", adaptStrategy, stringBuffer, sizeof(stringBuffer), &flg)); 1506 if (flg) PetscCall(DMForestSetAdaptivityStrategy(dm, (DMForestAdaptivityStrategy)stringBuffer)); 1507 PetscCall(DMForestGetGradeFactor(dm, &grade)); 1508 PetscCall(PetscOptionsBoundedInt("-dm_forest_grade_factor", "grade factor between neighboring cells", "DMForestSetGradeFactor", grade, &grade, &flg, 0)); 1509 if (flg) PetscCall(DMForestSetGradeFactor(dm, grade)); 1510 PetscCall(DMForestGetCellWeightFactor(dm, &weightsFactor)); 1511 PetscCall(PetscOptionsReal("-dm_forest_cell_weight_factor", "multiplying weight factor for cell refinement", "DMForestSetCellWeightFactor", weightsFactor, &weightsFactor, &flg)); 1512 if (flg) PetscCall(DMForestSetCellWeightFactor(dm, weightsFactor)); 1513 PetscOptionsHeadEnd(); 1514 PetscFunctionReturn(PETSC_SUCCESS); 1515 } 1516 1517 PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 1518 { 1519 PetscFunctionBegin; 1520 if (subdm) PetscCall(DMClone(dm, subdm)); 1521 PetscCall(DMCreateSectionSubDM(dm, numFields, fields, is, subdm)); 1522 PetscFunctionReturn(PETSC_SUCCESS); 1523 } 1524 1525 PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined) 1526 { 1527 DMLabel refine; 1528 DM fineDM; 1529 1530 PetscFunctionBegin; 1531 PetscCall(DMGetFineDM(dm, &fineDM)); 1532 if (fineDM) { 1533 PetscCall(PetscObjectReference((PetscObject)fineDM)); 1534 *dmRefined = fineDM; 1535 PetscFunctionReturn(PETSC_SUCCESS); 1536 } 1537 PetscCall(DMForestTemplate(dm, comm, dmRefined)); 1538 PetscCall(DMGetLabel(dm, "refine", &refine)); 1539 if (!refine) { 1540 PetscCall(DMLabelCreate(PETSC_COMM_SELF, "refine", &refine)); 1541 PetscCall(DMLabelSetDefaultValue(refine, DM_ADAPT_REFINE)); 1542 } else PetscCall(PetscObjectReference((PetscObject)refine)); 1543 PetscCall(DMForestSetAdaptivityLabel(*dmRefined, refine)); 1544 PetscCall(DMLabelDestroy(&refine)); 1545 PetscFunctionReturn(PETSC_SUCCESS); 1546 } 1547 1548 PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened) 1549 { 1550 DMLabel coarsen; 1551 DM coarseDM; 1552 1553 PetscFunctionBegin; 1554 { 1555 PetscMPIInt mpiComparison; 1556 MPI_Comm dmcomm = PetscObjectComm((PetscObject)dm); 1557 1558 PetscCallMPI(MPI_Comm_compare(comm, dmcomm, &mpiComparison)); 1559 PetscCheck(mpiComparison == MPI_IDENT || mpiComparison == MPI_CONGRUENT, dmcomm, PETSC_ERR_SUP, "No support for different communicators yet"); 1560 } 1561 PetscCall(DMGetCoarseDM(dm, &coarseDM)); 1562 if (coarseDM) { 1563 PetscCall(PetscObjectReference((PetscObject)coarseDM)); 1564 *dmCoarsened = coarseDM; 1565 PetscFunctionReturn(PETSC_SUCCESS); 1566 } 1567 PetscCall(DMForestTemplate(dm, comm, dmCoarsened)); 1568 PetscCall(DMForestSetAdaptivityPurpose(*dmCoarsened, DM_ADAPT_COARSEN)); 1569 PetscCall(DMGetLabel(dm, "coarsen", &coarsen)); 1570 if (!coarsen) { 1571 PetscCall(DMLabelCreate(PETSC_COMM_SELF, "coarsen", &coarsen)); 1572 PetscCall(DMLabelSetDefaultValue(coarsen, DM_ADAPT_COARSEN)); 1573 } else PetscCall(PetscObjectReference((PetscObject)coarsen)); 1574 PetscCall(DMForestSetAdaptivityLabel(*dmCoarsened, coarsen)); 1575 PetscCall(DMLabelDestroy(&coarsen)); 1576 PetscFunctionReturn(PETSC_SUCCESS); 1577 } 1578 1579 PetscErrorCode DMAdaptLabel_Forest(DM dm, PETSC_UNUSED Vec metric, DMLabel label, PETSC_UNUSED DMLabel rgLabel, DM *adaptedDM) 1580 { 1581 PetscBool success; 1582 1583 PetscFunctionBegin; 1584 PetscCall(DMForestTemplate(dm, PetscObjectComm((PetscObject)dm), adaptedDM)); 1585 PetscCall(DMForestSetAdaptivityLabel(*adaptedDM, label)); 1586 PetscCall(DMSetUp(*adaptedDM)); 1587 PetscCall(DMForestGetAdaptivitySuccess(*adaptedDM, &success)); 1588 if (!success) { 1589 PetscCall(DMDestroy(adaptedDM)); 1590 *adaptedDM = NULL; 1591 } 1592 PetscFunctionReturn(PETSC_SUCCESS); 1593 } 1594 1595 static PetscErrorCode DMInitialize_Forest(DM dm) 1596 { 1597 PetscFunctionBegin; 1598 PetscCall(PetscMemzero(dm->ops, sizeof(*(dm->ops)))); 1599 1600 dm->ops->clone = DMClone_Forest; 1601 dm->ops->setfromoptions = DMSetFromOptions_Forest; 1602 dm->ops->destroy = DMDestroy_Forest; 1603 dm->ops->createsubdm = DMCreateSubDM_Forest; 1604 dm->ops->refine = DMRefine_Forest; 1605 dm->ops->coarsen = DMCoarsen_Forest; 1606 PetscFunctionReturn(PETSC_SUCCESS); 1607 } 1608 1609 /*MC 1610 1611 DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh. Forests usually have a base `DM` 1612 (see `DMForestGetBaseDM()`), from which it is refined. The refinement and partitioning of forests is considered 1613 immutable after `DMSetUp()` is called. To adapt a mesh, one should call `DMForestTemplate()` to create a new mesh that 1614 will default to being identical to it, specify how that mesh should differ, and then calling `DMSetUp()` on the new 1615 mesh. 1616 1617 To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the 1618 previous mesh whose values indicate which cells should be refined (`DM_ADAPT_REFINE`) or coarsened (`DM_ADAPT_COARSEN`) 1619 and how (subtypes are free to allow additional values for things like anisotropic refinement). The label should be 1620 given to the *new* mesh with `DMForestSetAdaptivityLabel()`. 1621 1622 Level: advanced 1623 1624 .seealso: `DMType`, `DM`, `DMCreate()`, `DMSetType()`, `DMForestGetBaseDM()`, `DMForestSetBaseDM()`, `DMForestTemplate()`, `DMForestSetAdaptivityLabel()` 1625 M*/ 1626 1627 PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm) 1628 { 1629 DM_Forest *forest; 1630 1631 PetscFunctionBegin; 1632 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1633 PetscCall(PetscNew(&forest)); 1634 dm->dim = 0; 1635 dm->data = forest; 1636 forest->refct = 1; 1637 forest->data = NULL; 1638 forest->topology = NULL; 1639 forest->adapt = NULL; 1640 forest->base = NULL; 1641 forest->adaptPurpose = DM_ADAPT_DETERMINE; 1642 forest->adjDim = PETSC_DEFAULT; 1643 forest->overlap = PETSC_DEFAULT; 1644 forest->minRefinement = PETSC_DEFAULT; 1645 forest->maxRefinement = PETSC_DEFAULT; 1646 forest->initRefinement = PETSC_DEFAULT; 1647 forest->cStart = PETSC_DETERMINE; 1648 forest->cEnd = PETSC_DETERMINE; 1649 forest->cellSF = NULL; 1650 forest->adaptLabel = NULL; 1651 forest->gradeFactor = 2; 1652 forest->cellWeights = NULL; 1653 forest->cellWeightsCopyMode = PETSC_USE_POINTER; 1654 forest->weightsFactor = 1.; 1655 forest->weightCapacity = 1.; 1656 PetscCall(DMForestSetAdaptivityStrategy(dm, DMFORESTADAPTALL)); 1657 PetscCall(DMInitialize_Forest(dm)); 1658 PetscFunctionReturn(PETSC_SUCCESS); 1659 } 1660