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