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