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