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