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