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