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