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