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