19be51f97SToby Isaac #include <petsc/private/dmforestimpl.h> /*I "petscdmforest.h" I*/ 29be51f97SToby Isaac #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 3a1b0c543SToby Isaac #include <petsc/private/dmlabelimpl.h> /*I "petscdmlabel.h" I*/ 4ef19d27cSToby Isaac #include <petscsf.h> 5db4d5e8cSToby Isaac 627d4645fSToby Isaac PetscBool DMForestPackageInitialized = PETSC_FALSE; 727d4645fSToby Isaac 827d4645fSToby Isaac typedef struct _DMForestTypeLink *DMForestTypeLink; 927d4645fSToby Isaac 109371c9d4SSatish Balay struct _DMForestTypeLink { 1127d4645fSToby Isaac char *name; 1227d4645fSToby Isaac DMForestTypeLink next; 1327d4645fSToby Isaac }; 1427d4645fSToby Isaac 1527d4645fSToby Isaac DMForestTypeLink DMForestTypeList; 1627d4645fSToby Isaac 17d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMForestPackageFinalize(void) 18d71ae5a4SJacob Faibussowitsch { 1927d4645fSToby Isaac DMForestTypeLink oldLink, link = DMForestTypeList; 2027d4645fSToby Isaac 2127d4645fSToby Isaac PetscFunctionBegin; 2227d4645fSToby Isaac while (link) { 2327d4645fSToby Isaac oldLink = link; 249566063dSJacob Faibussowitsch PetscCall(PetscFree(oldLink->name)); 2527d4645fSToby Isaac link = oldLink->next; 269566063dSJacob Faibussowitsch PetscCall(PetscFree(oldLink)); 2727d4645fSToby Isaac } 283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2927d4645fSToby Isaac } 3027d4645fSToby Isaac 31d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMForestPackageInitialize(void) 32d71ae5a4SJacob Faibussowitsch { 3327d4645fSToby Isaac PetscFunctionBegin; 343ba16761SJacob Faibussowitsch if (DMForestPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS); 3527d4645fSToby Isaac DMForestPackageInitialized = PETSC_TRUE; 36f885a11aSToby Isaac 379566063dSJacob Faibussowitsch PetscCall(DMForestRegisterType(DMFOREST)); 389566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(DMForestPackageFinalize)); 393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4027d4645fSToby Isaac } 4127d4645fSToby Isaac 429be51f97SToby Isaac /*@C 43dce8aebaSBarry Smith DMForestRegisterType - Registers a `DMType` as a subtype of `DMFOREST` (so that `DMIsForest()` will be correct) 449be51f97SToby Isaac 459be51f97SToby Isaac Not Collective 469be51f97SToby Isaac 4760225df5SJacob Faibussowitsch Input Parameter: 489be51f97SToby Isaac . name - the name of the type 499be51f97SToby Isaac 509be51f97SToby Isaac Level: advanced 519be51f97SToby Isaac 52db781477SPatrick Sanan .seealso: `DMFOREST`, `DMIsForest()` 539be51f97SToby Isaac @*/ 54d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestRegisterType(DMType name) 55d71ae5a4SJacob Faibussowitsch { 5627d4645fSToby Isaac DMForestTypeLink link; 5727d4645fSToby Isaac 5827d4645fSToby Isaac PetscFunctionBegin; 599566063dSJacob Faibussowitsch PetscCall(DMForestPackageInitialize()); 609566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 619566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name, &link->name)); 6227d4645fSToby Isaac link->next = DMForestTypeList; 6327d4645fSToby Isaac DMForestTypeList = link; 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6527d4645fSToby Isaac } 6627d4645fSToby Isaac 679be51f97SToby Isaac /*@ 689be51f97SToby Isaac DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes 699be51f97SToby Isaac 709be51f97SToby Isaac Not Collective 719be51f97SToby Isaac 7260225df5SJacob Faibussowitsch Input Parameter: 739be51f97SToby Isaac . dm - the DM object 749be51f97SToby Isaac 7560225df5SJacob Faibussowitsch Output Parameter: 769be51f97SToby Isaac . isForest - whether dm is a subtype of DMFOREST 779be51f97SToby Isaac 789be51f97SToby Isaac Level: intermediate 799be51f97SToby Isaac 80db781477SPatrick Sanan .seealso: `DMFOREST`, `DMForestRegisterType()` 819be51f97SToby Isaac @*/ 82d71ae5a4SJacob Faibussowitsch PetscErrorCode DMIsForest(DM dm, PetscBool *isForest) 83d71ae5a4SJacob Faibussowitsch { 8427d4645fSToby Isaac DMForestTypeLink link = DMForestTypeList; 8527d4645fSToby Isaac 8627d4645fSToby Isaac PetscFunctionBegin; 8727d4645fSToby Isaac while (link) { 8827d4645fSToby Isaac PetscBool sameType; 899566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm, link->name, &sameType)); 9027d4645fSToby Isaac if (sameType) { 9127d4645fSToby Isaac *isForest = PETSC_TRUE; 923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9327d4645fSToby Isaac } 9427d4645fSToby Isaac link = link->next; 9527d4645fSToby Isaac } 9627d4645fSToby Isaac *isForest = PETSC_FALSE; 973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9827d4645fSToby Isaac } 9927d4645fSToby Isaac 1009be51f97SToby Isaac /*@ 101a4e35b19SJacob Faibussowitsch DMForestTemplate - Create a new `DM` that will be adapted from a source `DM`. 1029be51f97SToby Isaac 10320f4b53cSBarry Smith Collective 1049be51f97SToby Isaac 1059be51f97SToby Isaac Input Parameters: 10620f4b53cSBarry Smith + dm - the source `DM` object 10720f4b53cSBarry Smith - comm - the communicator for the new `DM` (this communicator is currently ignored, but is present so that `DMForestTemplate()` can be used within `DMCoarsen()`) 1089be51f97SToby Isaac 1099be51f97SToby Isaac Output Parameter: 11020f4b53cSBarry Smith . tdm - the new `DM` object 1119be51f97SToby Isaac 1129be51f97SToby Isaac Level: intermediate 1139be51f97SToby Isaac 114a4e35b19SJacob Faibussowitsch Notes: 115a4e35b19SJacob Faibussowitsch The new `DM` reproduces the configuration of the source, but is not yet setup, so that the 116a4e35b19SJacob Faibussowitsch user can then define only the ways that the new `DM` should differ (by, e.g., refinement or 117a4e35b19SJacob Faibussowitsch repartitioning). The source `DM` is also set as the adaptivity source `DM` of the new `DM` 118a4e35b19SJacob Faibussowitsch (see `DMForestSetAdaptivityForest()`). 119a4e35b19SJacob Faibussowitsch 12020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetAdaptivityForest()` 1219be51f97SToby Isaac @*/ 122d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm) 123d71ae5a4SJacob Faibussowitsch { 124a0452a8eSToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 12520e8089bSToby Isaac DMType type; 126a0452a8eSToby Isaac DM base; 127a0452a8eSToby Isaac DMForestTopology topology; 12805e99e11SStefano Zampini MatType mtype; 129a0452a8eSToby Isaac PetscInt dim, overlap, ref, factor; 130a0452a8eSToby Isaac DMForestAdaptivityStrategy strat; 131795844e7SToby Isaac void *ctx; 13249fc9a2fSToby Isaac PetscErrorCode (*map)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void *); 1333e58adeeSToby Isaac void *mapCtx; 134a0452a8eSToby Isaac 135a0452a8eSToby Isaac PetscFunctionBegin; 136a0452a8eSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1379566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), tdm)); 1389566063dSJacob Faibussowitsch PetscCall(DMGetType(dm, &type)); 1399566063dSJacob Faibussowitsch PetscCall(DMSetType(*tdm, type)); 1409566063dSJacob Faibussowitsch PetscCall(DMForestGetBaseDM(dm, &base)); 1419566063dSJacob Faibussowitsch PetscCall(DMForestSetBaseDM(*tdm, base)); 1429566063dSJacob Faibussowitsch PetscCall(DMForestGetTopology(dm, &topology)); 1439566063dSJacob Faibussowitsch PetscCall(DMForestSetTopology(*tdm, topology)); 1449566063dSJacob Faibussowitsch PetscCall(DMForestGetAdjacencyDimension(dm, &dim)); 1459566063dSJacob Faibussowitsch PetscCall(DMForestSetAdjacencyDimension(*tdm, dim)); 1469566063dSJacob Faibussowitsch PetscCall(DMForestGetPartitionOverlap(dm, &overlap)); 1479566063dSJacob Faibussowitsch PetscCall(DMForestSetPartitionOverlap(*tdm, overlap)); 1489566063dSJacob Faibussowitsch PetscCall(DMForestGetMinimumRefinement(dm, &ref)); 1499566063dSJacob Faibussowitsch PetscCall(DMForestSetMinimumRefinement(*tdm, ref)); 1509566063dSJacob Faibussowitsch PetscCall(DMForestGetMaximumRefinement(dm, &ref)); 1519566063dSJacob Faibussowitsch PetscCall(DMForestSetMaximumRefinement(*tdm, ref)); 1529566063dSJacob Faibussowitsch PetscCall(DMForestGetAdaptivityStrategy(dm, &strat)); 1539566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityStrategy(*tdm, strat)); 1549566063dSJacob Faibussowitsch PetscCall(DMForestGetGradeFactor(dm, &factor)); 1559566063dSJacob Faibussowitsch PetscCall(DMForestSetGradeFactor(*tdm, factor)); 1569566063dSJacob Faibussowitsch PetscCall(DMForestGetBaseCoordinateMapping(dm, &map, &mapCtx)); 1579566063dSJacob Faibussowitsch PetscCall(DMForestSetBaseCoordinateMapping(*tdm, map, mapCtx)); 1581baa6e33SBarry Smith if (forest->ftemplate) PetscCall((*forest->ftemplate)(dm, *tdm)); 1599566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityForest(*tdm, dm)); 1609566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm, *tdm)); 1619566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(dm, &ctx)); 1629566063dSJacob Faibussowitsch PetscCall(DMSetApplicationContext(*tdm, &ctx)); 16390b157c4SStefano Zampini { 1644fb89dddSMatthew G. Knepley const PetscReal *maxCell, *L, *Lstart; 165795844e7SToby Isaac 1664fb89dddSMatthew G. Knepley PetscCall(DMGetPeriodicity(dm, &maxCell, &Lstart, &L)); 1674fb89dddSMatthew G. Knepley PetscCall(DMSetPeriodicity(*tdm, maxCell, Lstart, L)); 168795844e7SToby Isaac } 1699566063dSJacob Faibussowitsch PetscCall(DMGetMatType(dm, &mtype)); 1709566063dSJacob Faibussowitsch PetscCall(DMSetMatType(*tdm, mtype)); 1713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 172a0452a8eSToby Isaac } 173a0452a8eSToby Isaac 17401d9d024SToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm); 17501d9d024SToby Isaac 176d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm) 177d71ae5a4SJacob Faibussowitsch { 178db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 179db4d5e8cSToby Isaac const char *type; 180db4d5e8cSToby Isaac 181db4d5e8cSToby Isaac PetscFunctionBegin; 182db4d5e8cSToby Isaac forest->refct++; 183db4d5e8cSToby Isaac (*newdm)->data = forest; 1849566063dSJacob Faibussowitsch PetscCall(PetscObjectGetType((PetscObject)dm, &type)); 1859566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject)*newdm, type)); 1869566063dSJacob Faibussowitsch PetscCall(DMInitialize_Forest(*newdm)); 1873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 188db4d5e8cSToby Isaac } 189db4d5e8cSToby Isaac 190d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDestroy_Forest(DM dm) 191d71ae5a4SJacob Faibussowitsch { 192db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 193db4d5e8cSToby Isaac 194db4d5e8cSToby Isaac PetscFunctionBegin; 19597ca8683SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMConvert_plex_p4est_C", NULL)); 19697ca8683SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMConvert_p4est_plex_C", NULL)); 19797ca8683SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMConvert_plex_p8est_C", NULL)); 19897ca8683SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMConvert_p8est_plex_C", NULL)); 19997ca8683SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMSetUpGLVisViewer_C", NULL)); 20097ca8683SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMCreateNeumannOverlap_C", NULL)); 20197ca8683SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)dm, "DMPlexGetOverlap_C", NULL)); 20297ca8683SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)dm, "MatComputeNeumannOverlap_C", NULL)); 2033ba16761SJacob Faibussowitsch if (--forest->refct > 0) PetscFunctionReturn(PETSC_SUCCESS); 2049566063dSJacob Faibussowitsch if (forest->destroy) PetscCall((*forest->destroy)(dm)); 2059566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&forest->cellSF)); 2069566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&forest->preCoarseToFine)); 2079566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&forest->coarseToPreFine)); 2089566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&forest->adaptLabel)); 2099566063dSJacob Faibussowitsch PetscCall(PetscFree(forest->adaptStrategy)); 2109566063dSJacob Faibussowitsch PetscCall(DMDestroy(&forest->base)); 2119566063dSJacob Faibussowitsch PetscCall(DMDestroy(&forest->adapt)); 2129566063dSJacob Faibussowitsch PetscCall(PetscFree(forest->topology)); 2139566063dSJacob Faibussowitsch PetscCall(PetscFree(forest)); 2143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 215db4d5e8cSToby Isaac } 216db4d5e8cSToby Isaac 2179be51f97SToby Isaac /*@C 21820f4b53cSBarry Smith DMForestSetTopology - Set the topology of a `DMFOREST` during the pre-setup phase. The topology is a string (e.g. 21920f4b53cSBarry Smith "cube", "shell") and can be interpreted by subtypes of `DMFOREST`) to construct the base DM of a forest during 22020f4b53cSBarry Smith `DMSetUp()`. 2219be51f97SToby Isaac 22220f4b53cSBarry Smith Logically collectiv 2239be51f97SToby Isaac 22460225df5SJacob Faibussowitsch Input Parameters: 2259be51f97SToby Isaac + dm - the forest 2269be51f97SToby Isaac - topology - the topology of the forest 2279be51f97SToby Isaac 2289be51f97SToby Isaac Level: intermediate 2299be51f97SToby Isaac 23020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetTopology()`, `DMForestSetBaseDM()` 2319be51f97SToby Isaac @*/ 232d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology) 233d71ae5a4SJacob Faibussowitsch { 234db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 235db4d5e8cSToby Isaac 236db4d5e8cSToby Isaac PetscFunctionBegin; 237db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 23828b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the topology after setup"); 2399566063dSJacob Faibussowitsch PetscCall(PetscFree(forest->topology)); 2409566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy((const char *)topology, (char **)&forest->topology)); 2413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 242db4d5e8cSToby Isaac } 243db4d5e8cSToby Isaac 2449be51f97SToby Isaac /*@C 24520f4b53cSBarry Smith DMForestGetTopology - Get a string describing the topology of a `DMFOREST`. 2469be51f97SToby Isaac 24720f4b53cSBarry Smith Not Collective 2489be51f97SToby Isaac 24960225df5SJacob Faibussowitsch Input Parameter: 2509be51f97SToby Isaac . dm - the forest 2519be51f97SToby Isaac 25260225df5SJacob Faibussowitsch Output Parameter: 2539be51f97SToby Isaac . topology - the topology of the forest (e.g., 'cube', 'shell') 2549be51f97SToby Isaac 2559be51f97SToby Isaac Level: intermediate 2569be51f97SToby Isaac 25720f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetTopology()` 2589be51f97SToby Isaac @*/ 259d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology) 260d71ae5a4SJacob Faibussowitsch { 261dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 262dd8e54a2SToby Isaac 263dd8e54a2SToby Isaac PetscFunctionBegin; 264dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2654f572ea9SToby Isaac PetscAssertPointer(topology, 2); 266dd8e54a2SToby Isaac *topology = forest->topology; 2673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 268dd8e54a2SToby Isaac } 269dd8e54a2SToby Isaac 2709be51f97SToby Isaac /*@ 271a4e35b19SJacob Faibussowitsch DMForestSetBaseDM - During the pre-setup phase, set the `DM` that defines the base mesh of a 272a4e35b19SJacob Faibussowitsch `DMFOREST` forest. 2739be51f97SToby Isaac 27420f4b53cSBarry Smith Logically Collective 2759be51f97SToby Isaac 2769be51f97SToby Isaac Input Parameters: 2779be51f97SToby Isaac + dm - the forest 27820f4b53cSBarry Smith - base - the base `DM` of the forest 279765b024eSBarry Smith 2809be51f97SToby Isaac Level: intermediate 2819be51f97SToby Isaac 282a4e35b19SJacob Faibussowitsch Notes: 283a4e35b19SJacob Faibussowitsch The forest will be hierarchically refined from the base, and all refinements/coarsenings of 284a4e35b19SJacob Faibussowitsch the forest will share its base. In general, two forest must share a base to be comparable, 285a4e35b19SJacob Faibussowitsch to do things like construct interpolators. 286a4e35b19SJacob Faibussowitsch 28720f4b53cSBarry Smith Currently the base `DM` must be a `DMPLEX` 28820f4b53cSBarry Smith 28920f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetBaseDM()` 2909be51f97SToby Isaac @*/ 291d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetBaseDM(DM dm, DM base) 292d71ae5a4SJacob Faibussowitsch { 293dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 294dd8e54a2SToby Isaac PetscInt dim, dimEmbed; 295dd8e54a2SToby Isaac 296dd8e54a2SToby Isaac PetscFunctionBegin; 297dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 29828b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the base after setup"); 2999566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)base)); 3009566063dSJacob Faibussowitsch PetscCall(DMDestroy(&forest->base)); 301dd8e54a2SToby Isaac forest->base = base; 302a0452a8eSToby Isaac if (base) { 3034fb89dddSMatthew G. Knepley const PetscReal *maxCell, *Lstart, *L; 30428dfcf7cSStefano Zampini 305a0452a8eSToby Isaac PetscValidHeaderSpecific(base, DM_CLASSID, 2); 3069566063dSJacob Faibussowitsch PetscCall(DMGetDimension(base, &dim)); 3079566063dSJacob Faibussowitsch PetscCall(DMSetDimension(dm, dim)); 3089566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(base, &dimEmbed)); 3099566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(dm, dimEmbed)); 3104fb89dddSMatthew G. Knepley PetscCall(DMGetPeriodicity(base, &maxCell, &Lstart, &L)); 3114fb89dddSMatthew G. Knepley PetscCall(DMSetPeriodicity(dm, maxCell, Lstart, L)); 3124fb89dddSMatthew G. Knepley } else PetscCall(DMSetPeriodicity(dm, NULL, NULL, NULL)); 3133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 314dd8e54a2SToby Isaac } 315dd8e54a2SToby Isaac 3169be51f97SToby Isaac /*@ 317a4e35b19SJacob Faibussowitsch DMForestGetBaseDM - Get the base `DM` of a `DMFOREST` 3189be51f97SToby Isaac 31920f4b53cSBarry Smith Not Collective 3209be51f97SToby Isaac 3219be51f97SToby Isaac Input Parameter: 3229be51f97SToby Isaac . dm - the forest 3239be51f97SToby Isaac 3249be51f97SToby Isaac Output Parameter: 325a4e35b19SJacob Faibussowitsch . base - the base `DM` of the forest 326a4e35b19SJacob Faibussowitsch 327a4e35b19SJacob Faibussowitsch Level: intermediate 3289be51f97SToby Isaac 329367003a6SStefano Zampini Notes: 330367003a6SStefano Zampini After DMSetUp(), the base DM will be redundantly distributed across MPI processes 331367003a6SStefano Zampini 332a4e35b19SJacob Faibussowitsch The forest will be hierarchically refined from the base, and all refinements/coarsenings of 333a4e35b19SJacob Faibussowitsch the forest will share its base. In general, two forest must share a base to be comparable, 334a4e35b19SJacob Faibussowitsch to do things like construct interpolators. 3359be51f97SToby Isaac 33620f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetBaseDM()` 3379be51f97SToby Isaac @*/ 338d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetBaseDM(DM dm, DM *base) 339d71ae5a4SJacob Faibussowitsch { 340dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 341dd8e54a2SToby Isaac 342dd8e54a2SToby Isaac PetscFunctionBegin; 343dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3444f572ea9SToby Isaac PetscAssertPointer(base, 2); 345dd8e54a2SToby Isaac *base = forest->base; 3463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 347dd8e54a2SToby Isaac } 348dd8e54a2SToby Isaac 349d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetBaseCoordinateMapping(DM dm, PetscErrorCode (*func)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void *), void *ctx) 350d71ae5a4SJacob Faibussowitsch { 351cf38a08cSToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 352cf38a08cSToby Isaac 353cf38a08cSToby Isaac PetscFunctionBegin; 354cf38a08cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 355cf38a08cSToby Isaac forest->mapcoordinates = func; 356cf38a08cSToby Isaac forest->mapcoordinatesctx = ctx; 3573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 358cf38a08cSToby Isaac } 359cf38a08cSToby Isaac 360d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetBaseCoordinateMapping(DM dm, PetscErrorCode (**func)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void *), void *ctx) 361d71ae5a4SJacob Faibussowitsch { 362cf38a08cSToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 363cf38a08cSToby Isaac 364cf38a08cSToby Isaac PetscFunctionBegin; 365cf38a08cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 366cf38a08cSToby Isaac if (func) *func = forest->mapcoordinates; 367cf38a08cSToby Isaac if (ctx) *((void **)ctx) = forest->mapcoordinatesctx; 3683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 369cf38a08cSToby Isaac } 370cf38a08cSToby Isaac 3719be51f97SToby Isaac /*@ 372a4e35b19SJacob Faibussowitsch DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the 373a4e35b19SJacob Faibussowitsch current forest will be adapted (e.g., the current forest will be 374a4e35b19SJacob Faibussowitsch refined/coarsened/repartitioned from it) in `DMSetUp()`. 3759be51f97SToby Isaac 37620f4b53cSBarry Smith Logically Collective 3779be51f97SToby Isaac 378d8d19677SJose E. Roman Input Parameters: 3799be51f97SToby Isaac + dm - the new forest, which will be constructed from adapt 3809be51f97SToby Isaac - adapt - the old forest 3819be51f97SToby Isaac 3829be51f97SToby Isaac Level: intermediate 3839be51f97SToby Isaac 38420f4b53cSBarry Smith Note: 385a4e35b19SJacob Faibussowitsch Usually not needed by users directly, `DMForestTemplate()` constructs a new forest to be 386a4e35b19SJacob Faibussowitsch adapted from an old forest and calls this routine. 387a4e35b19SJacob Faibussowitsch 388a4e35b19SJacob Faibussowitsch This can be called after setup with `adapt` = `NULL`, which will clear all internal data 389a4e35b19SJacob Faibussowitsch related to the adaptivity forest from `dm`. This way, repeatedly adapting does not leave 390a4e35b19SJacob Faibussowitsch stale `DM` objects in memory. 39120f4b53cSBarry Smith 39220f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetAdaptivityForest()`, `DMForestSetAdaptivityPurpose()` 3939be51f97SToby Isaac @*/ 394d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityForest(DM dm, DM adapt) 395d71ae5a4SJacob Faibussowitsch { 396dffe73a3SToby Isaac DM_Forest *forest, *adaptForest, *oldAdaptForest; 397dffe73a3SToby Isaac DM oldAdapt; 398456cc5b7SMatthew G. Knepley PetscBool isForest; 399dd8e54a2SToby Isaac 400dd8e54a2SToby Isaac PetscFunctionBegin; 401dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4021fd50544SStefano Zampini if (adapt) PetscValidHeaderSpecific(adapt, DM_CLASSID, 2); 4039566063dSJacob Faibussowitsch PetscCall(DMIsForest(dm, &isForest)); 4043ba16761SJacob Faibussowitsch if (!isForest) PetscFunctionReturn(PETSC_SUCCESS); 4051dca8a05SBarry Smith PetscCheck(adapt == NULL || !dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the adaptation forest after setup"); 406ba936b91SToby Isaac forest = (DM_Forest *)dm->data; 4079566063dSJacob Faibussowitsch PetscCall(DMForestGetAdaptivityForest(dm, &oldAdapt)); 408193eb951SToby Isaac adaptForest = (DM_Forest *)(adapt ? adapt->data : NULL); 409193eb951SToby Isaac oldAdaptForest = (DM_Forest *)(oldAdapt ? oldAdapt->data : NULL); 410dffe73a3SToby Isaac if (adaptForest != oldAdaptForest) { 4119566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&forest->preCoarseToFine)); 4129566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&forest->coarseToPreFine)); 4139566063dSJacob Faibussowitsch if (forest->clearadaptivityforest) PetscCall((*forest->clearadaptivityforest)(dm)); 414dffe73a3SToby Isaac } 41526d9498aSToby Isaac switch (forest->adaptPurpose) { 416cd3c525cSToby Isaac case DM_ADAPT_DETERMINE: 4179566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)adapt)); 418*f4f49eeaSPierre Jolivet PetscCall(DMDestroy(&forest->adapt)); 419ba936b91SToby Isaac forest->adapt = adapt; 42026d9498aSToby Isaac break; 421d71ae5a4SJacob Faibussowitsch case DM_ADAPT_REFINE: 422d71ae5a4SJacob Faibussowitsch PetscCall(DMSetCoarseDM(dm, adapt)); 423d71ae5a4SJacob Faibussowitsch break; 424a1b0c543SToby Isaac case DM_ADAPT_COARSEN: 425d71ae5a4SJacob Faibussowitsch case DM_ADAPT_COARSEN_LAST: 426d71ae5a4SJacob Faibussowitsch PetscCall(DMSetFineDM(dm, adapt)); 427d71ae5a4SJacob Faibussowitsch break; 428d71ae5a4SJacob Faibussowitsch default: 429d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid adaptivity purpose"); 43026d9498aSToby Isaac } 4313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 432dd8e54a2SToby Isaac } 433dd8e54a2SToby Isaac 4349be51f97SToby Isaac /*@ 4359be51f97SToby Isaac DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted. 4369be51f97SToby Isaac 43720f4b53cSBarry Smith Not Collective 4389be51f97SToby Isaac 4399be51f97SToby Isaac Input Parameter: 4409be51f97SToby Isaac . dm - the forest 4419be51f97SToby Isaac 4429be51f97SToby Isaac Output Parameter: 44320f4b53cSBarry Smith . adapt - the forest from which `dm` is/was adapted 4449be51f97SToby Isaac 4459be51f97SToby Isaac Level: intermediate 4469be51f97SToby Isaac 44720f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetAdaptivityForest()`, `DMForestSetAdaptivityPurpose()` 4489be51f97SToby Isaac @*/ 449d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt) 450d71ae5a4SJacob Faibussowitsch { 451ba936b91SToby Isaac DM_Forest *forest; 452dd8e54a2SToby Isaac 453dd8e54a2SToby Isaac PetscFunctionBegin; 454dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 455ba936b91SToby Isaac forest = (DM_Forest *)dm->data; 45626d9498aSToby Isaac switch (forest->adaptPurpose) { 457d71ae5a4SJacob Faibussowitsch case DM_ADAPT_DETERMINE: 458d71ae5a4SJacob Faibussowitsch *adapt = forest->adapt; 459d71ae5a4SJacob Faibussowitsch break; 460d71ae5a4SJacob Faibussowitsch case DM_ADAPT_REFINE: 461d71ae5a4SJacob Faibussowitsch PetscCall(DMGetCoarseDM(dm, adapt)); 462d71ae5a4SJacob Faibussowitsch break; 463a1b0c543SToby Isaac case DM_ADAPT_COARSEN: 464d71ae5a4SJacob Faibussowitsch case DM_ADAPT_COARSEN_LAST: 465d71ae5a4SJacob Faibussowitsch PetscCall(DMGetFineDM(dm, adapt)); 466d71ae5a4SJacob Faibussowitsch break; 467d71ae5a4SJacob Faibussowitsch default: 468d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "invalid adaptivity purpose"); 46926d9498aSToby Isaac } 4703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 47126d9498aSToby Isaac } 47226d9498aSToby Isaac 4739be51f97SToby Isaac /*@ 47420f4b53cSBarry Smith DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current `DM` is being adapted from its 47520f4b53cSBarry Smith source (set with `DMForestSetAdaptivityForest()`) for the purpose of refinement (`DM_ADAPT_REFINE`), coarsening 476a4e35b19SJacob Faibussowitsch (`DM_ADAPT_COARSEN`), or undefined (`DM_ADAPT_DETERMINE`). 4779be51f97SToby Isaac 47820f4b53cSBarry Smith Logically Collective 4799be51f97SToby Isaac 4809be51f97SToby Isaac Input Parameters: 4819be51f97SToby Isaac + dm - the forest 482bf2d5fbbSStefano Zampini - purpose - the adaptivity purpose 4839be51f97SToby Isaac 4849be51f97SToby Isaac Level: advanced 4859be51f97SToby Isaac 486a4e35b19SJacob Faibussowitsch Notes: 487a4e35b19SJacob Faibussowitsch This only matters for reference counting during `DMDestroy()`. Cyclic references 488a4e35b19SJacob Faibussowitsch can be found between `DM`s only if the cyclic reference is due to a fine/coarse relationship 489a4e35b19SJacob Faibussowitsch (see `DMSetFineDM()`/`DMSetCoarseDM()`). If the purpose is not refinement or coarsening, and 490a4e35b19SJacob Faibussowitsch the user does not maintain a reference to the post-adaptation forest (i.e., the one created 491a4e35b19SJacob Faibussowitsch by `DMForestTemplate()`), this can cause a memory leak. This method is used by subtypes 492a4e35b19SJacob Faibussowitsch of `DMFOREST` when automatically constructing mesh hierarchies. 493a4e35b19SJacob Faibussowitsch 49420f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestTemplate()`, `DMForestSetAdaptivityForest()`, `DMForestGetAdaptivityForest()`, `DMAdaptFlag` 4959be51f97SToby Isaac @*/ 496d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMAdaptFlag purpose) 497d71ae5a4SJacob Faibussowitsch { 49826d9498aSToby Isaac DM_Forest *forest; 49926d9498aSToby Isaac 50026d9498aSToby Isaac PetscFunctionBegin; 50126d9498aSToby Isaac forest = (DM_Forest *)dm->data; 50228b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the adaptation forest after setup"); 50326d9498aSToby Isaac if (purpose != forest->adaptPurpose) { 50426d9498aSToby Isaac DM adapt; 50526d9498aSToby Isaac 5069566063dSJacob Faibussowitsch PetscCall(DMForestGetAdaptivityForest(dm, &adapt)); 5079566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)adapt)); 5089566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityForest(dm, NULL)); 509f885a11aSToby Isaac 51026d9498aSToby Isaac forest->adaptPurpose = purpose; 511f885a11aSToby Isaac 5129566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityForest(dm, adapt)); 5139566063dSJacob Faibussowitsch PetscCall(DMDestroy(&adapt)); 51426d9498aSToby Isaac } 5153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 516dd8e54a2SToby Isaac } 517dd8e54a2SToby Isaac 51856c0450aSToby Isaac /*@ 51920f4b53cSBarry Smith DMForestGetAdaptivityPurpose - Get whether the current `DM` is being adapted from its source (set with 52020f4b53cSBarry Smith `DMForestSetAdaptivityForest()`) for the purpose of refinement (`DM_ADAPT_REFINE`), coarsening (`DM_ADAPT_COARSEN`), 52120f4b53cSBarry Smith coarsening only the last level (`DM_ADAPT_COARSEN_LAST`) or undefined (`DM_ADAPT_DETERMINE`). 52256c0450aSToby Isaac 52320f4b53cSBarry Smith Not Collective 52456c0450aSToby Isaac 52556c0450aSToby Isaac Input Parameter: 52656c0450aSToby Isaac . dm - the forest 52756c0450aSToby Isaac 52856c0450aSToby Isaac Output Parameter: 529bf2d5fbbSStefano Zampini . purpose - the adaptivity purpose 53056c0450aSToby Isaac 53156c0450aSToby Isaac Level: advanced 53256c0450aSToby Isaac 533a4e35b19SJacob Faibussowitsch Notes: 534a4e35b19SJacob Faibussowitsch This only matters for reference counting: during `DMDestroy()`. Cyclic references 535a4e35b19SJacob Faibussowitsch can be found between `DM`s only if the cyclic reference is due to a fine/coarse relationship 536a4e35b19SJacob Faibussowitsch (See `DMSetFineDM()`/`DMSetCoarseDM()`). If the purpose is not refinement or coarsening, and 537a4e35b19SJacob Faibussowitsch the user does not maintain a reference to the post-adaptation forest (i.e., the one created 538a4e35b19SJacob Faibussowitsch by `DMForestTemplate()`), this can cause a memory leak. This method is used by subtypes 539a4e35b19SJacob Faibussowitsch of `DMFOREST` when automatically constructing mesh hierarchies. 540a4e35b19SJacob Faibussowitsch 54120f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestTemplate()`, `DMForestSetAdaptivityForest()`, `DMForestGetAdaptivityForest()`, `DMAdaptFlag` 54256c0450aSToby Isaac @*/ 543d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMAdaptFlag *purpose) 544d71ae5a4SJacob Faibussowitsch { 54556c0450aSToby Isaac DM_Forest *forest; 54656c0450aSToby Isaac 54756c0450aSToby Isaac PetscFunctionBegin; 54856c0450aSToby Isaac forest = (DM_Forest *)dm->data; 54956c0450aSToby Isaac *purpose = forest->adaptPurpose; 5503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 55156c0450aSToby Isaac } 55256c0450aSToby Isaac 5539be51f97SToby Isaac /*@ 5549be51f97SToby Isaac DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine 5559be51f97SToby Isaac cell adjacency (for the purposes of partitioning and overlap). 5569be51f97SToby Isaac 55720f4b53cSBarry Smith Logically Collective 5589be51f97SToby Isaac 5599be51f97SToby Isaac Input Parameters: 5609be51f97SToby Isaac + dm - the forest 5619be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency) 5629be51f97SToby Isaac 5639be51f97SToby Isaac Level: intermediate 5649be51f97SToby Isaac 56520f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()`, `DMForestSetPartitionOverlap()` 5669be51f97SToby Isaac @*/ 567d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim) 568d71ae5a4SJacob Faibussowitsch { 569dd8e54a2SToby Isaac PetscInt dim; 570dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 571dd8e54a2SToby Isaac 572dd8e54a2SToby Isaac PetscFunctionBegin; 573dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 57428b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the adjacency dimension after setup"); 57563a3b9bcSJacob Faibussowitsch PetscCheck(adjDim >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "adjacency dim cannot be < 0: %" PetscInt_FMT, adjDim); 5769566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 57763a3b9bcSJacob Faibussowitsch PetscCheck(adjDim <= dim, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "adjacency dim cannot be > %" PetscInt_FMT ": %" PetscInt_FMT, dim, adjDim); 578dd8e54a2SToby Isaac forest->adjDim = adjDim; 5793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 580dd8e54a2SToby Isaac } 581dd8e54a2SToby Isaac 5829be51f97SToby Isaac /*@ 58320f4b53cSBarry Smith DMForestSetAdjacencyCodimension - Like `DMForestSetAdjacencyDimension()`, but specified as a co-dimension (so that, 5849be51f97SToby Isaac e.g., adjacency based on facets can be specified by codimension 1 in all cases) 5859be51f97SToby Isaac 58620f4b53cSBarry Smith Logically Collective 5879be51f97SToby Isaac 5889be51f97SToby Isaac Input Parameters: 5899be51f97SToby Isaac + dm - the forest 59020f4b53cSBarry Smith - adjCodim - default is the dimension of the forest (see `DMGetDimension()`), since this is the codimension of vertices 5919be51f97SToby Isaac 5929be51f97SToby Isaac Level: intermediate 5939be51f97SToby Isaac 59420f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetAdjacencyCodimension()`, `DMForestSetAdjacencyDimension()` 5959be51f97SToby Isaac @*/ 596d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim) 597d71ae5a4SJacob Faibussowitsch { 598dd8e54a2SToby Isaac PetscInt dim; 599dd8e54a2SToby Isaac 600dd8e54a2SToby Isaac PetscFunctionBegin; 601dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 6029566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 6039566063dSJacob Faibussowitsch PetscCall(DMForestSetAdjacencyDimension(dm, dim - adjCodim)); 6043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 605dd8e54a2SToby Isaac } 606dd8e54a2SToby Isaac 6079be51f97SToby Isaac /*@ 6089be51f97SToby Isaac DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the 6099be51f97SToby Isaac purposes of partitioning and overlap). 6109be51f97SToby Isaac 61120f4b53cSBarry Smith Not Collective 6129be51f97SToby Isaac 6139be51f97SToby Isaac Input Parameter: 6149be51f97SToby Isaac . dm - the forest 6159be51f97SToby Isaac 6169be51f97SToby Isaac Output Parameter: 6179be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency) 6189be51f97SToby Isaac 6199be51f97SToby Isaac Level: intermediate 6209be51f97SToby Isaac 62120f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetAdjacencyDimension()`, `DMForestGetAdjacencyCodimension()`, `DMForestSetPartitionOverlap()` 6229be51f97SToby Isaac @*/ 623d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim) 624d71ae5a4SJacob Faibussowitsch { 625dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 626dd8e54a2SToby Isaac 627dd8e54a2SToby Isaac PetscFunctionBegin; 628dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 6294f572ea9SToby Isaac PetscAssertPointer(adjDim, 2); 630dd8e54a2SToby Isaac *adjDim = forest->adjDim; 6313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 632dd8e54a2SToby Isaac } 633dd8e54a2SToby Isaac 6349be51f97SToby Isaac /*@ 63520f4b53cSBarry Smith DMForestGetAdjacencyCodimension - Like `DMForestGetAdjacencyDimension()`, but specified as a co-dimension (so that, 6369be51f97SToby Isaac e.g., adjacency based on facets can be specified by codimension 1 in all cases) 6379be51f97SToby Isaac 63820f4b53cSBarry Smith Not Collective 6399be51f97SToby Isaac 6409be51f97SToby Isaac Input Parameter: 6419be51f97SToby Isaac . dm - the forest 6429be51f97SToby Isaac 6439be51f97SToby Isaac Output Parameter: 64420f4b53cSBarry Smith . adjCodim - default isthe dimension of the forest (see `DMGetDimension()`), since this is the codimension of vertices 6459be51f97SToby Isaac 6469be51f97SToby Isaac Level: intermediate 6479be51f97SToby Isaac 64820f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetAdjacencyCodimension()`, `DMForestGetAdjacencyDimension()` 6499be51f97SToby Isaac @*/ 650d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim) 651d71ae5a4SJacob Faibussowitsch { 652dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 653dd8e54a2SToby Isaac PetscInt dim; 654dd8e54a2SToby Isaac 655dd8e54a2SToby Isaac PetscFunctionBegin; 656dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 6574f572ea9SToby Isaac PetscAssertPointer(adjCodim, 2); 6589566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm, &dim)); 659dd8e54a2SToby Isaac *adjCodim = dim - forest->adjDim; 6603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 661dd8e54a2SToby Isaac } 662dd8e54a2SToby Isaac 6639be51f97SToby Isaac /*@ 6649be51f97SToby Isaac DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel 6659be51f97SToby Isaac partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding 6669be51f97SToby Isaac adjacent cells 6679be51f97SToby Isaac 66820f4b53cSBarry Smith Logically Collective 6699be51f97SToby Isaac 6709be51f97SToby Isaac Input Parameters: 6719be51f97SToby Isaac + dm - the forest 6729be51f97SToby Isaac - overlap - default 0 6739be51f97SToby Isaac 6749be51f97SToby Isaac Level: intermediate 6759be51f97SToby Isaac 67620f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetPartitionOverlap()`, `DMForestSetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()` 6779be51f97SToby Isaac @*/ 678d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap) 679d71ae5a4SJacob Faibussowitsch { 680dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 681dd8e54a2SToby Isaac 682dd8e54a2SToby Isaac PetscFunctionBegin; 683dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 68428b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the overlap after setup"); 68563a3b9bcSJacob Faibussowitsch PetscCheck(overlap >= 0, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "overlap cannot be < 0: %" PetscInt_FMT, overlap); 686dd8e54a2SToby Isaac forest->overlap = overlap; 6873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 688dd8e54a2SToby Isaac } 689dd8e54a2SToby Isaac 6909be51f97SToby Isaac /*@ 6919be51f97SToby Isaac DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values 6929be51f97SToby Isaac > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells 6939be51f97SToby Isaac 69420f4b53cSBarry Smith Not Collective 6959be51f97SToby Isaac 6969be51f97SToby Isaac Input Parameter: 6979be51f97SToby Isaac . dm - the forest 6989be51f97SToby Isaac 6999be51f97SToby Isaac Output Parameter: 7009be51f97SToby Isaac . overlap - default 0 7019be51f97SToby Isaac 7029be51f97SToby Isaac Level: intermediate 7039be51f97SToby Isaac 70442747ad1SJacob Faibussowitsch .seealso: `DM`, `DMFOREST`, `DMForestSetAdjacencyDimension()`, `DMForestSetAdjacencyCodimension()` 7059be51f97SToby Isaac @*/ 706d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap) 707d71ae5a4SJacob Faibussowitsch { 708dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 709dd8e54a2SToby Isaac 710dd8e54a2SToby Isaac PetscFunctionBegin; 711dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 7124f572ea9SToby Isaac PetscAssertPointer(overlap, 2); 713dd8e54a2SToby Isaac *overlap = forest->overlap; 7143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 715dd8e54a2SToby Isaac } 716dd8e54a2SToby Isaac 7179be51f97SToby Isaac /*@ 7189be51f97SToby Isaac DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base 71920f4b53cSBarry Smith `DM`, see `DMForestGetBaseDM()`) allowed in the forest. If the forest is being created by coarsening a previous forest 72020f4b53cSBarry Smith (see `DMForestGetAdaptivityForest()`) this limits the amount of coarsening. 7219be51f97SToby Isaac 72220f4b53cSBarry Smith Logically Collective 7239be51f97SToby Isaac 7249be51f97SToby Isaac Input Parameters: 7259be51f97SToby Isaac + dm - the forest 72620f4b53cSBarry Smith - minRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`) 7279be51f97SToby Isaac 7289be51f97SToby Isaac Level: intermediate 7299be51f97SToby Isaac 73020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestSetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()` 7319be51f97SToby Isaac @*/ 732d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement) 733d71ae5a4SJacob Faibussowitsch { 734dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 735dd8e54a2SToby Isaac 736dd8e54a2SToby Isaac PetscFunctionBegin; 737dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 73828b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the minimum refinement after setup"); 739dd8e54a2SToby Isaac forest->minRefinement = minRefinement; 7403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 741dd8e54a2SToby Isaac } 742dd8e54a2SToby Isaac 7439be51f97SToby Isaac /*@ 74420f4b53cSBarry Smith DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base `DM`, see 74520f4b53cSBarry Smith `DMForestGetBaseDM()`) allowed in the forest. If the forest is being created by coarsening a previous forest (see 74620f4b53cSBarry Smith `DMForestGetAdaptivityForest()`), this limits the amount of coarsening. 7479be51f97SToby Isaac 74820f4b53cSBarry Smith Not Collective 7499be51f97SToby Isaac 7509be51f97SToby Isaac Input Parameter: 7519be51f97SToby Isaac . dm - the forest 7529be51f97SToby Isaac 7539be51f97SToby Isaac Output Parameter: 75420f4b53cSBarry Smith . minRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`) 7559be51f97SToby Isaac 7569be51f97SToby Isaac Level: intermediate 7579be51f97SToby Isaac 75820f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetMinimumRefinement()`, `DMForestGetMaximumRefinement()`, `DMForestGetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()` 7599be51f97SToby Isaac @*/ 760d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement) 761d71ae5a4SJacob Faibussowitsch { 762dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 763dd8e54a2SToby Isaac 764dd8e54a2SToby Isaac PetscFunctionBegin; 765dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 7664f572ea9SToby Isaac PetscAssertPointer(minRefinement, 2); 767dd8e54a2SToby Isaac *minRefinement = forest->minRefinement; 7683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 769dd8e54a2SToby Isaac } 770dd8e54a2SToby Isaac 7719be51f97SToby Isaac /*@ 7729be51f97SToby Isaac DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base 77320f4b53cSBarry Smith `DM`, see `DMForestGetBaseDM()`) allowed in the forest. 7749be51f97SToby Isaac 77520f4b53cSBarry Smith Logically Collective 7769be51f97SToby Isaac 7779be51f97SToby Isaac Input Parameters: 7789be51f97SToby Isaac + dm - the forest 77960225df5SJacob Faibussowitsch - initRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`) 7809be51f97SToby Isaac 7819be51f97SToby Isaac Level: intermediate 7829be51f97SToby Isaac 78320f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestGetBaseDM()` 7849be51f97SToby Isaac @*/ 785d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement) 786d71ae5a4SJacob Faibussowitsch { 78756ba9f64SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 78856ba9f64SToby Isaac 78956ba9f64SToby Isaac PetscFunctionBegin; 79056ba9f64SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 79128b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the initial refinement after setup"); 79256ba9f64SToby Isaac forest->initRefinement = initRefinement; 7933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 79456ba9f64SToby Isaac } 79556ba9f64SToby Isaac 7969be51f97SToby Isaac /*@ 79720f4b53cSBarry Smith DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base `DM`, see 79820f4b53cSBarry Smith `DMForestGetBaseDM()`) allowed in the forest. 7999be51f97SToby Isaac 80020f4b53cSBarry Smith Not Collective 8019be51f97SToby Isaac 8029be51f97SToby Isaac Input Parameter: 8039be51f97SToby Isaac . dm - the forest 8049be51f97SToby Isaac 80501d2d390SJose E. Roman Output Parameter: 80620f4b53cSBarry Smith . initRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`) 8079be51f97SToby Isaac 8089be51f97SToby Isaac Level: intermediate 8099be51f97SToby Isaac 81020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetMinimumRefinement()`, `DMForestSetMaximumRefinement()`, `DMForestGetBaseDM()` 8119be51f97SToby Isaac @*/ 812d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement) 813d71ae5a4SJacob Faibussowitsch { 81456ba9f64SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 81556ba9f64SToby Isaac 81656ba9f64SToby Isaac PetscFunctionBegin; 81756ba9f64SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8184f572ea9SToby Isaac PetscAssertPointer(initRefinement, 2); 81956ba9f64SToby Isaac *initRefinement = forest->initRefinement; 8203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 82156ba9f64SToby Isaac } 82256ba9f64SToby Isaac 8239be51f97SToby Isaac /*@ 8249be51f97SToby Isaac DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base 82520f4b53cSBarry Smith `DM`, see `DMForestGetBaseDM()`) allowed in the forest. If the forest is being created by refining a previous forest 82620f4b53cSBarry Smith (see `DMForestGetAdaptivityForest()`), this limits the amount of refinement. 8279be51f97SToby Isaac 82820f4b53cSBarry Smith Logically Collective 8299be51f97SToby Isaac 8309be51f97SToby Isaac Input Parameters: 8319be51f97SToby Isaac + dm - the forest 83220f4b53cSBarry Smith - maxRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`) 8339be51f97SToby Isaac 8349be51f97SToby Isaac Level: intermediate 8359be51f97SToby Isaac 83642747ad1SJacob Faibussowitsch .seealso: `DM`, `DMFOREST`, `DMForestGetMinimumRefinement()`, `DMForestSetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityDM()` 8379be51f97SToby Isaac @*/ 838d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement) 839d71ae5a4SJacob Faibussowitsch { 840dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 841dd8e54a2SToby Isaac 842dd8e54a2SToby Isaac PetscFunctionBegin; 843dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 84428b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the maximum refinement after setup"); 845c7eeac06SToby Isaac forest->maxRefinement = maxRefinement; 8463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 847dd8e54a2SToby Isaac } 848dd8e54a2SToby Isaac 8499be51f97SToby Isaac /*@ 85020f4b53cSBarry Smith DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base `DM`, see 85120f4b53cSBarry Smith `DMForestGetBaseDM()`) allowed in the forest. If the forest is being created by refining a previous forest (see 852f826b5fcSPierre Jolivet `DMForestGetAdaptivityForest()`), this limits the amount of refinement. 8539be51f97SToby Isaac 85420f4b53cSBarry Smith Not Collective 8559be51f97SToby Isaac 8569be51f97SToby Isaac Input Parameter: 8579be51f97SToby Isaac . dm - the forest 8589be51f97SToby Isaac 8599be51f97SToby Isaac Output Parameter: 86020f4b53cSBarry Smith . maxRefinement - default `PETSC_DEFAULT` (interpreted by the subtype of `DMFOREST`) 8619be51f97SToby Isaac 8629be51f97SToby Isaac Level: intermediate 8639be51f97SToby Isaac 86420f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetMaximumRefinement()`, `DMForestGetMinimumRefinement()`, `DMForestGetInitialRefinement()`, `DMForestGetBaseDM()`, `DMForestGetAdaptivityForest()` 8659be51f97SToby Isaac @*/ 866d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement) 867d71ae5a4SJacob Faibussowitsch { 868dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 869dd8e54a2SToby Isaac 870dd8e54a2SToby Isaac PetscFunctionBegin; 871dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8724f572ea9SToby Isaac PetscAssertPointer(maxRefinement, 2); 873c7eeac06SToby Isaac *maxRefinement = forest->maxRefinement; 8743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 875dd8e54a2SToby Isaac } 876c7eeac06SToby Isaac 8779be51f97SToby Isaac /*@C 8789be51f97SToby Isaac DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes. 8799be51f97SToby Isaac 88020f4b53cSBarry Smith Logically Collective 8819be51f97SToby Isaac 8829be51f97SToby Isaac Input Parameters: 8839be51f97SToby Isaac + dm - the forest 88420f4b53cSBarry Smith - adaptStrategy - default `DMFORESTADAPTALL` 8859be51f97SToby Isaac 8869be51f97SToby Isaac Level: advanced 8879be51f97SToby Isaac 88820f4b53cSBarry Smith Notes: 88920f4b53cSBarry Smith Subtypes of `DMFOREST` may define their own strategies. Two default strategies are `DMFORESTADAPTALL`, which indicates that all processes must agree 89020f4b53cSBarry Smith for a refinement/coarsening flag to be valid, and `DMFORESTADAPTANY`, which indicates that only one process needs to 89120f4b53cSBarry Smith specify refinement/coarsening. 89220f4b53cSBarry Smith 89320f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetAdaptivityStrategy()`, `DMFORESTADAPTALL`, `DMFORESTADAPTANY` 8949be51f97SToby Isaac @*/ 895d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy) 896d71ae5a4SJacob Faibussowitsch { 897c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 898c7eeac06SToby Isaac 899c7eeac06SToby Isaac PetscFunctionBegin; 900c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 9019566063dSJacob Faibussowitsch PetscCall(PetscFree(forest->adaptStrategy)); 9029566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy((const char *)adaptStrategy, (char **)&forest->adaptStrategy)); 9033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 904c7eeac06SToby Isaac } 905c7eeac06SToby Isaac 9069be51f97SToby Isaac /*@C 90760225df5SJacob Faibussowitsch DMForestGetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes. 9089be51f97SToby Isaac 90920f4b53cSBarry Smith Not Collective 9109be51f97SToby Isaac 9119be51f97SToby Isaac Input Parameter: 9129be51f97SToby Isaac . dm - the forest 9139be51f97SToby Isaac 9149be51f97SToby Isaac Output Parameter: 91520f4b53cSBarry Smith . adaptStrategy - the adaptivity strategy (default `DMFORESTADAPTALL`) 9169be51f97SToby Isaac 9179be51f97SToby Isaac Level: advanced 9189be51f97SToby Isaac 91920f4b53cSBarry Smith Note: 92020f4b53cSBarry Smith Subtypes 92120f4b53cSBarry Smith of `DMFOREST` may define their own strategies. Two default strategies are `DMFORESTADAPTALL`, which indicates that all 92220f4b53cSBarry Smith processes must agree for a refinement/coarsening flag to be valid, and `DMFORESTADAPTANY`, which indicates that only 92320f4b53cSBarry Smith one process needs to specify refinement/coarsening. 92420f4b53cSBarry Smith 92520f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMFORESTADAPTALL`, `DMFORESTADAPTANY`, `DMForestSetAdaptivityStrategy()` 9269be51f97SToby Isaac @*/ 927d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy) 928d71ae5a4SJacob Faibussowitsch { 929c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 930c7eeac06SToby Isaac 931c7eeac06SToby Isaac PetscFunctionBegin; 932c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 9334f572ea9SToby Isaac PetscAssertPointer(adaptStrategy, 2); 934c7eeac06SToby Isaac *adaptStrategy = forest->adaptStrategy; 9353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 936c7eeac06SToby Isaac } 937c7eeac06SToby Isaac 9382a133e43SToby Isaac /*@ 9392a133e43SToby Isaac DMForestGetAdaptivitySuccess - Return whether the requested adaptation (refinement, coarsening, repartitioning, 94020f4b53cSBarry Smith etc.) was successful. 9412a133e43SToby Isaac 94220f4b53cSBarry Smith Collective 9432a133e43SToby Isaac 9442a133e43SToby Isaac Input Parameter: 9452a133e43SToby Isaac . dm - the post-adaptation forest 9462a133e43SToby Isaac 9472a133e43SToby Isaac Output Parameter: 94820f4b53cSBarry Smith . success - `PETSC_TRUE` if the post-adaptation forest is different from the pre-adaptation forest. 9492a133e43SToby Isaac 9502a133e43SToby Isaac Level: intermediate 9512a133e43SToby Isaac 95220f4b53cSBarry Smith Notes: 953d72c4dc2SStefano Zampini `PETSC_FALSE` indicates that the post-adaptation forest is the same as the pre-adaptation 95420f4b53cSBarry Smith forest. A requested adaptation may have been unsuccessful if, for example, the requested refinement would have 95520f4b53cSBarry Smith exceeded the maximum refinement level. 95620f4b53cSBarry Smith 95720f4b53cSBarry Smith .seealso: `DM`, `DMFOREST` 9582a133e43SToby Isaac @*/ 959d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success) 960d71ae5a4SJacob Faibussowitsch { 9612a133e43SToby Isaac DM_Forest *forest; 9622a133e43SToby Isaac 9632a133e43SToby Isaac PetscFunctionBegin; 9642a133e43SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 96528b400f6SJacob Faibussowitsch PetscCheck(dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "DMSetUp() has not been called yet."); 9662a133e43SToby Isaac forest = (DM_Forest *)dm->data; 967*f4f49eeaSPierre Jolivet PetscCall(forest->getadaptivitysuccess(dm, success)); 9683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9692a133e43SToby Isaac } 9702a133e43SToby Isaac 971bf9b5d84SToby Isaac /*@ 97220f4b53cSBarry Smith DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer `PetscSF`s should be computed 97320f4b53cSBarry Smith relating the cells of the pre-adaptation forest to the post-adaptiation forest. 974bf9b5d84SToby Isaac 97520f4b53cSBarry Smith Logically Collective 976bf9b5d84SToby Isaac 977bf9b5d84SToby Isaac Input Parameters: 978bf9b5d84SToby Isaac + dm - the post-adaptation forest 97920f4b53cSBarry Smith - computeSF - default `PETSC_TRUE` 980bf9b5d84SToby Isaac 981bf9b5d84SToby Isaac Level: advanced 982bf9b5d84SToby Isaac 98320f4b53cSBarry Smith Note: 98420f4b53cSBarry Smith After `DMSetUp()` is called, the transfer `PetscSF`s can be accessed with `DMForestGetAdaptivitySF()`. 98520f4b53cSBarry Smith 98620f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetComputeAdaptivitySF()`, `DMForestGetAdaptivitySF()` 987bf9b5d84SToby Isaac @*/ 988d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF) 989d71ae5a4SJacob Faibussowitsch { 990bf9b5d84SToby Isaac DM_Forest *forest; 991bf9b5d84SToby Isaac 992bf9b5d84SToby Isaac PetscFunctionBegin; 993bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 99428b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot compute adaptivity PetscSFs after setup is called"); 995bf9b5d84SToby Isaac forest = (DM_Forest *)dm->data; 996bf9b5d84SToby Isaac forest->computeAdaptSF = computeSF; 9973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 998bf9b5d84SToby Isaac } 999bf9b5d84SToby Isaac 1000d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time) 1001d71ae5a4SJacob Faibussowitsch { 100280b27e07SToby Isaac DM_Forest *forest; 100380b27e07SToby Isaac 100480b27e07SToby Isaac PetscFunctionBegin; 100580b27e07SToby Isaac PetscValidHeaderSpecific(dmIn, DM_CLASSID, 1); 100680b27e07SToby Isaac PetscValidHeaderSpecific(vecIn, VEC_CLASSID, 2); 100780b27e07SToby Isaac PetscValidHeaderSpecific(dmOut, DM_CLASSID, 3); 100880b27e07SToby Isaac PetscValidHeaderSpecific(vecOut, VEC_CLASSID, 4); 100980b27e07SToby Isaac forest = (DM_Forest *)dmIn->data; 101028b400f6SJacob Faibussowitsch PetscCheck(forest->transfervec, PetscObjectComm((PetscObject)dmIn), PETSC_ERR_SUP, "DMForestTransferVec() not implemented"); 1011*f4f49eeaSPierre Jolivet PetscCall(forest->transfervec(dmIn, vecIn, dmOut, vecOut, useBCs, time)); 10123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 101380b27e07SToby Isaac } 101480b27e07SToby Isaac 1015d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestTransferVecFromBase(DM dm, Vec vecIn, Vec vecOut) 1016d71ae5a4SJacob Faibussowitsch { 1017ac34a06fSStefano Zampini DM_Forest *forest; 1018ac34a06fSStefano Zampini 1019ac34a06fSStefano Zampini PetscFunctionBegin; 1020ac34a06fSStefano Zampini PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1021ac34a06fSStefano Zampini PetscValidHeaderSpecific(vecIn, VEC_CLASSID, 2); 1022ac34a06fSStefano Zampini PetscValidHeaderSpecific(vecOut, VEC_CLASSID, 3); 1023ac34a06fSStefano Zampini forest = (DM_Forest *)dm->data; 102428b400f6SJacob Faibussowitsch PetscCheck(forest->transfervecfrombase, PetscObjectComm((PetscObject)dm), PETSC_ERR_SUP, "DMForestTransferVecFromBase() not implemented"); 1025*f4f49eeaSPierre Jolivet PetscCall(forest->transfervecfrombase(dm, vecIn, vecOut)); 10263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1027ac34a06fSStefano Zampini } 1028ac34a06fSStefano Zampini 1029bf9b5d84SToby Isaac /*@ 103020f4b53cSBarry Smith DMForestGetComputeAdaptivitySF - Get whether transfer `PetscSF`s should be computed relating the cells of the 103120f4b53cSBarry Smith pre-adaptation forest to the post-adaptiation forest. After `DMSetUp()` is called, these transfer PetscSFs can be 103220f4b53cSBarry Smith accessed with `DMForestGetAdaptivitySF()`. 1033bf9b5d84SToby Isaac 103420f4b53cSBarry Smith Not Collective 1035bf9b5d84SToby Isaac 1036bf9b5d84SToby Isaac Input Parameter: 1037bf9b5d84SToby Isaac . dm - the post-adaptation forest 1038bf9b5d84SToby Isaac 1039bf9b5d84SToby Isaac Output Parameter: 104020f4b53cSBarry Smith . computeSF - default `PETSC_TRUE` 1041bf9b5d84SToby Isaac 1042bf9b5d84SToby Isaac Level: advanced 1043bf9b5d84SToby Isaac 104420f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetComputeAdaptivitySF()`, `DMForestGetAdaptivitySF()` 1045bf9b5d84SToby Isaac @*/ 1046d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF) 1047d71ae5a4SJacob Faibussowitsch { 1048bf9b5d84SToby Isaac DM_Forest *forest; 1049bf9b5d84SToby Isaac 1050bf9b5d84SToby Isaac PetscFunctionBegin; 1051bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1052bf9b5d84SToby Isaac forest = (DM_Forest *)dm->data; 1053bf9b5d84SToby Isaac *computeSF = forest->computeAdaptSF; 10543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1055bf9b5d84SToby Isaac } 1056bf9b5d84SToby Isaac 1057bf9b5d84SToby Isaac /*@ 1058a4e35b19SJacob Faibussowitsch DMForestGetAdaptivitySF - Get `PetscSF`s that relate the pre-adaptation forest to the 1059a4e35b19SJacob Faibussowitsch post-adaptation forest. 1060bf9b5d84SToby Isaac 106120f4b53cSBarry Smith Not Collective 1062bf9b5d84SToby Isaac 1063bf9b5d84SToby Isaac Input Parameter: 106420f4b53cSBarry Smith . dm - the post-adaptation forest 1065bf9b5d84SToby Isaac 106620f4b53cSBarry Smith Output Parameters: 106720f4b53cSBarry Smith + preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post- 106820f4b53cSBarry Smith - coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre- 1069bf9b5d84SToby Isaac 1070bf9b5d84SToby Isaac Level: advanced 1071bf9b5d84SToby Isaac 1072a4e35b19SJacob Faibussowitsch Notes: 1073a4e35b19SJacob Faibussowitsch Adaptation can be any combination of refinement, coarsening, repartition, and change of 1074a4e35b19SJacob Faibussowitsch overlap, so there may be some cells of the pre-adaptation that are parents of post-adaptation 1075a4e35b19SJacob Faibussowitsch cells, and vice versa. Therefore there are two `PetscSF`s: one that relates pre-adaptation 1076a4e35b19SJacob Faibussowitsch coarse cells to post-adaptation fine cells, and one that relates pre-adaptation fine cells to 1077a4e35b19SJacob Faibussowitsch post-adaptation coarse cells. 1078a4e35b19SJacob Faibussowitsch 107920f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetComputeAdaptivitySF()`, `DMForestSetComputeAdaptivitySF()` 1080bf9b5d84SToby Isaac @*/ 1081d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine) 1082d71ae5a4SJacob Faibussowitsch { 1083bf9b5d84SToby Isaac DM_Forest *forest; 1084bf9b5d84SToby Isaac 1085bf9b5d84SToby Isaac PetscFunctionBegin; 1086bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 10879566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 1088bf9b5d84SToby Isaac forest = (DM_Forest *)dm->data; 1089f885a11aSToby Isaac if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine; 1090f885a11aSToby Isaac if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine; 10913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1092bf9b5d84SToby Isaac } 1093bf9b5d84SToby Isaac 10949be51f97SToby Isaac /*@ 1095a4e35b19SJacob Faibussowitsch DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the 1096a4e35b19SJacob Faibussowitsch mesh, e.g. give 2 to indicate that the diameter of neighboring cells should differ by at most 1097a4e35b19SJacob Faibussowitsch a factor of 2. 10989be51f97SToby Isaac 109920f4b53cSBarry Smith Logically Collective 11009be51f97SToby Isaac 11019be51f97SToby Isaac Input Parameters: 11029be51f97SToby Isaac + dm - the forest 11039be51f97SToby Isaac - grade - the grading factor 11049be51f97SToby Isaac 11059be51f97SToby Isaac Level: advanced 11069be51f97SToby Isaac 1107a4e35b19SJacob Faibussowitsch Note: 1108a4e35b19SJacob Faibussowitsch Subtypes of `DMFOREST` may only support one particular choice of grading factor. 1109a4e35b19SJacob Faibussowitsch 111020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetGradeFactor()` 11119be51f97SToby Isaac @*/ 1112d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade) 1113d71ae5a4SJacob Faibussowitsch { 1114c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 1115c7eeac06SToby Isaac 1116c7eeac06SToby Isaac PetscFunctionBegin; 1117c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 111828b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the grade factor after setup"); 1119c7eeac06SToby Isaac forest->gradeFactor = grade; 11203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1121c7eeac06SToby Isaac } 1122c7eeac06SToby Isaac 11239be51f97SToby Isaac /*@ 11249be51f97SToby Isaac DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of 112520f4b53cSBarry Smith neighboring cells should differ by at most a factor of 2. Subtypes of `DMFOREST` may only support one particular 11269be51f97SToby Isaac choice of grading factor. 11279be51f97SToby Isaac 112820f4b53cSBarry Smith Not Collective 11299be51f97SToby Isaac 11309be51f97SToby Isaac Input Parameter: 11319be51f97SToby Isaac . dm - the forest 11329be51f97SToby Isaac 11339be51f97SToby Isaac Output Parameter: 11349be51f97SToby Isaac . grade - the grading factor 11359be51f97SToby Isaac 11369be51f97SToby Isaac Level: advanced 11379be51f97SToby Isaac 113820f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetGradeFactor()` 11399be51f97SToby Isaac @*/ 1140d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade) 1141d71ae5a4SJacob Faibussowitsch { 1142c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 1143c7eeac06SToby Isaac 1144c7eeac06SToby Isaac PetscFunctionBegin; 1145c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 11464f572ea9SToby Isaac PetscAssertPointer(grade, 2); 1147c7eeac06SToby Isaac *grade = forest->gradeFactor; 11483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1149c7eeac06SToby Isaac } 1150c7eeac06SToby Isaac 11519be51f97SToby Isaac /*@ 11529be51f97SToby Isaac DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes 1153a4e35b19SJacob Faibussowitsch the cell weight (see `DMForestSetCellWeights()`) when calculating partitions. 11549be51f97SToby Isaac 115520f4b53cSBarry Smith Logically Collective 11569be51f97SToby Isaac 11579be51f97SToby Isaac Input Parameters: 11589be51f97SToby Isaac + dm - the forest 115960225df5SJacob Faibussowitsch - weightsFactor - default 1. 11609be51f97SToby Isaac 11619be51f97SToby Isaac Level: advanced 11629be51f97SToby Isaac 1163a4e35b19SJacob Faibussowitsch Note: 1164a4e35b19SJacob Faibussowitsch The final weight of a cell will be (cellWeight) * (weightFactor^refinementLevel). A factor 1165a4e35b19SJacob Faibussowitsch of 1 indicates that the weight of a cell does not depend on its level; a factor of 2, for 1166a4e35b19SJacob Faibussowitsch example, might be appropriate for sub-cycling time-stepping methods, when the computation 1167a4e35b19SJacob Faibussowitsch associated with a cell is multiplied by a factor of 2 for each additional level of 1168a4e35b19SJacob Faibussowitsch refinement. 1169a4e35b19SJacob Faibussowitsch 117020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetCellWeightFactor()`, `DMForestSetCellWeights()` 11719be51f97SToby Isaac @*/ 1172d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor) 1173d71ae5a4SJacob Faibussowitsch { 1174c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 1175c7eeac06SToby Isaac 1176c7eeac06SToby Isaac PetscFunctionBegin; 1177c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 117828b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the weights factor after setup"); 1179c7eeac06SToby Isaac forest->weightsFactor = weightsFactor; 11803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1181c7eeac06SToby Isaac } 1182c7eeac06SToby Isaac 11839be51f97SToby Isaac /*@ 11849be51f97SToby Isaac DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see 1185a4e35b19SJacob Faibussowitsch `DMForestSetCellWeights()`) when calculating partitions. 11869be51f97SToby Isaac 118720f4b53cSBarry Smith Not Collective 11889be51f97SToby Isaac 11899be51f97SToby Isaac Input Parameter: 11909be51f97SToby Isaac . dm - the forest 11919be51f97SToby Isaac 11929be51f97SToby Isaac Output Parameter: 119360225df5SJacob Faibussowitsch . weightsFactor - default 1. 11949be51f97SToby Isaac 11959be51f97SToby Isaac Level: advanced 11969be51f97SToby Isaac 1197a4e35b19SJacob Faibussowitsch Note: 1198a4e35b19SJacob Faibussowitsch The final weight of a cell will be (cellWeight) * (weightFactor^refinementLevel). A factor 1199a4e35b19SJacob Faibussowitsch of 1 indicates that the weight of a cell does not depend on its level; a factor of 2, for 1200a4e35b19SJacob Faibussowitsch example, might be appropriate for sub-cycling time-stepping methods, when the computation 1201a4e35b19SJacob Faibussowitsch associated with a cell is multiplied by a factor of 2 for each additional level of 1202a4e35b19SJacob Faibussowitsch refinement. 1203a4e35b19SJacob Faibussowitsch 120420f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetCellWeightFactor()`, `DMForestSetCellWeights()` 12059be51f97SToby Isaac @*/ 1206d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor) 1207d71ae5a4SJacob Faibussowitsch { 1208c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 1209c7eeac06SToby Isaac 1210c7eeac06SToby Isaac PetscFunctionBegin; 1211c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12124f572ea9SToby Isaac PetscAssertPointer(weightsFactor, 2); 1213c7eeac06SToby Isaac *weightsFactor = forest->weightsFactor; 12143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1215c7eeac06SToby Isaac } 1216c7eeac06SToby Isaac 12179be51f97SToby Isaac /*@ 12189be51f97SToby Isaac DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process 12199be51f97SToby Isaac 122020f4b53cSBarry Smith Not Collective 12219be51f97SToby Isaac 12229be51f97SToby Isaac Input Parameter: 12239be51f97SToby Isaac . dm - the forest 12249be51f97SToby Isaac 12259be51f97SToby Isaac Output Parameters: 12269be51f97SToby Isaac + cStart - the first cell on this process 12279be51f97SToby Isaac - cEnd - one after the final cell on this process 12289be51f97SToby Isaac 12291a244344SSatish Balay Level: intermediate 12309be51f97SToby Isaac 123120f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetCellSF()` 12329be51f97SToby Isaac @*/ 1233d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd) 1234d71ae5a4SJacob Faibussowitsch { 1235c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 1236c7eeac06SToby Isaac 1237c7eeac06SToby Isaac PetscFunctionBegin; 1238c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12394f572ea9SToby Isaac PetscAssertPointer(cStart, 2); 12404f572ea9SToby Isaac PetscAssertPointer(cEnd, 3); 124148a46eb9SPierre Jolivet if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) PetscCall(forest->createcellchart(dm, &forest->cStart, &forest->cEnd)); 1242c7eeac06SToby Isaac *cStart = forest->cStart; 1243c7eeac06SToby Isaac *cEnd = forest->cEnd; 12443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1245c7eeac06SToby Isaac } 1246c7eeac06SToby Isaac 12479be51f97SToby Isaac /*@ 124820f4b53cSBarry Smith DMForestGetCellSF - After the setup phase, get the `PetscSF` for overlapping cells between processes 12499be51f97SToby Isaac 125020f4b53cSBarry Smith Not Collective 12519be51f97SToby Isaac 12529be51f97SToby Isaac Input Parameter: 12539be51f97SToby Isaac . dm - the forest 12549be51f97SToby Isaac 12559be51f97SToby Isaac Output Parameter: 125620f4b53cSBarry Smith . cellSF - the `PetscSF` 12579be51f97SToby Isaac 12581a244344SSatish Balay Level: intermediate 12599be51f97SToby Isaac 126020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetCellChart()` 12619be51f97SToby Isaac @*/ 1262d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF) 1263d71ae5a4SJacob Faibussowitsch { 1264c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 1265c7eeac06SToby Isaac 1266c7eeac06SToby Isaac PetscFunctionBegin; 1267c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12684f572ea9SToby Isaac PetscAssertPointer(cellSF, 2); 126948a46eb9SPierre Jolivet if ((!forest->cellSF) && forest->createcellsf) PetscCall(forest->createcellsf(dm, &forest->cellSF)); 1270c7eeac06SToby Isaac *cellSF = forest->cellSF; 12713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1272c7eeac06SToby Isaac } 1273c7eeac06SToby Isaac 12749be51f97SToby Isaac /*@C 12759be51f97SToby Isaac DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see 1276a4e35b19SJacob Faibussowitsch `DMForestGetAdaptivityForest()`) that holds the adaptation flags (refinement, coarsening, or some combination). 12779be51f97SToby Isaac 127820f4b53cSBarry Smith Logically Collective 12799be51f97SToby Isaac 12809be51f97SToby Isaac Input Parameters: 128160225df5SJacob Faibussowitsch + dm - the forest 128260225df5SJacob Faibussowitsch - adaptLabel - the label in the pre-adaptation forest 12839be51f97SToby Isaac 12849be51f97SToby Isaac Level: intermediate 12859be51f97SToby Isaac 1286a4e35b19SJacob Faibussowitsch Note: 1287a4e35b19SJacob Faibussowitsch The interpretation of the label values is up to the subtype of `DMFOREST`, but 1288a4e35b19SJacob Faibussowitsch `DM_ADAPT_DETERMINE`, `DM_ADAPT_KEEP`, `DM_ADAPT_REFINE`, and `DM_ADAPT_COARSEN` have been 1289a4e35b19SJacob Faibussowitsch reserved as choices that should be accepted by all subtypes. 1290a4e35b19SJacob Faibussowitsch 129120f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetAdaptivityLabel()` 12929be51f97SToby Isaac @*/ 1293d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetAdaptivityLabel(DM dm, DMLabel adaptLabel) 1294d71ae5a4SJacob Faibussowitsch { 1295c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 1296c7eeac06SToby Isaac 1297c7eeac06SToby Isaac PetscFunctionBegin; 1298c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 12991fd50544SStefano Zampini if (adaptLabel) PetscValidHeaderSpecific(adaptLabel, DMLABEL_CLASSID, 2); 13009566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)adaptLabel)); 13019566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&forest->adaptLabel)); 13021fd50544SStefano Zampini forest->adaptLabel = adaptLabel; 13033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1304c7eeac06SToby Isaac } 1305c7eeac06SToby Isaac 13069be51f97SToby Isaac /*@C 130720f4b53cSBarry Smith DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see `DMForestGetAdaptivityForest()`) that 1308a4e35b19SJacob Faibussowitsch holds the adaptation flags (refinement, coarsening, or some combination). 13099be51f97SToby Isaac 131020f4b53cSBarry Smith Not Collective 13119be51f97SToby Isaac 13129be51f97SToby Isaac Input Parameter: 13139be51f97SToby Isaac . dm - the forest 13149be51f97SToby Isaac 13159be51f97SToby Isaac Output Parameter: 13169be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest 13179be51f97SToby Isaac 13189be51f97SToby Isaac Level: intermediate 13199be51f97SToby Isaac 1320a4e35b19SJacob Faibussowitsch Note: 1321a4e35b19SJacob Faibussowitsch The interpretation of the label values is up to the subtype of `DMFOREST`, but 1322a4e35b19SJacob Faibussowitsch `DM_ADAPT_DETERMINE`, `DM_ADAPT_KEEP`, `DM_ADAPT_REFINE`, and `DM_ADAPT_COARSEN` have been 1323a4e35b19SJacob Faibussowitsch reserved as choices that should be accepted by all subtypes. 1324a4e35b19SJacob Faibussowitsch 132520f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetAdaptivityLabel()` 13269be51f97SToby Isaac @*/ 1327d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetAdaptivityLabel(DM dm, DMLabel *adaptLabel) 1328d71ae5a4SJacob Faibussowitsch { 1329c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 1330c7eeac06SToby Isaac 1331c7eeac06SToby Isaac PetscFunctionBegin; 1332c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1333ba936b91SToby Isaac *adaptLabel = forest->adaptLabel; 13343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1335c7eeac06SToby Isaac } 1336c7eeac06SToby Isaac 13379be51f97SToby Isaac /*@ 133820f4b53cSBarry Smith DMForestSetCellWeights - Set the weights assigned to each of the cells (see `DMForestGetCellChart()`) of the current 133920f4b53cSBarry Smith process: weights are used to determine parallel partitioning. 13409be51f97SToby Isaac 134120f4b53cSBarry Smith Logically Collective 13429be51f97SToby Isaac 13439be51f97SToby Isaac Input Parameters: 13449be51f97SToby Isaac + dm - the forest 134520f4b53cSBarry Smith . weights - the array of weights (see `DMForestSetWeightCapacity()`) for all cells, or `NULL` to indicate each cell has weight 1. 13469be51f97SToby Isaac - copyMode - how weights should reference weights 13479be51f97SToby Isaac 13489be51f97SToby Isaac Level: advanced 13499be51f97SToby Isaac 135020f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestGetCellWeights()`, `DMForestSetWeightCapacity()` 13519be51f97SToby Isaac @*/ 1352d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode) 1353d71ae5a4SJacob Faibussowitsch { 1354c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 1355c7eeac06SToby Isaac PetscInt cStart, cEnd; 1356c7eeac06SToby Isaac 1357c7eeac06SToby Isaac PetscFunctionBegin; 1358c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 13599566063dSJacob Faibussowitsch PetscCall(DMForestGetCellChart(dm, &cStart, &cEnd)); 136063a3b9bcSJacob Faibussowitsch PetscCheck(cEnd >= cStart, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "cell chart [%" PetscInt_FMT ",%" PetscInt_FMT ") is not valid", cStart, cEnd); 1361c7eeac06SToby Isaac if (copyMode == PETSC_COPY_VALUES) { 136248a46eb9SPierre Jolivet if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) PetscCall(PetscMalloc1(cEnd - cStart, &forest->cellWeights)); 13639566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(forest->cellWeights, weights, cEnd - cStart)); 1364c7eeac06SToby Isaac forest->cellWeightsCopyMode = PETSC_OWN_POINTER; 13653ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1366c7eeac06SToby Isaac } 136748a46eb9SPierre Jolivet if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) PetscCall(PetscFree(forest->cellWeights)); 1368c7eeac06SToby Isaac forest->cellWeights = weights; 1369c7eeac06SToby Isaac forest->cellWeightsCopyMode = copyMode; 13703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1371c7eeac06SToby Isaac } 1372c7eeac06SToby Isaac 13739be51f97SToby Isaac /*@ 137420f4b53cSBarry Smith DMForestGetCellWeights - Get the weights assigned to each of the cells (see `DMForestGetCellChart()`) of the current 137520f4b53cSBarry Smith process: weights are used to determine parallel partitioning. 13769be51f97SToby Isaac 137720f4b53cSBarry Smith Not Collective 13789be51f97SToby Isaac 13799be51f97SToby Isaac Input Parameter: 13809be51f97SToby Isaac . dm - the forest 13819be51f97SToby Isaac 13829be51f97SToby Isaac Output Parameter: 138320f4b53cSBarry Smith . weights - the array of weights for all cells, or `NULL` to indicate each cell has weight 1. 13849be51f97SToby Isaac 13859be51f97SToby Isaac Level: advanced 13869be51f97SToby Isaac 138720f4b53cSBarry Smith .seealso: `DM`, `DMFOREST`, `DMForestSetCellWeights()`, `DMForestSetWeightCapacity()` 13889be51f97SToby Isaac @*/ 1389d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights) 1390d71ae5a4SJacob Faibussowitsch { 1391c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 1392c7eeac06SToby Isaac 1393c7eeac06SToby Isaac PetscFunctionBegin; 1394c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 13954f572ea9SToby Isaac PetscAssertPointer(weights, 2); 1396c7eeac06SToby Isaac *weights = forest->cellWeights; 13973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1398c7eeac06SToby Isaac } 1399c7eeac06SToby Isaac 14009be51f97SToby Isaac /*@ 14019be51f97SToby Isaac DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning 1402a4e35b19SJacob Faibussowitsch a pre-adaptation forest (see `DMForestGetAdaptivityForest()`). 14039be51f97SToby Isaac 140420f4b53cSBarry Smith Logically Collective 14059be51f97SToby Isaac 140660225df5SJacob Faibussowitsch Input Parameters: 14079be51f97SToby Isaac + dm - the forest 14089be51f97SToby Isaac - capacity - this process's capacity 14099be51f97SToby Isaac 14109be51f97SToby Isaac Level: advanced 14119be51f97SToby Isaac 1412a4e35b19SJacob Faibussowitsch Note: 1413a4e35b19SJacob Faibussowitsch After partitioning, the ratio of the weight of each process's cells to the process's capacity 1414a4e35b19SJacob Faibussowitsch will be roughly equal for all processes. A capacity of 0 indicates that the current process 1415a4e35b19SJacob Faibussowitsch should not have any cells after repartitioning. 1416a4e35b19SJacob Faibussowitsch 1417a4e35b19SJacob Faibussowitsch .seealso: `DM`, `DMFOREST`, `DMForestGetWeightCapacity()`, `DMForestSetCellWeights()`, `DMForestSetCellWeightFactor()` 14189be51f97SToby Isaac @*/ 1419d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity) 1420d71ae5a4SJacob Faibussowitsch { 1421c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 1422c7eeac06SToby Isaac 1423c7eeac06SToby Isaac PetscFunctionBegin; 1424c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 142528b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE, "Cannot change the weight capacity after setup"); 142663a3b9bcSJacob Faibussowitsch PetscCheck(capacity >= 0., PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_OUTOFRANGE, "Cannot have negative weight capacity; %g", (double)capacity); 1427c7eeac06SToby Isaac forest->weightCapacity = capacity; 14283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1429c7eeac06SToby Isaac } 1430c7eeac06SToby Isaac 14319be51f97SToby Isaac /*@ 14329be51f97SToby Isaac DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see 1433a4e35b19SJacob Faibussowitsch `DMForestGetAdaptivityForest()`). 14349be51f97SToby Isaac 143520f4b53cSBarry Smith Not Collective 14369be51f97SToby Isaac 143760225df5SJacob Faibussowitsch Input Parameter: 14389be51f97SToby Isaac . dm - the forest 14399be51f97SToby Isaac 144060225df5SJacob Faibussowitsch Output Parameter: 14419be51f97SToby Isaac . capacity - this process's capacity 14429be51f97SToby Isaac 14439be51f97SToby Isaac Level: advanced 14449be51f97SToby Isaac 1445a4e35b19SJacob Faibussowitsch Note: 1446a4e35b19SJacob Faibussowitsch After partitioning, the ratio of the weight of each process's cells to the process's capacity 1447a4e35b19SJacob Faibussowitsch will be roughly equal for all processes. A capacity of 0 indicates that the current process 1448a4e35b19SJacob Faibussowitsch should not have any cells after repartitioning. 1449a4e35b19SJacob Faibussowitsch 1450a4e35b19SJacob Faibussowitsch .seealso: `DM`, `DMFOREST`, `DMForestSetWeightCapacity()`, `DMForestSetCellWeights()`, `DMForestSetCellWeightFactor()` 14519be51f97SToby Isaac @*/ 1452d71ae5a4SJacob Faibussowitsch PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity) 1453d71ae5a4SJacob Faibussowitsch { 1454c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *)dm->data; 1455c7eeac06SToby Isaac 1456c7eeac06SToby Isaac PetscFunctionBegin; 1457c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 14584f572ea9SToby Isaac PetscAssertPointer(capacity, 2); 1459c7eeac06SToby Isaac *capacity = forest->weightCapacity; 14603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1461c7eeac06SToby Isaac } 1462c7eeac06SToby Isaac 1463d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(DM dm, PetscOptionItems *PetscOptionsObject) 1464d71ae5a4SJacob Faibussowitsch { 146556ba9f64SToby Isaac PetscBool flg, flg1, flg2, flg3, flg4; 1466dd8e54a2SToby Isaac DMForestTopology oldTopo; 1467c7eeac06SToby Isaac char stringBuffer[256]; 1468dd8e54a2SToby Isaac PetscViewer viewer; 1469dd8e54a2SToby Isaac PetscViewerFormat format; 147056ba9f64SToby Isaac PetscInt adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade; 1471c7eeac06SToby Isaac PetscReal weightsFactor; 1472c7eeac06SToby Isaac DMForestAdaptivityStrategy adaptStrategy; 1473db4d5e8cSToby Isaac 1474db4d5e8cSToby Isaac PetscFunctionBegin; 14759566063dSJacob Faibussowitsch PetscCall(DMForestGetTopology(dm, &oldTopo)); 1476d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "DMForest Options"); 14779566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-dm_forest_topology", "the topology of the forest's base mesh", "DMForestSetTopology", oldTopo, stringBuffer, sizeof(stringBuffer), &flg1)); 14789566063dSJacob Faibussowitsch PetscCall(PetscOptionsViewer("-dm_forest_base_dm", "load the base DM from a viewer specification", "DMForestSetBaseDM", &viewer, &format, &flg2)); 14799566063dSJacob Faibussowitsch PetscCall(PetscOptionsViewer("-dm_forest_coarse_forest", "load the coarse forest from a viewer specification", "DMForestSetCoarseForest", &viewer, &format, &flg3)); 14809566063dSJacob Faibussowitsch PetscCall(PetscOptionsViewer("-dm_forest_fine_forest", "load the fine forest from a viewer specification", "DMForestSetFineForest", &viewer, &format, &flg4)); 14811dca8a05SBarry Smith PetscCheck((PetscInt)flg1 + (PetscInt)flg2 + (PetscInt)flg3 + (PetscInt)flg4 <= 1, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_INCOMP, "Specify only one of -dm_forest_{topology,base_dm,coarse_forest,fine_forest}"); 148256ba9f64SToby Isaac if (flg1) { 14839566063dSJacob Faibussowitsch PetscCall(DMForestSetTopology(dm, (DMForestTopology)stringBuffer)); 14849566063dSJacob Faibussowitsch PetscCall(DMForestSetBaseDM(dm, NULL)); 14859566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityForest(dm, NULL)); 148656ba9f64SToby Isaac } 148756ba9f64SToby Isaac if (flg2) { 1488dd8e54a2SToby Isaac DM base; 1489dd8e54a2SToby Isaac 14909566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &base)); 14919566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, format)); 14929566063dSJacob Faibussowitsch PetscCall(DMLoad(base, viewer)); 14939566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 14949566063dSJacob Faibussowitsch PetscCall(DMForestSetBaseDM(dm, base)); 14959566063dSJacob Faibussowitsch PetscCall(DMDestroy(&base)); 14969566063dSJacob Faibussowitsch PetscCall(DMForestSetTopology(dm, NULL)); 14979566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityForest(dm, NULL)); 1498dd8e54a2SToby Isaac } 149956ba9f64SToby Isaac if (flg3) { 1500dd8e54a2SToby Isaac DM coarse; 1501dd8e54a2SToby Isaac 15029566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &coarse)); 15039566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, format)); 15049566063dSJacob Faibussowitsch PetscCall(DMLoad(coarse, viewer)); 15059566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 15069566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityForest(dm, coarse)); 15079566063dSJacob Faibussowitsch PetscCall(DMDestroy(&coarse)); 15089566063dSJacob Faibussowitsch PetscCall(DMForestSetTopology(dm, NULL)); 15099566063dSJacob Faibussowitsch PetscCall(DMForestSetBaseDM(dm, NULL)); 1510dd8e54a2SToby Isaac } 151156ba9f64SToby Isaac if (flg4) { 1512dd8e54a2SToby Isaac DM fine; 1513dd8e54a2SToby Isaac 15149566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm), &fine)); 15159566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, format)); 15169566063dSJacob Faibussowitsch PetscCall(DMLoad(fine, viewer)); 15179566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 15189566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityForest(dm, fine)); 15199566063dSJacob Faibussowitsch PetscCall(DMDestroy(&fine)); 15209566063dSJacob Faibussowitsch PetscCall(DMForestSetTopology(dm, NULL)); 15219566063dSJacob Faibussowitsch PetscCall(DMForestSetBaseDM(dm, NULL)); 1522dd8e54a2SToby Isaac } 15239566063dSJacob Faibussowitsch PetscCall(DMForestGetAdjacencyDimension(dm, &adjDim)); 15249566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_forest_adjacency_dimension", "set the dimension of points that define adjacency in the forest", "DMForestSetAdjacencyDimension", adjDim, &adjDim, &flg, 0)); 1525dd8e54a2SToby Isaac if (flg) { 15269566063dSJacob Faibussowitsch PetscCall(DMForestSetAdjacencyDimension(dm, adjDim)); 1527f885a11aSToby Isaac } else { 15289566063dSJacob Faibussowitsch PetscCall(DMForestGetAdjacencyCodimension(dm, &adjCodim)); 15299566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_forest_adjacency_codimension", "set the codimension of points that define adjacency in the forest", "DMForestSetAdjacencyCodimension", adjCodim, &adjCodim, &flg, 1)); 15301baa6e33SBarry Smith if (flg) PetscCall(DMForestSetAdjacencyCodimension(dm, adjCodim)); 1531dd8e54a2SToby Isaac } 15329566063dSJacob Faibussowitsch PetscCall(DMForestGetPartitionOverlap(dm, &overlap)); 15339566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_forest_partition_overlap", "set the degree of partition overlap", "DMForestSetPartitionOverlap", overlap, &overlap, &flg, 0)); 15341baa6e33SBarry Smith if (flg) PetscCall(DMForestSetPartitionOverlap(dm, overlap)); 1535a6121fbdSMatthew G. Knepley #if 0 15369566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_refine","equivalent to -dm_forest_set_minimum_refinement and -dm_forest_set_initial_refinement with the same value",NULL,minRefinement,&minRefinement,&flg,0)); 1537a6121fbdSMatthew G. Knepley if (flg) { 15389566063dSJacob Faibussowitsch PetscCall(DMForestSetMinimumRefinement(dm,minRefinement)); 15399566063dSJacob Faibussowitsch PetscCall(DMForestSetInitialRefinement(dm,minRefinement)); 1540a6121fbdSMatthew G. Knepley } 15419566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_refine_hierarchy","equivalent to -dm_forest_set_minimum_refinement 0 and -dm_forest_set_initial_refinement",NULL,initRefinement,&initRefinement,&flg,0)); 1542a6121fbdSMatthew G. Knepley if (flg) { 15439566063dSJacob Faibussowitsch PetscCall(DMForestSetMinimumRefinement(dm,0)); 15449566063dSJacob Faibussowitsch PetscCall(DMForestSetInitialRefinement(dm,initRefinement)); 1545a6121fbdSMatthew G. Knepley } 1546a6121fbdSMatthew G. Knepley #endif 15479566063dSJacob Faibussowitsch PetscCall(DMForestGetMinimumRefinement(dm, &minRefinement)); 15489566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_forest_minimum_refinement", "set the minimum level of refinement in the forest", "DMForestSetMinimumRefinement", minRefinement, &minRefinement, &flg, 0)); 15491baa6e33SBarry Smith if (flg) PetscCall(DMForestSetMinimumRefinement(dm, minRefinement)); 15509566063dSJacob Faibussowitsch PetscCall(DMForestGetInitialRefinement(dm, &initRefinement)); 15519566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_forest_initial_refinement", "set the initial level of refinement in the forest", "DMForestSetInitialRefinement", initRefinement, &initRefinement, &flg, 0)); 15521baa6e33SBarry Smith if (flg) PetscCall(DMForestSetInitialRefinement(dm, initRefinement)); 15539566063dSJacob Faibussowitsch PetscCall(DMForestGetMaximumRefinement(dm, &maxRefinement)); 15549566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_forest_maximum_refinement", "set the maximum level of refinement in the forest", "DMForestSetMaximumRefinement", maxRefinement, &maxRefinement, &flg, 0)); 15551baa6e33SBarry Smith if (flg) PetscCall(DMForestSetMaximumRefinement(dm, maxRefinement)); 15569566063dSJacob Faibussowitsch PetscCall(DMForestGetAdaptivityStrategy(dm, &adaptStrategy)); 15579566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-dm_forest_adaptivity_strategy", "the forest's adaptivity-flag resolution strategy", "DMForestSetAdaptivityStrategy", adaptStrategy, stringBuffer, sizeof(stringBuffer), &flg)); 15581baa6e33SBarry Smith if (flg) PetscCall(DMForestSetAdaptivityStrategy(dm, (DMForestAdaptivityStrategy)stringBuffer)); 15599566063dSJacob Faibussowitsch PetscCall(DMForestGetGradeFactor(dm, &grade)); 15609566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_forest_grade_factor", "grade factor between neighboring cells", "DMForestSetGradeFactor", grade, &grade, &flg, 0)); 15611baa6e33SBarry Smith if (flg) PetscCall(DMForestSetGradeFactor(dm, grade)); 15629566063dSJacob Faibussowitsch PetscCall(DMForestGetCellWeightFactor(dm, &weightsFactor)); 15639566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-dm_forest_cell_weight_factor", "multiplying weight factor for cell refinement", "DMForestSetCellWeightFactor", weightsFactor, &weightsFactor, &flg)); 15641baa6e33SBarry Smith if (flg) PetscCall(DMForestSetCellWeightFactor(dm, weightsFactor)); 1565d0609cedSBarry Smith PetscOptionsHeadEnd(); 15663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1567db4d5e8cSToby Isaac } 1568db4d5e8cSToby Isaac 156966976f2fSJacob Faibussowitsch static PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 1570d71ae5a4SJacob Faibussowitsch { 1571d8984e3bSMatthew G. Knepley PetscFunctionBegin; 15729566063dSJacob Faibussowitsch if (subdm) PetscCall(DMClone(dm, subdm)); 1573dd072f5fSMatthew G. Knepley PetscCall(DMCreateSectionSubDM(dm, numFields, fields, NULL, NULL, is, subdm)); 15743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1575d8984e3bSMatthew G. Knepley } 1576d8984e3bSMatthew G. Knepley 157766976f2fSJacob Faibussowitsch static PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined) 1578d71ae5a4SJacob Faibussowitsch { 15795421bac9SToby Isaac DMLabel refine; 15805421bac9SToby Isaac DM fineDM; 15815421bac9SToby Isaac 15825421bac9SToby Isaac PetscFunctionBegin; 15839566063dSJacob Faibussowitsch PetscCall(DMGetFineDM(dm, &fineDM)); 15845421bac9SToby Isaac if (fineDM) { 15859566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fineDM)); 15865421bac9SToby Isaac *dmRefined = fineDM; 15873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15885421bac9SToby Isaac } 15899566063dSJacob Faibussowitsch PetscCall(DMForestTemplate(dm, comm, dmRefined)); 15909566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "refine", &refine)); 15915421bac9SToby Isaac if (!refine) { 15929566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "refine", &refine)); 15939566063dSJacob Faibussowitsch PetscCall(DMLabelSetDefaultValue(refine, DM_ADAPT_REFINE)); 15941baa6e33SBarry Smith } else PetscCall(PetscObjectReference((PetscObject)refine)); 15959566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityLabel(*dmRefined, refine)); 15969566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&refine)); 15973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15985421bac9SToby Isaac } 15995421bac9SToby Isaac 160066976f2fSJacob Faibussowitsch static PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened) 1601d71ae5a4SJacob Faibussowitsch { 16025421bac9SToby Isaac DMLabel coarsen; 16035421bac9SToby Isaac DM coarseDM; 16045421bac9SToby Isaac 16055421bac9SToby Isaac PetscFunctionBegin; 1606b42326f7SStefano Zampini if (comm != MPI_COMM_NULL) { 16074098eed7SToby Isaac PetscMPIInt mpiComparison; 16084098eed7SToby Isaac MPI_Comm dmcomm = PetscObjectComm((PetscObject)dm); 16094098eed7SToby Isaac 16109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(comm, dmcomm, &mpiComparison)); 16111dca8a05SBarry Smith PetscCheck(mpiComparison == MPI_IDENT || mpiComparison == MPI_CONGRUENT, dmcomm, PETSC_ERR_SUP, "No support for different communicators yet"); 16124098eed7SToby Isaac } 16139566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(dm, &coarseDM)); 16145421bac9SToby Isaac if (coarseDM) { 16159566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)coarseDM)); 16165421bac9SToby Isaac *dmCoarsened = coarseDM; 16173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16185421bac9SToby Isaac } 16199566063dSJacob Faibussowitsch PetscCall(DMForestTemplate(dm, comm, dmCoarsened)); 16209566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityPurpose(*dmCoarsened, DM_ADAPT_COARSEN)); 16219566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm, "coarsen", &coarsen)); 16225421bac9SToby Isaac if (!coarsen) { 16239566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "coarsen", &coarsen)); 16249566063dSJacob Faibussowitsch PetscCall(DMLabelSetDefaultValue(coarsen, DM_ADAPT_COARSEN)); 16251baa6e33SBarry Smith } else PetscCall(PetscObjectReference((PetscObject)coarsen)); 16269566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityLabel(*dmCoarsened, coarsen)); 16279566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&coarsen)); 16283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16295421bac9SToby Isaac } 16305421bac9SToby Isaac 1631d71ae5a4SJacob Faibussowitsch PetscErrorCode DMAdaptLabel_Forest(DM dm, PETSC_UNUSED Vec metric, DMLabel label, PETSC_UNUSED DMLabel rgLabel, DM *adaptedDM) 1632d71ae5a4SJacob Faibussowitsch { 163309350103SToby Isaac PetscBool success; 163409350103SToby Isaac 163509350103SToby Isaac PetscFunctionBegin; 16369566063dSJacob Faibussowitsch PetscCall(DMForestTemplate(dm, PetscObjectComm((PetscObject)dm), adaptedDM)); 16379566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityLabel(*adaptedDM, label)); 16389566063dSJacob Faibussowitsch PetscCall(DMSetUp(*adaptedDM)); 16399566063dSJacob Faibussowitsch PetscCall(DMForestGetAdaptivitySuccess(*adaptedDM, &success)); 164009350103SToby Isaac if (!success) { 16419566063dSJacob Faibussowitsch PetscCall(DMDestroy(adaptedDM)); 164209350103SToby Isaac *adaptedDM = NULL; 164309350103SToby Isaac } 16443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 164509350103SToby Isaac } 164609350103SToby Isaac 1647d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMInitialize_Forest(DM dm) 1648d71ae5a4SJacob Faibussowitsch { 1649d222f98bSToby Isaac PetscFunctionBegin; 1650*f4f49eeaSPierre Jolivet PetscCall(PetscMemzero(dm->ops, sizeof(*dm->ops))); 1651d222f98bSToby Isaac 1652d222f98bSToby Isaac dm->ops->clone = DMClone_Forest; 1653d222f98bSToby Isaac dm->ops->setfromoptions = DMSetFromOptions_Forest; 1654d222f98bSToby Isaac dm->ops->destroy = DMDestroy_Forest; 1655d8984e3bSMatthew G. Knepley dm->ops->createsubdm = DMCreateSubDM_Forest; 16565421bac9SToby Isaac dm->ops->refine = DMRefine_Forest; 16575421bac9SToby Isaac dm->ops->coarsen = DMCoarsen_Forest; 16583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1659d222f98bSToby Isaac } 1660d222f98bSToby Isaac 16619be51f97SToby Isaac /*MC 16629be51f97SToby Isaac 166320f4b53cSBarry Smith DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh. Forests usually have a base `DM` 166420f4b53cSBarry Smith (see `DMForestGetBaseDM()`), from which it is refined. The refinement and partitioning of forests is considered 166520f4b53cSBarry Smith immutable after `DMSetUp()` is called. To adapt a mesh, one should call `DMForestTemplate()` to create a new mesh that 166620f4b53cSBarry Smith will default to being identical to it, specify how that mesh should differ, and then calling `DMSetUp()` on the new 1667bae1f979SBarry Smith mesh. 1668bae1f979SBarry Smith 1669bae1f979SBarry Smith To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the 167020f4b53cSBarry Smith previous mesh whose values indicate which cells should be refined (`DM_ADAPT_REFINE`) or coarsened (`DM_ADAPT_COARSEN`) 1671bae1f979SBarry Smith and how (subtypes are free to allow additional values for things like anisotropic refinement). The label should be 167220f4b53cSBarry Smith given to the *new* mesh with `DMForestSetAdaptivityLabel()`. 16739be51f97SToby Isaac 16749be51f97SToby Isaac Level: advanced 16759be51f97SToby Isaac 167620f4b53cSBarry Smith .seealso: `DMType`, `DM`, `DMCreate()`, `DMSetType()`, `DMForestGetBaseDM()`, `DMForestSetBaseDM()`, `DMForestTemplate()`, `DMForestSetAdaptivityLabel()` 16779be51f97SToby Isaac M*/ 16789be51f97SToby Isaac 1679d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm) 1680d71ae5a4SJacob Faibussowitsch { 1681db4d5e8cSToby Isaac DM_Forest *forest; 1682db4d5e8cSToby Isaac 1683db4d5e8cSToby Isaac PetscFunctionBegin; 1684db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 16854dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&forest)); 1686db4d5e8cSToby Isaac dm->dim = 0; 1687db4d5e8cSToby Isaac dm->data = forest; 1688db4d5e8cSToby Isaac forest->refct = 1; 1689db4d5e8cSToby Isaac forest->data = NULL; 1690db4d5e8cSToby Isaac forest->topology = NULL; 1691cd3c525cSToby Isaac forest->adapt = NULL; 1692db4d5e8cSToby Isaac forest->base = NULL; 16936a87ffbfSToby Isaac forest->adaptPurpose = DM_ADAPT_DETERMINE; 1694db4d5e8cSToby Isaac forest->adjDim = PETSC_DEFAULT; 1695db4d5e8cSToby Isaac forest->overlap = PETSC_DEFAULT; 1696db4d5e8cSToby Isaac forest->minRefinement = PETSC_DEFAULT; 1697db4d5e8cSToby Isaac forest->maxRefinement = PETSC_DEFAULT; 169856ba9f64SToby Isaac forest->initRefinement = PETSC_DEFAULT; 1699c7eeac06SToby Isaac forest->cStart = PETSC_DETERMINE; 1700c7eeac06SToby Isaac forest->cEnd = PETSC_DETERMINE; 1701cd3c525cSToby Isaac forest->cellSF = NULL; 1702ebdf65a2SToby Isaac forest->adaptLabel = NULL; 1703db4d5e8cSToby Isaac forest->gradeFactor = 2; 1704db4d5e8cSToby Isaac forest->cellWeights = NULL; 1705db4d5e8cSToby Isaac forest->cellWeightsCopyMode = PETSC_USE_POINTER; 1706db4d5e8cSToby Isaac forest->weightsFactor = 1.; 1707db4d5e8cSToby Isaac forest->weightCapacity = 1.; 17089566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityStrategy(dm, DMFORESTADAPTALL)); 17099566063dSJacob Faibussowitsch PetscCall(DMInitialize_Forest(dm)); 17103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1711db4d5e8cSToby Isaac } 1712