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