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 1027d4645fSToby Isaac struct _DMForestTypeLink 1127d4645fSToby Isaac { 1227d4645fSToby Isaac char *name; 1327d4645fSToby Isaac DMForestTypeLink next; 1427d4645fSToby Isaac }; 1527d4645fSToby Isaac 1627d4645fSToby Isaac DMForestTypeLink DMForestTypeList; 1727d4645fSToby Isaac 1827d4645fSToby Isaac static PetscErrorCode DMForestPackageFinalize(void) 1927d4645fSToby Isaac { 2027d4645fSToby Isaac DMForestTypeLink oldLink, link = DMForestTypeList; 2127d4645fSToby Isaac 2227d4645fSToby Isaac PetscFunctionBegin; 2327d4645fSToby Isaac while (link) { 2427d4645fSToby Isaac oldLink = link; 259566063dSJacob Faibussowitsch PetscCall(PetscFree(oldLink->name)); 2627d4645fSToby Isaac link = oldLink->next; 279566063dSJacob Faibussowitsch PetscCall(PetscFree(oldLink)); 2827d4645fSToby Isaac } 2927d4645fSToby Isaac PetscFunctionReturn(0); 3027d4645fSToby Isaac } 3127d4645fSToby Isaac 3227d4645fSToby Isaac static PetscErrorCode DMForestPackageInitialize(void) 3327d4645fSToby Isaac { 3427d4645fSToby Isaac PetscFunctionBegin; 3527d4645fSToby Isaac if (DMForestPackageInitialized) PetscFunctionReturn(0); 3627d4645fSToby Isaac DMForestPackageInitialized = PETSC_TRUE; 37f885a11aSToby Isaac 389566063dSJacob Faibussowitsch PetscCall(DMForestRegisterType(DMFOREST)); 399566063dSJacob Faibussowitsch PetscCall(PetscRegisterFinalize(DMForestPackageFinalize)); 4027d4645fSToby Isaac PetscFunctionReturn(0); 4127d4645fSToby Isaac } 4227d4645fSToby Isaac 439be51f97SToby Isaac /*@C 449be51f97SToby Isaac DMForestRegisterType - Registers a DMType as a subtype of DMFOREST (so that DMIsForest() will be correct) 459be51f97SToby Isaac 469be51f97SToby Isaac Not Collective 479be51f97SToby Isaac 489be51f97SToby Isaac Input parameter: 499be51f97SToby Isaac . name - the name of the type 509be51f97SToby Isaac 519be51f97SToby Isaac Level: advanced 529be51f97SToby Isaac 539be51f97SToby Isaac .seealso: DMFOREST, DMIsForest() 549be51f97SToby Isaac @*/ 5527d4645fSToby Isaac PetscErrorCode DMForestRegisterType(DMType name) 5627d4645fSToby Isaac { 5727d4645fSToby Isaac DMForestTypeLink link; 5827d4645fSToby Isaac 5927d4645fSToby Isaac PetscFunctionBegin; 609566063dSJacob Faibussowitsch PetscCall(DMForestPackageInitialize()); 619566063dSJacob Faibussowitsch PetscCall(PetscNew(&link)); 629566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(name,&link->name)); 6327d4645fSToby Isaac link->next = DMForestTypeList; 6427d4645fSToby Isaac DMForestTypeList = link; 6527d4645fSToby Isaac PetscFunctionReturn(0); 6627d4645fSToby Isaac } 6727d4645fSToby Isaac 689be51f97SToby Isaac /*@ 699be51f97SToby Isaac DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes 709be51f97SToby Isaac 719be51f97SToby Isaac Not Collective 729be51f97SToby Isaac 739be51f97SToby Isaac Input parameter: 749be51f97SToby Isaac . dm - the DM object 759be51f97SToby Isaac 769be51f97SToby Isaac Output parameter: 779be51f97SToby Isaac . isForest - whether dm is a subtype of DMFOREST 789be51f97SToby Isaac 799be51f97SToby Isaac Level: intermediate 809be51f97SToby Isaac 819be51f97SToby Isaac .seealso: DMFOREST, DMForestRegisterType() 829be51f97SToby Isaac @*/ 8327d4645fSToby Isaac PetscErrorCode DMIsForest(DM dm, PetscBool *isForest) 8427d4645fSToby Isaac { 8527d4645fSToby Isaac DMForestTypeLink link = DMForestTypeList; 8627d4645fSToby Isaac 8727d4645fSToby Isaac PetscFunctionBegin; 8827d4645fSToby Isaac while (link) { 8927d4645fSToby Isaac PetscBool sameType; 909566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)dm,link->name,&sameType)); 9127d4645fSToby Isaac if (sameType) { 9227d4645fSToby Isaac *isForest = PETSC_TRUE; 9327d4645fSToby Isaac PetscFunctionReturn(0); 9427d4645fSToby Isaac } 9527d4645fSToby Isaac link = link->next; 9627d4645fSToby Isaac } 9727d4645fSToby Isaac *isForest = PETSC_FALSE; 9827d4645fSToby Isaac PetscFunctionReturn(0); 9927d4645fSToby Isaac } 10027d4645fSToby Isaac 1019be51f97SToby Isaac /*@ 1029be51f97SToby Isaac DMForestTemplate - Create a new DM that will be adapted from a source DM. The new DM reproduces the configuration 1039be51f97SToby Isaac of the source, but is not yet setup, so that the user can then define only the ways that the new DM should differ 1049be51f97SToby Isaac (by, e.g., refinement or repartitioning). The source DM is also set as the adaptivity source DM of the new DM (see 1059be51f97SToby Isaac DMForestSetAdaptivityForest()). 1069be51f97SToby Isaac 1079be51f97SToby Isaac Collective on dm 1089be51f97SToby Isaac 1099be51f97SToby Isaac Input Parameters: 1109be51f97SToby Isaac + dm - the source DM object 1119be51f97SToby Isaac - comm - the communicator for the new DM (this communicator is currently ignored, but is present so that DMForestTemplate() can be used within DMCoarsen()) 1129be51f97SToby Isaac 1139be51f97SToby Isaac Output Parameter: 1149be51f97SToby Isaac . tdm - the new DM object 1159be51f97SToby Isaac 1169be51f97SToby Isaac Level: intermediate 1179be51f97SToby Isaac 1189be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest() 1199be51f97SToby Isaac @*/ 1209be51f97SToby Isaac PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm) 121a0452a8eSToby Isaac { 122a0452a8eSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 12320e8089bSToby Isaac DMType type; 124a0452a8eSToby Isaac DM base; 125a0452a8eSToby Isaac DMForestTopology topology; 12605e99e11SStefano Zampini MatType mtype; 127a0452a8eSToby Isaac PetscInt dim, overlap, ref, factor; 128a0452a8eSToby Isaac DMForestAdaptivityStrategy strat; 129795844e7SToby Isaac void *ctx; 13049fc9a2fSToby Isaac PetscErrorCode (*map)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void*); 1313e58adeeSToby Isaac void *mapCtx; 132a0452a8eSToby Isaac 133a0452a8eSToby Isaac PetscFunctionBegin; 134a0452a8eSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1359566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm),tdm)); 1369566063dSJacob Faibussowitsch PetscCall(DMGetType(dm,&type)); 1379566063dSJacob Faibussowitsch PetscCall(DMSetType(*tdm,type)); 1389566063dSJacob Faibussowitsch PetscCall(DMForestGetBaseDM(dm,&base)); 1399566063dSJacob Faibussowitsch PetscCall(DMForestSetBaseDM(*tdm,base)); 1409566063dSJacob Faibussowitsch PetscCall(DMForestGetTopology(dm,&topology)); 1419566063dSJacob Faibussowitsch PetscCall(DMForestSetTopology(*tdm,topology)); 1429566063dSJacob Faibussowitsch PetscCall(DMForestGetAdjacencyDimension(dm,&dim)); 1439566063dSJacob Faibussowitsch PetscCall(DMForestSetAdjacencyDimension(*tdm,dim)); 1449566063dSJacob Faibussowitsch PetscCall(DMForestGetPartitionOverlap(dm,&overlap)); 1459566063dSJacob Faibussowitsch PetscCall(DMForestSetPartitionOverlap(*tdm,overlap)); 1469566063dSJacob Faibussowitsch PetscCall(DMForestGetMinimumRefinement(dm,&ref)); 1479566063dSJacob Faibussowitsch PetscCall(DMForestSetMinimumRefinement(*tdm,ref)); 1489566063dSJacob Faibussowitsch PetscCall(DMForestGetMaximumRefinement(dm,&ref)); 1499566063dSJacob Faibussowitsch PetscCall(DMForestSetMaximumRefinement(*tdm,ref)); 1509566063dSJacob Faibussowitsch PetscCall(DMForestGetAdaptivityStrategy(dm,&strat)); 1519566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityStrategy(*tdm,strat)); 1529566063dSJacob Faibussowitsch PetscCall(DMForestGetGradeFactor(dm,&factor)); 1539566063dSJacob Faibussowitsch PetscCall(DMForestSetGradeFactor(*tdm,factor)); 1549566063dSJacob Faibussowitsch PetscCall(DMForestGetBaseCoordinateMapping(dm,&map,&mapCtx)); 1559566063dSJacob Faibussowitsch PetscCall(DMForestSetBaseCoordinateMapping(*tdm,map,mapCtx)); 156a0452a8eSToby Isaac if (forest->ftemplate) { 1579566063dSJacob Faibussowitsch PetscCall((*forest->ftemplate)(dm, *tdm)); 158a0452a8eSToby Isaac } 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 { 16490b157c4SStefano Zampini PetscBool isper; 165795844e7SToby Isaac const PetscReal *maxCell, *L; 166795844e7SToby Isaac const DMBoundaryType *bd; 167795844e7SToby Isaac 1689566063dSJacob Faibussowitsch PetscCall(DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd)); 1699566063dSJacob Faibussowitsch PetscCall(DMSetPeriodicity(*tdm,isper,maxCell,L,bd)); 170795844e7SToby Isaac } 1719566063dSJacob Faibussowitsch PetscCall(DMGetMatType(dm,&mtype)); 1729566063dSJacob Faibussowitsch PetscCall(DMSetMatType(*tdm,mtype)); 173a0452a8eSToby Isaac PetscFunctionReturn(0); 174a0452a8eSToby Isaac } 175a0452a8eSToby Isaac 17601d9d024SToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm); 17701d9d024SToby Isaac 178db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm) 179db4d5e8cSToby Isaac { 180db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 181db4d5e8cSToby Isaac const char *type; 182db4d5e8cSToby Isaac 183db4d5e8cSToby Isaac PetscFunctionBegin; 184db4d5e8cSToby Isaac forest->refct++; 185db4d5e8cSToby Isaac (*newdm)->data = forest; 1869566063dSJacob Faibussowitsch PetscCall(PetscObjectGetType((PetscObject) dm, &type)); 1879566063dSJacob Faibussowitsch PetscCall(PetscObjectChangeTypeName((PetscObject) *newdm, type)); 1889566063dSJacob Faibussowitsch PetscCall(DMInitialize_Forest(*newdm)); 189db4d5e8cSToby Isaac PetscFunctionReturn(0); 190db4d5e8cSToby Isaac } 191db4d5e8cSToby Isaac 192d222f98bSToby Isaac static PetscErrorCode DMDestroy_Forest(DM dm) 193db4d5e8cSToby Isaac { 194db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 195db4d5e8cSToby Isaac 196db4d5e8cSToby Isaac PetscFunctionBegin; 197db4d5e8cSToby Isaac if (--forest->refct > 0) PetscFunctionReturn(0); 1989566063dSJacob Faibussowitsch if (forest->destroy) PetscCall((*forest->destroy)(dm)); 1999566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&forest->cellSF)); 2009566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&forest->preCoarseToFine)); 2019566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&forest->coarseToPreFine)); 2029566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&forest->adaptLabel)); 2039566063dSJacob Faibussowitsch PetscCall(PetscFree(forest->adaptStrategy)); 2049566063dSJacob Faibussowitsch PetscCall(DMDestroy(&forest->base)); 2059566063dSJacob Faibussowitsch PetscCall(DMDestroy(&forest->adapt)); 2069566063dSJacob Faibussowitsch PetscCall(PetscFree(forest->topology)); 2079566063dSJacob Faibussowitsch PetscCall(PetscFree(forest)); 208db4d5e8cSToby Isaac PetscFunctionReturn(0); 209db4d5e8cSToby Isaac } 210db4d5e8cSToby Isaac 2119be51f97SToby Isaac /*@C 2129be51f97SToby Isaac DMForestSetTopology - Set the topology of a DMForest during the pre-setup phase. The topology is a string (e.g. 21368d54884SBarry Smith "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest during 2149be51f97SToby Isaac DMSetUp(). 2159be51f97SToby Isaac 2169be51f97SToby Isaac Logically collective on dm 2179be51f97SToby Isaac 2189be51f97SToby Isaac Input parameters: 2199be51f97SToby Isaac + dm - the forest 2209be51f97SToby Isaac - topology - the topology of the forest 2219be51f97SToby Isaac 2229be51f97SToby Isaac Level: intermediate 2239be51f97SToby Isaac 224fd292e60Sprj- .seealso: DMForestGetTopology(), DMForestSetBaseDM() 2259be51f97SToby Isaac @*/ 226dd8e54a2SToby Isaac PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology) 227db4d5e8cSToby Isaac { 228db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 229db4d5e8cSToby Isaac 230db4d5e8cSToby Isaac PetscFunctionBegin; 231db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 23228b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the topology after setup"); 2339566063dSJacob Faibussowitsch PetscCall(PetscFree(forest->topology)); 2349566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy((const char*)topology,(char**) &forest->topology)); 235db4d5e8cSToby Isaac PetscFunctionReturn(0); 236db4d5e8cSToby Isaac } 237db4d5e8cSToby Isaac 2389be51f97SToby Isaac /*@C 2399be51f97SToby Isaac DMForestGetTopology - Get a string describing the topology of a DMForest. 2409be51f97SToby Isaac 2419be51f97SToby Isaac Not collective 2429be51f97SToby Isaac 2439be51f97SToby Isaac Input parameter: 2449be51f97SToby Isaac . dm - the forest 2459be51f97SToby Isaac 2469be51f97SToby Isaac Output parameter: 2479be51f97SToby Isaac . topology - the topology of the forest (e.g., 'cube', 'shell') 2489be51f97SToby Isaac 2499be51f97SToby Isaac Level: intermediate 2509be51f97SToby Isaac 2519be51f97SToby Isaac .seealso: DMForestSetTopology() 2529be51f97SToby Isaac @*/ 253dd8e54a2SToby Isaac PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology) 254dd8e54a2SToby Isaac { 255dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 256dd8e54a2SToby Isaac 257dd8e54a2SToby Isaac PetscFunctionBegin; 258dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 259dd8e54a2SToby Isaac PetscValidPointer(topology,2); 260dd8e54a2SToby Isaac *topology = forest->topology; 261dd8e54a2SToby Isaac PetscFunctionReturn(0); 262dd8e54a2SToby Isaac } 263dd8e54a2SToby Isaac 2649be51f97SToby Isaac /*@ 2659be51f97SToby Isaac DMForestSetBaseDM - During the pre-setup phase, set the DM that defines the base mesh of a DMForest forest. The 2669be51f97SToby Isaac forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its 267765b024eSBarry Smith base. In general, two forest must share a base to be comparable, to do things like construct interpolators. 2689be51f97SToby Isaac 2699be51f97SToby Isaac Logically collective on dm 2709be51f97SToby Isaac 2719be51f97SToby Isaac Input Parameters: 2729be51f97SToby Isaac + dm - the forest 2739be51f97SToby Isaac - base - the base DM of the forest 2749be51f97SToby Isaac 27595452b02SPatrick Sanan Notes: 27695452b02SPatrick Sanan Currently the base DM must be a DMPLEX 277765b024eSBarry Smith 2789be51f97SToby Isaac Level: intermediate 2799be51f97SToby Isaac 280fd292e60Sprj- .seealso: DMForestGetBaseDM() 2819be51f97SToby Isaac @*/ 282dd8e54a2SToby Isaac PetscErrorCode DMForestSetBaseDM(DM dm, DM base) 283dd8e54a2SToby Isaac { 284dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 285dd8e54a2SToby Isaac PetscInt dim, dimEmbed; 286dd8e54a2SToby Isaac 287dd8e54a2SToby Isaac PetscFunctionBegin; 288dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 28928b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the base after setup"); 2909566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)base)); 2919566063dSJacob Faibussowitsch PetscCall(DMDestroy(&forest->base)); 292dd8e54a2SToby Isaac forest->base = base; 293a0452a8eSToby Isaac if (base) { 29428dfcf7cSStefano Zampini PetscBool isper; 29528dfcf7cSStefano Zampini const PetscReal *maxCell, *L; 29628dfcf7cSStefano Zampini const DMBoundaryType *bd; 29728dfcf7cSStefano Zampini 298a0452a8eSToby Isaac PetscValidHeaderSpecific(base, DM_CLASSID, 2); 2999566063dSJacob Faibussowitsch PetscCall(DMGetDimension(base,&dim)); 3009566063dSJacob Faibussowitsch PetscCall(DMSetDimension(dm,dim)); 3019566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDim(base,&dimEmbed)); 3029566063dSJacob Faibussowitsch PetscCall(DMSetCoordinateDim(dm,dimEmbed)); 3039566063dSJacob Faibussowitsch PetscCall(DMGetPeriodicity(base,&isper,&maxCell,&L,&bd)); 3049566063dSJacob Faibussowitsch PetscCall(DMSetPeriodicity(dm,isper,maxCell,L,bd)); 30528dfcf7cSStefano Zampini } else { 3069566063dSJacob Faibussowitsch PetscCall(DMSetPeriodicity(dm,PETSC_FALSE,NULL,NULL,NULL)); 307a0452a8eSToby Isaac } 308dd8e54a2SToby Isaac PetscFunctionReturn(0); 309dd8e54a2SToby Isaac } 310dd8e54a2SToby Isaac 3119be51f97SToby Isaac /*@ 3129be51f97SToby Isaac DMForestGetBaseDM - Get the base DM of a DMForest forest. The forest will be hierarchically refined from the base, 31368d54884SBarry Smith and all refinements/coarsenings of the forest will share its base. In general, two forest must share a base to be 3149be51f97SToby Isaac comparable, to do things like construct interpolators. 3159be51f97SToby Isaac 3169be51f97SToby Isaac Not collective 3179be51f97SToby Isaac 3189be51f97SToby Isaac Input Parameter: 3199be51f97SToby Isaac . dm - the forest 3209be51f97SToby Isaac 3219be51f97SToby Isaac Output Parameter: 3229be51f97SToby Isaac . base - the base DM of the forest 3239be51f97SToby Isaac 324367003a6SStefano Zampini Notes: 325367003a6SStefano Zampini After DMSetUp(), the base DM will be redundantly distributed across MPI processes 326367003a6SStefano Zampini 3279be51f97SToby Isaac Level: intermediate 3289be51f97SToby Isaac 329fd292e60Sprj- .seealso: DMForestSetBaseDM() 3309be51f97SToby Isaac @*/ 331dd8e54a2SToby Isaac PetscErrorCode DMForestGetBaseDM(DM dm, DM *base) 332dd8e54a2SToby Isaac { 333dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 334dd8e54a2SToby Isaac 335dd8e54a2SToby Isaac PetscFunctionBegin; 336dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 337dd8e54a2SToby Isaac PetscValidPointer(base, 2); 338dd8e54a2SToby Isaac *base = forest->base; 339dd8e54a2SToby Isaac PetscFunctionReturn(0); 340dd8e54a2SToby Isaac } 341dd8e54a2SToby Isaac 34299478f86SToby Isaac PetscErrorCode DMForestSetBaseCoordinateMapping(DM dm, PetscErrorCode (*func)(DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx) 343cf38a08cSToby Isaac { 344cf38a08cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 345cf38a08cSToby Isaac 346cf38a08cSToby Isaac PetscFunctionBegin; 347cf38a08cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 348cf38a08cSToby Isaac forest->mapcoordinates = func; 349cf38a08cSToby Isaac forest->mapcoordinatesctx = ctx; 350cf38a08cSToby Isaac PetscFunctionReturn(0); 351cf38a08cSToby Isaac } 352cf38a08cSToby Isaac 35399478f86SToby Isaac PetscErrorCode DMForestGetBaseCoordinateMapping(DM dm, PetscErrorCode (**func) (DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx) 354cf38a08cSToby Isaac { 355cf38a08cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 356cf38a08cSToby Isaac 357cf38a08cSToby Isaac PetscFunctionBegin; 358cf38a08cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 359cf38a08cSToby Isaac if (func) *func = forest->mapcoordinates; 360cf38a08cSToby Isaac if (ctx) *((void**) ctx) = forest->mapcoordinatesctx; 361cf38a08cSToby Isaac PetscFunctionReturn(0); 362cf38a08cSToby Isaac } 363cf38a08cSToby Isaac 3649be51f97SToby Isaac /*@ 3659be51f97SToby Isaac DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be 3669be51f97SToby Isaac adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) im DMSetUp(). Usually not needed 3679be51f97SToby Isaac by users directly: DMForestTemplate() constructs a new forest to be adapted from an old forest and calls this 3689be51f97SToby Isaac routine. 3699be51f97SToby Isaac 370dffe73a3SToby Isaac Note that this can be called after setup with adapt = NULL, which will clear all internal data related to the 371dffe73a3SToby Isaac adaptivity forest from dm. This way, repeatedly adapting does not leave stale DM objects in memory. 372dffe73a3SToby Isaac 3739be51f97SToby Isaac Logically collective on dm 3749be51f97SToby Isaac 375d8d19677SJose E. Roman Input Parameters: 3769be51f97SToby Isaac + dm - the new forest, which will be constructed from adapt 3779be51f97SToby Isaac - adapt - the old forest 3789be51f97SToby Isaac 3799be51f97SToby Isaac Level: intermediate 3809be51f97SToby Isaac 3819be51f97SToby Isaac .seealso: DMForestGetAdaptivityForest(), DMForestSetAdaptivityPurpose() 3829be51f97SToby Isaac @*/ 383ba936b91SToby Isaac PetscErrorCode DMForestSetAdaptivityForest(DM dm,DM adapt) 384dd8e54a2SToby Isaac { 385dffe73a3SToby Isaac DM_Forest *forest, *adaptForest, *oldAdaptForest; 386dffe73a3SToby Isaac DM oldAdapt; 387456cc5b7SMatthew G. Knepley PetscBool isForest; 388dd8e54a2SToby Isaac 389dd8e54a2SToby Isaac PetscFunctionBegin; 390dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3911fd50544SStefano Zampini if (adapt) PetscValidHeaderSpecific(adapt, DM_CLASSID, 2); 3929566063dSJacob Faibussowitsch PetscCall(DMIsForest(dm, &isForest)); 393456cc5b7SMatthew G. Knepley if (!isForest) PetscFunctionReturn(0); 394*1dca8a05SBarry Smith PetscCheck(adapt == NULL || !dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup"); 395ba936b91SToby Isaac forest = (DM_Forest*) dm->data; 3969566063dSJacob Faibussowitsch PetscCall(DMForestGetAdaptivityForest(dm,&oldAdapt)); 397193eb951SToby Isaac adaptForest = (DM_Forest*) (adapt ? adapt->data : NULL); 398193eb951SToby Isaac oldAdaptForest = (DM_Forest*) (oldAdapt ? oldAdapt->data : NULL); 399dffe73a3SToby Isaac if (adaptForest != oldAdaptForest) { 4009566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&forest->preCoarseToFine)); 4019566063dSJacob Faibussowitsch PetscCall(PetscSFDestroy(&forest->coarseToPreFine)); 4029566063dSJacob Faibussowitsch if (forest->clearadaptivityforest) PetscCall((*forest->clearadaptivityforest)(dm)); 403dffe73a3SToby Isaac } 40426d9498aSToby Isaac switch (forest->adaptPurpose) { 405cd3c525cSToby Isaac case DM_ADAPT_DETERMINE: 4069566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)adapt)); 4079566063dSJacob Faibussowitsch PetscCall(DMDestroy(&(forest->adapt))); 408ba936b91SToby Isaac forest->adapt = adapt; 40926d9498aSToby Isaac break; 410a1b0c543SToby Isaac case DM_ADAPT_REFINE: 4119566063dSJacob Faibussowitsch PetscCall(DMSetCoarseDM(dm,adapt)); 41226d9498aSToby Isaac break; 413a1b0c543SToby Isaac case DM_ADAPT_COARSEN: 414bf2d5fbbSStefano Zampini case DM_ADAPT_COARSEN_LAST: 4159566063dSJacob Faibussowitsch PetscCall(DMSetFineDM(dm,adapt)); 41626d9498aSToby Isaac break; 41726d9498aSToby Isaac default: 41826d9498aSToby Isaac SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose"); 41926d9498aSToby Isaac } 420dd8e54a2SToby Isaac PetscFunctionReturn(0); 421dd8e54a2SToby Isaac } 422dd8e54a2SToby Isaac 4239be51f97SToby Isaac /*@ 4249be51f97SToby Isaac DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted. 4259be51f97SToby Isaac 4269be51f97SToby Isaac Not collective 4279be51f97SToby Isaac 4289be51f97SToby Isaac Input Parameter: 4299be51f97SToby Isaac . dm - the forest 4309be51f97SToby Isaac 4319be51f97SToby Isaac Output Parameter: 4329be51f97SToby Isaac . adapt - the forest from which dm is/was adapted 4339be51f97SToby Isaac 4349be51f97SToby Isaac Level: intermediate 4359be51f97SToby Isaac 4369be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose() 4379be51f97SToby Isaac @*/ 438ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt) 439dd8e54a2SToby Isaac { 440ba936b91SToby Isaac DM_Forest *forest; 441dd8e54a2SToby Isaac 442dd8e54a2SToby Isaac PetscFunctionBegin; 443dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 444ba936b91SToby Isaac forest = (DM_Forest*) dm->data; 44526d9498aSToby Isaac switch (forest->adaptPurpose) { 446cd3c525cSToby Isaac case DM_ADAPT_DETERMINE: 447ba936b91SToby Isaac *adapt = forest->adapt; 44826d9498aSToby Isaac break; 449a1b0c543SToby Isaac case DM_ADAPT_REFINE: 4509566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(dm,adapt)); 45126d9498aSToby Isaac break; 452a1b0c543SToby Isaac case DM_ADAPT_COARSEN: 453bf2d5fbbSStefano Zampini case DM_ADAPT_COARSEN_LAST: 4549566063dSJacob Faibussowitsch PetscCall(DMGetFineDM(dm,adapt)); 45526d9498aSToby Isaac break; 45626d9498aSToby Isaac default: 45726d9498aSToby Isaac SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose"); 45826d9498aSToby Isaac } 45926d9498aSToby Isaac PetscFunctionReturn(0); 46026d9498aSToby Isaac } 46126d9498aSToby Isaac 4629be51f97SToby Isaac /*@ 4639be51f97SToby Isaac DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its 464a1b0c543SToby Isaac source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening 465cd3c525cSToby Isaac (DM_ADAPT_COARSEN), or undefined (DM_ADAPT_DETERMINE). This only matters for the purposes of reference counting: 4669be51f97SToby Isaac during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse 4679be51f97SToby Isaac relationship (see DMSetFineDM()/DMSetCoarseDM()). If the purpose is not refinement or coarsening, and the user does 4689be51f97SToby Isaac not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can 4699be51f97SToby Isaac cause a memory leak. This method is used by subtypes of DMForest when automatically constructing mesh hierarchies. 4709be51f97SToby Isaac 4719be51f97SToby Isaac Logically collective on dm 4729be51f97SToby Isaac 4739be51f97SToby Isaac Input Parameters: 4749be51f97SToby Isaac + dm - the forest 475bf2d5fbbSStefano Zampini - purpose - the adaptivity purpose 4769be51f97SToby Isaac 4779be51f97SToby Isaac Level: advanced 4789be51f97SToby Isaac 479bf2d5fbbSStefano Zampini .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest(), DMAdaptFlag 4809be51f97SToby Isaac @*/ 481a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMAdaptFlag purpose) 48226d9498aSToby Isaac { 48326d9498aSToby Isaac DM_Forest *forest; 48426d9498aSToby Isaac 48526d9498aSToby Isaac PetscFunctionBegin; 48626d9498aSToby Isaac forest = (DM_Forest*) dm->data; 48728b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup"); 48826d9498aSToby Isaac if (purpose != forest->adaptPurpose) { 48926d9498aSToby Isaac DM adapt; 49026d9498aSToby Isaac 4919566063dSJacob Faibussowitsch PetscCall(DMForestGetAdaptivityForest(dm,&adapt)); 4929566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)adapt)); 4939566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityForest(dm,NULL)); 494f885a11aSToby Isaac 49526d9498aSToby Isaac forest->adaptPurpose = purpose; 496f885a11aSToby Isaac 4979566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityForest(dm,adapt)); 4989566063dSJacob Faibussowitsch PetscCall(DMDestroy(&adapt)); 49926d9498aSToby Isaac } 500dd8e54a2SToby Isaac PetscFunctionReturn(0); 501dd8e54a2SToby Isaac } 502dd8e54a2SToby Isaac 50356c0450aSToby Isaac /*@ 50456c0450aSToby Isaac DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with 505bf2d5fbbSStefano Zampini DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening (DM_ADAPT_COARSEN), 506bf2d5fbbSStefano Zampini coarsening only the last level (DM_ADAPT_COARSEN_LAST) or undefined (DM_ADAPT_DETERMINE). 507bf2d5fbbSStefano Zampini This only matters for the purposes of reference counting: during DMDestroy(), cyclic 50856c0450aSToby Isaac references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see 50956c0450aSToby Isaac DMSetFineDM()/DMSetCoarseDM()). If the purpose is not refinement or coarsening, and the user does not maintain a 51056c0450aSToby Isaac reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory 51156c0450aSToby Isaac leak. This method is used by subtypes of DMForest when automatically constructing mesh hierarchies. 51256c0450aSToby Isaac 51356c0450aSToby Isaac Not collective 51456c0450aSToby Isaac 51556c0450aSToby Isaac Input Parameter: 51656c0450aSToby Isaac . dm - the forest 51756c0450aSToby Isaac 51856c0450aSToby Isaac Output Parameter: 519bf2d5fbbSStefano Zampini . purpose - the adaptivity purpose 52056c0450aSToby Isaac 52156c0450aSToby Isaac Level: advanced 52256c0450aSToby Isaac 523bf2d5fbbSStefano Zampini .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest(), DMAdaptFlag 52456c0450aSToby Isaac @*/ 525a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMAdaptFlag *purpose) 52656c0450aSToby Isaac { 52756c0450aSToby Isaac DM_Forest *forest; 52856c0450aSToby Isaac 52956c0450aSToby Isaac PetscFunctionBegin; 53056c0450aSToby Isaac forest = (DM_Forest*) dm->data; 53156c0450aSToby Isaac *purpose = forest->adaptPurpose; 53256c0450aSToby Isaac PetscFunctionReturn(0); 53356c0450aSToby Isaac } 53456c0450aSToby Isaac 5359be51f97SToby Isaac /*@ 5369be51f97SToby Isaac DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine 5379be51f97SToby Isaac cell adjacency (for the purposes of partitioning and overlap). 5389be51f97SToby Isaac 5399be51f97SToby Isaac Logically collective on dm 5409be51f97SToby Isaac 5419be51f97SToby Isaac Input Parameters: 5429be51f97SToby Isaac + dm - the forest 5439be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency) 5449be51f97SToby Isaac 5459be51f97SToby Isaac Level: intermediate 5469be51f97SToby Isaac 5479be51f97SToby Isaac .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap() 5489be51f97SToby Isaac @*/ 549dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim) 550dd8e54a2SToby Isaac { 551dd8e54a2SToby Isaac PetscInt dim; 552dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 553dd8e54a2SToby Isaac 554dd8e54a2SToby Isaac PetscFunctionBegin; 555dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 55628b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup"); 55763a3b9bcSJacob Faibussowitsch PetscCheck(adjDim >= 0,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %" PetscInt_FMT, adjDim); 5589566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm,&dim)); 55963a3b9bcSJacob Faibussowitsch PetscCheck(adjDim <= dim,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %" PetscInt_FMT ": %" PetscInt_FMT, dim, adjDim); 560dd8e54a2SToby Isaac forest->adjDim = adjDim; 561dd8e54a2SToby Isaac PetscFunctionReturn(0); 562dd8e54a2SToby Isaac } 563dd8e54a2SToby Isaac 5649be51f97SToby Isaac /*@ 5659be51f97SToby Isaac DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that, 5669be51f97SToby Isaac e.g., adjacency based on facets can be specified by codimension 1 in all cases) 5679be51f97SToby Isaac 5689be51f97SToby Isaac Logically collective on dm 5699be51f97SToby Isaac 5709be51f97SToby Isaac Input Parameters: 5719be51f97SToby Isaac + dm - the forest 5729be51f97SToby Isaac - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices 5739be51f97SToby Isaac 5749be51f97SToby Isaac Level: intermediate 5759be51f97SToby Isaac 5769be51f97SToby Isaac .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension() 5779be51f97SToby Isaac @*/ 578dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim) 579dd8e54a2SToby Isaac { 580dd8e54a2SToby Isaac PetscInt dim; 581dd8e54a2SToby Isaac 582dd8e54a2SToby Isaac PetscFunctionBegin; 583dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5849566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm,&dim)); 5859566063dSJacob Faibussowitsch PetscCall(DMForestSetAdjacencyDimension(dm,dim-adjCodim)); 586dd8e54a2SToby Isaac PetscFunctionReturn(0); 587dd8e54a2SToby Isaac } 588dd8e54a2SToby Isaac 5899be51f97SToby Isaac /*@ 5909be51f97SToby Isaac DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the 5919be51f97SToby Isaac purposes of partitioning and overlap). 5929be51f97SToby Isaac 5939be51f97SToby Isaac Not collective 5949be51f97SToby Isaac 5959be51f97SToby Isaac Input Parameter: 5969be51f97SToby Isaac . dm - the forest 5979be51f97SToby Isaac 5989be51f97SToby Isaac Output Parameter: 5999be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency) 6009be51f97SToby Isaac 6019be51f97SToby Isaac Level: intermediate 6029be51f97SToby Isaac 6039be51f97SToby Isaac .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap() 6049be51f97SToby Isaac @*/ 605dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim) 606dd8e54a2SToby Isaac { 607dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 608dd8e54a2SToby Isaac 609dd8e54a2SToby Isaac PetscFunctionBegin; 610dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 611dd8e54a2SToby Isaac PetscValidIntPointer(adjDim,2); 612dd8e54a2SToby Isaac *adjDim = forest->adjDim; 613dd8e54a2SToby Isaac PetscFunctionReturn(0); 614dd8e54a2SToby Isaac } 615dd8e54a2SToby Isaac 6169be51f97SToby Isaac /*@ 6179be51f97SToby Isaac DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that, 6189be51f97SToby Isaac e.g., adjacency based on facets can be specified by codimension 1 in all cases) 6199be51f97SToby Isaac 6209be51f97SToby Isaac Not collective 6219be51f97SToby Isaac 6229be51f97SToby Isaac Input Parameter: 6239be51f97SToby Isaac . dm - the forest 6249be51f97SToby Isaac 6259be51f97SToby Isaac Output Parameter: 6269be51f97SToby Isaac . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices 6279be51f97SToby Isaac 6289be51f97SToby Isaac Level: intermediate 6299be51f97SToby Isaac 6309be51f97SToby Isaac .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension() 6319be51f97SToby Isaac @*/ 632dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim) 633dd8e54a2SToby Isaac { 634dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 635dd8e54a2SToby Isaac PetscInt dim; 636dd8e54a2SToby Isaac 637dd8e54a2SToby Isaac PetscFunctionBegin; 638dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 639dd8e54a2SToby Isaac PetscValidIntPointer(adjCodim,2); 6409566063dSJacob Faibussowitsch PetscCall(DMGetDimension(dm,&dim)); 641dd8e54a2SToby Isaac *adjCodim = dim - forest->adjDim; 642dd8e54a2SToby Isaac PetscFunctionReturn(0); 643dd8e54a2SToby Isaac } 644dd8e54a2SToby Isaac 6459be51f97SToby Isaac /*@ 6469be51f97SToby Isaac DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel 6479be51f97SToby Isaac partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding 6489be51f97SToby Isaac adjacent cells 6499be51f97SToby Isaac 6509be51f97SToby Isaac Logically collective on dm 6519be51f97SToby Isaac 6529be51f97SToby Isaac Input Parameters: 6539be51f97SToby Isaac + dm - the forest 6549be51f97SToby Isaac - overlap - default 0 6559be51f97SToby Isaac 6569be51f97SToby Isaac Level: intermediate 6579be51f97SToby Isaac 6589be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension() 6599be51f97SToby Isaac @*/ 660dd8e54a2SToby Isaac PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap) 661dd8e54a2SToby Isaac { 662dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 663dd8e54a2SToby Isaac 664dd8e54a2SToby Isaac PetscFunctionBegin; 665dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 66628b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup"); 66763a3b9bcSJacob Faibussowitsch PetscCheck(overlap >= 0,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %" PetscInt_FMT, overlap); 668dd8e54a2SToby Isaac forest->overlap = overlap; 669dd8e54a2SToby Isaac PetscFunctionReturn(0); 670dd8e54a2SToby Isaac } 671dd8e54a2SToby Isaac 6729be51f97SToby Isaac /*@ 6739be51f97SToby Isaac DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values 6749be51f97SToby Isaac > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells 6759be51f97SToby Isaac 6769be51f97SToby Isaac Not collective 6779be51f97SToby Isaac 6789be51f97SToby Isaac Input Parameter: 6799be51f97SToby Isaac . dm - the forest 6809be51f97SToby Isaac 6819be51f97SToby Isaac Output Parameter: 6829be51f97SToby Isaac . overlap - default 0 6839be51f97SToby Isaac 6849be51f97SToby Isaac Level: intermediate 6859be51f97SToby Isaac 6869be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension() 6879be51f97SToby Isaac @*/ 688dd8e54a2SToby Isaac PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap) 689dd8e54a2SToby Isaac { 690dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 691dd8e54a2SToby Isaac 692dd8e54a2SToby Isaac PetscFunctionBegin; 693dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 694dd8e54a2SToby Isaac PetscValidIntPointer(overlap,2); 695dd8e54a2SToby Isaac *overlap = forest->overlap; 696dd8e54a2SToby Isaac PetscFunctionReturn(0); 697dd8e54a2SToby Isaac } 698dd8e54a2SToby Isaac 6999be51f97SToby Isaac /*@ 7009be51f97SToby Isaac DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base 7019be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. If the forest is being created by coarsening a previous forest 7029be51f97SToby Isaac (see DMForestGetAdaptivityForest()) this limits the amount of coarsening. 7039be51f97SToby Isaac 7049be51f97SToby Isaac Logically collective on dm 7059be51f97SToby Isaac 7069be51f97SToby Isaac Input Parameters: 7079be51f97SToby Isaac + dm - the forest 7089be51f97SToby Isaac - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7099be51f97SToby Isaac 7109be51f97SToby Isaac Level: intermediate 7119be51f97SToby Isaac 7129be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 7139be51f97SToby Isaac @*/ 714dd8e54a2SToby Isaac PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement) 715dd8e54a2SToby Isaac { 716dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 717dd8e54a2SToby Isaac 718dd8e54a2SToby Isaac PetscFunctionBegin; 719dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 72028b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup"); 721dd8e54a2SToby Isaac forest->minRefinement = minRefinement; 722dd8e54a2SToby Isaac PetscFunctionReturn(0); 723dd8e54a2SToby Isaac } 724dd8e54a2SToby Isaac 7259be51f97SToby Isaac /*@ 7269be51f97SToby Isaac DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see 7279be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. If the forest is being created by coarsening a previous forest (see 7289be51f97SToby Isaac DMForestGetAdaptivityForest()), this limits the amount of coarsening. 7299be51f97SToby Isaac 7309be51f97SToby Isaac Not collective 7319be51f97SToby Isaac 7329be51f97SToby Isaac Input Parameter: 7339be51f97SToby Isaac . dm - the forest 7349be51f97SToby Isaac 7359be51f97SToby Isaac Output Parameter: 7369be51f97SToby Isaac . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7379be51f97SToby Isaac 7389be51f97SToby Isaac Level: intermediate 7399be51f97SToby Isaac 7409be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestGetMaximumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 7419be51f97SToby Isaac @*/ 742dd8e54a2SToby Isaac PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement) 743dd8e54a2SToby Isaac { 744dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 745dd8e54a2SToby Isaac 746dd8e54a2SToby Isaac PetscFunctionBegin; 747dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 748dd8e54a2SToby Isaac PetscValidIntPointer(minRefinement,2); 749dd8e54a2SToby Isaac *minRefinement = forest->minRefinement; 750dd8e54a2SToby Isaac PetscFunctionReturn(0); 751dd8e54a2SToby Isaac } 752dd8e54a2SToby Isaac 7539be51f97SToby Isaac /*@ 7549be51f97SToby Isaac DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base 7559be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. 7569be51f97SToby Isaac 7579be51f97SToby Isaac Logically collective on dm 7589be51f97SToby Isaac 7599be51f97SToby Isaac Input Parameters: 7609be51f97SToby Isaac + dm - the forest 7619be51f97SToby Isaac - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7629be51f97SToby Isaac 7639be51f97SToby Isaac Level: intermediate 7649be51f97SToby Isaac 7659be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM() 7669be51f97SToby Isaac @*/ 76756ba9f64SToby Isaac PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement) 76856ba9f64SToby Isaac { 76956ba9f64SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 77056ba9f64SToby Isaac 77156ba9f64SToby Isaac PetscFunctionBegin; 77256ba9f64SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 77328b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup"); 77456ba9f64SToby Isaac forest->initRefinement = initRefinement; 77556ba9f64SToby Isaac PetscFunctionReturn(0); 77656ba9f64SToby Isaac } 77756ba9f64SToby Isaac 7789be51f97SToby Isaac /*@ 7799be51f97SToby Isaac DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see 7809be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. 7819be51f97SToby Isaac 7829be51f97SToby Isaac Not collective 7839be51f97SToby Isaac 7849be51f97SToby Isaac Input Parameter: 7859be51f97SToby Isaac . dm - the forest 7869be51f97SToby Isaac 78701d2d390SJose E. Roman Output Parameter: 788bf2d5fbbSStefano Zampini . initRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7899be51f97SToby Isaac 7909be51f97SToby Isaac Level: intermediate 7919be51f97SToby Isaac 7929be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM() 7939be51f97SToby Isaac @*/ 79456ba9f64SToby Isaac PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement) 79556ba9f64SToby Isaac { 79656ba9f64SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 79756ba9f64SToby Isaac 79856ba9f64SToby Isaac PetscFunctionBegin; 79956ba9f64SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 80056ba9f64SToby Isaac PetscValidIntPointer(initRefinement,2); 80156ba9f64SToby Isaac *initRefinement = forest->initRefinement; 80256ba9f64SToby Isaac PetscFunctionReturn(0); 80356ba9f64SToby Isaac } 80456ba9f64SToby Isaac 8059be51f97SToby Isaac /*@ 8069be51f97SToby Isaac DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base 8079be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. If the forest is being created by refining a previous forest 8089be51f97SToby Isaac (see DMForestGetAdaptivityForest()), this limits the amount of refinement. 8099be51f97SToby Isaac 8109be51f97SToby Isaac Logically collective on dm 8119be51f97SToby Isaac 8129be51f97SToby Isaac Input Parameters: 8139be51f97SToby Isaac + dm - the forest 8149be51f97SToby Isaac - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 8159be51f97SToby Isaac 8169be51f97SToby Isaac Level: intermediate 8179be51f97SToby Isaac 8189be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityDM() 8199be51f97SToby Isaac @*/ 820c7eeac06SToby Isaac PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement) 821dd8e54a2SToby Isaac { 822dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 823dd8e54a2SToby Isaac 824dd8e54a2SToby Isaac PetscFunctionBegin; 825dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 82628b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup"); 827c7eeac06SToby Isaac forest->maxRefinement = maxRefinement; 828dd8e54a2SToby Isaac PetscFunctionReturn(0); 829dd8e54a2SToby Isaac } 830dd8e54a2SToby Isaac 8319be51f97SToby Isaac /*@ 8329be51f97SToby Isaac DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see 8339be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. If the forest is being created by refining a previous forest (see 8349be51f97SToby Isaac DMForestGetAdaptivityForest()), this limits the amount of refinement. 8359be51f97SToby Isaac 8369be51f97SToby Isaac Not collective 8379be51f97SToby Isaac 8389be51f97SToby Isaac Input Parameter: 8399be51f97SToby Isaac . dm - the forest 8409be51f97SToby Isaac 8419be51f97SToby Isaac Output Parameter: 8429be51f97SToby Isaac . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 8439be51f97SToby Isaac 8449be51f97SToby Isaac Level: intermediate 8459be51f97SToby Isaac 8469be51f97SToby Isaac .seealso: DMForestSetMaximumRefinement(), DMForestGetMinimumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 8479be51f97SToby Isaac @*/ 848c7eeac06SToby Isaac PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement) 849dd8e54a2SToby Isaac { 850dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 851dd8e54a2SToby Isaac 852dd8e54a2SToby Isaac PetscFunctionBegin; 853dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 854c7eeac06SToby Isaac PetscValidIntPointer(maxRefinement,2); 855c7eeac06SToby Isaac *maxRefinement = forest->maxRefinement; 856dd8e54a2SToby Isaac PetscFunctionReturn(0); 857dd8e54a2SToby Isaac } 858c7eeac06SToby Isaac 8599be51f97SToby Isaac /*@C 8609be51f97SToby Isaac DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes. 8619be51f97SToby Isaac Subtypes of DMForest may define their own strategies. Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree 8629be51f97SToby Isaac for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to 8639be51f97SToby Isaac specify refinement/coarsening. 8649be51f97SToby Isaac 8659be51f97SToby Isaac Logically collective on dm 8669be51f97SToby Isaac 8679be51f97SToby Isaac Input Parameters: 8689be51f97SToby Isaac + dm - the forest 8699be51f97SToby Isaac - adaptStrategy - default DMFORESTADAPTALL 8709be51f97SToby Isaac 8719be51f97SToby Isaac Level: advanced 8729be51f97SToby Isaac 8739be51f97SToby Isaac .seealso: DMForestGetAdaptivityStrategy() 8749be51f97SToby Isaac @*/ 875c7eeac06SToby Isaac PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy) 876c7eeac06SToby Isaac { 877c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 878c7eeac06SToby Isaac 879c7eeac06SToby Isaac PetscFunctionBegin; 880c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 8819566063dSJacob Faibussowitsch PetscCall(PetscFree(forest->adaptStrategy)); 8829566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy((const char*) adaptStrategy,(char**)&forest->adaptStrategy)); 883c7eeac06SToby Isaac PetscFunctionReturn(0); 884c7eeac06SToby Isaac } 885c7eeac06SToby Isaac 8869be51f97SToby Isaac /*@C 8879be51f97SToby Isaac DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes. Subtypes 8889be51f97SToby Isaac of DMForest may define their own strategies. Two default strategies are DMFORESTADAPTALL, which indicates that all 8899be51f97SToby Isaac processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only 8909be51f97SToby Isaac one process needs to specify refinement/coarsening. 8919be51f97SToby Isaac 8929be51f97SToby Isaac Not collective 8939be51f97SToby Isaac 8949be51f97SToby Isaac Input Parameter: 8959be51f97SToby Isaac . dm - the forest 8969be51f97SToby Isaac 8979be51f97SToby Isaac Output Parameter: 8989be51f97SToby Isaac . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL) 8999be51f97SToby Isaac 9009be51f97SToby Isaac Level: advanced 9019be51f97SToby Isaac 9029be51f97SToby Isaac .seealso: DMForestSetAdaptivityStrategy() 9039be51f97SToby Isaac @*/ 904c7eeac06SToby Isaac PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy) 905c7eeac06SToby Isaac { 906c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 907c7eeac06SToby Isaac 908c7eeac06SToby Isaac PetscFunctionBegin; 909c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 910c7eeac06SToby Isaac PetscValidPointer(adaptStrategy,2); 911c7eeac06SToby Isaac *adaptStrategy = forest->adaptStrategy; 912c7eeac06SToby Isaac PetscFunctionReturn(0); 913c7eeac06SToby Isaac } 914c7eeac06SToby Isaac 9152a133e43SToby Isaac /*@ 9162a133e43SToby Isaac DMForestGetAdaptivitySuccess - Return whether the requested adaptation (refinement, coarsening, repartitioning, 9172a133e43SToby Isaac etc.) was successful. PETSC_FALSE indicates that the post-adaptation forest is the same as the pre-adpatation 9182a133e43SToby Isaac forest. A requested adaptation may have been unsuccessful if, for example, the requested refinement would have 9192a133e43SToby Isaac exceeded the maximum refinement level. 9202a133e43SToby Isaac 9212a133e43SToby Isaac Collective on dm 9222a133e43SToby Isaac 9232a133e43SToby Isaac Input Parameter: 9242a133e43SToby Isaac 9252a133e43SToby Isaac . dm - the post-adaptation forest 9262a133e43SToby Isaac 9272a133e43SToby Isaac Output Parameter: 9282a133e43SToby Isaac 9292a133e43SToby Isaac . success - PETSC_TRUE if the post-adaptation forest is different from the pre-adaptation forest. 9302a133e43SToby Isaac 9312a133e43SToby Isaac Level: intermediate 9322a133e43SToby Isaac 9332a133e43SToby Isaac .see 9342a133e43SToby Isaac @*/ 9352a133e43SToby Isaac PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success) 9362a133e43SToby Isaac { 9372a133e43SToby Isaac DM_Forest *forest; 9382a133e43SToby Isaac 9392a133e43SToby Isaac PetscFunctionBegin; 9402a133e43SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 94128b400f6SJacob Faibussowitsch PetscCheck(dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() has not been called yet."); 9422a133e43SToby Isaac forest = (DM_Forest *) dm->data; 9439566063dSJacob Faibussowitsch PetscCall((forest->getadaptivitysuccess)(dm,success)); 9442a133e43SToby Isaac PetscFunctionReturn(0); 9452a133e43SToby Isaac } 9462a133e43SToby Isaac 947bf9b5d84SToby Isaac /*@ 948bf9b5d84SToby Isaac DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed 949bf9b5d84SToby Isaac relating the cells of the pre-adaptation forest to the post-adaptiation forest. After DMSetUp() is called, these transfer PetscSFs can be accessed with DMForestGetAdaptivitySF(). 950bf9b5d84SToby Isaac 951bf9b5d84SToby Isaac Logically collective on dm 952bf9b5d84SToby Isaac 953bf9b5d84SToby Isaac Input Parameters: 954bf9b5d84SToby Isaac + dm - the post-adaptation forest 955bf9b5d84SToby Isaac - computeSF - default PETSC_TRUE 956bf9b5d84SToby Isaac 957bf9b5d84SToby Isaac Level: advanced 958bf9b5d84SToby Isaac 959bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF() 960bf9b5d84SToby Isaac @*/ 961bf9b5d84SToby Isaac PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF) 962bf9b5d84SToby Isaac { 963bf9b5d84SToby Isaac DM_Forest *forest; 964bf9b5d84SToby Isaac 965bf9b5d84SToby Isaac PetscFunctionBegin; 966bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 96728b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called"); 968bf9b5d84SToby Isaac forest = (DM_Forest*) dm->data; 969bf9b5d84SToby Isaac forest->computeAdaptSF = computeSF; 970bf9b5d84SToby Isaac PetscFunctionReturn(0); 971bf9b5d84SToby Isaac } 972bf9b5d84SToby Isaac 9730eb7e1eaSToby Isaac PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time) 97480b27e07SToby Isaac { 97580b27e07SToby Isaac DM_Forest *forest; 97680b27e07SToby Isaac 97780b27e07SToby Isaac PetscFunctionBegin; 97880b27e07SToby Isaac PetscValidHeaderSpecific(dmIn ,DM_CLASSID ,1); 97980b27e07SToby Isaac PetscValidHeaderSpecific(vecIn ,VEC_CLASSID ,2); 98080b27e07SToby Isaac PetscValidHeaderSpecific(dmOut ,DM_CLASSID ,3); 98180b27e07SToby Isaac PetscValidHeaderSpecific(vecOut ,VEC_CLASSID ,4); 98280b27e07SToby Isaac forest = (DM_Forest *) dmIn->data; 98328b400f6SJacob Faibussowitsch PetscCheck(forest->transfervec,PetscObjectComm((PetscObject)dmIn),PETSC_ERR_SUP,"DMForestTransferVec() not implemented"); 9849566063dSJacob Faibussowitsch PetscCall((forest->transfervec)(dmIn,vecIn,dmOut,vecOut,useBCs,time)); 98580b27e07SToby Isaac PetscFunctionReturn(0); 98680b27e07SToby Isaac } 98780b27e07SToby Isaac 988ac34a06fSStefano Zampini PetscErrorCode DMForestTransferVecFromBase(DM dm, Vec vecIn, Vec vecOut) 989ac34a06fSStefano Zampini { 990ac34a06fSStefano Zampini DM_Forest *forest; 991ac34a06fSStefano Zampini 992ac34a06fSStefano Zampini PetscFunctionBegin; 993ac34a06fSStefano Zampini PetscValidHeaderSpecific(dm ,DM_CLASSID ,1); 994ac34a06fSStefano Zampini PetscValidHeaderSpecific(vecIn ,VEC_CLASSID ,2); 995ac34a06fSStefano Zampini PetscValidHeaderSpecific(vecOut ,VEC_CLASSID ,3); 996ac34a06fSStefano Zampini forest = (DM_Forest *) dm->data; 99728b400f6SJacob Faibussowitsch PetscCheck(forest->transfervecfrombase,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"DMForestTransferVecFromBase() not implemented"); 9989566063dSJacob Faibussowitsch PetscCall((forest->transfervecfrombase)(dm,vecIn,vecOut)); 999ac34a06fSStefano Zampini PetscFunctionReturn(0); 1000ac34a06fSStefano Zampini } 1001ac34a06fSStefano Zampini 1002bf9b5d84SToby Isaac /*@ 1003bf9b5d84SToby Isaac DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the 1004bf9b5d84SToby Isaac pre-adaptation forest to the post-adaptiation forest. After DMSetUp() is called, these transfer PetscSFs can be 1005bf9b5d84SToby Isaac accessed with DMForestGetAdaptivitySF(). 1006bf9b5d84SToby Isaac 1007bf9b5d84SToby Isaac Not collective 1008bf9b5d84SToby Isaac 1009bf9b5d84SToby Isaac Input Parameter: 1010bf9b5d84SToby Isaac . dm - the post-adaptation forest 1011bf9b5d84SToby Isaac 1012bf9b5d84SToby Isaac Output Parameter: 1013bf9b5d84SToby Isaac . computeSF - default PETSC_TRUE 1014bf9b5d84SToby Isaac 1015bf9b5d84SToby Isaac Level: advanced 1016bf9b5d84SToby Isaac 1017bf9b5d84SToby Isaac .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF() 1018bf9b5d84SToby Isaac @*/ 1019bf9b5d84SToby Isaac PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF) 1020bf9b5d84SToby Isaac { 1021bf9b5d84SToby Isaac DM_Forest *forest; 1022bf9b5d84SToby Isaac 1023bf9b5d84SToby Isaac PetscFunctionBegin; 1024bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1025bf9b5d84SToby Isaac forest = (DM_Forest*) dm->data; 1026bf9b5d84SToby Isaac *computeSF = forest->computeAdaptSF; 1027bf9b5d84SToby Isaac PetscFunctionReturn(0); 1028bf9b5d84SToby Isaac } 1029bf9b5d84SToby Isaac 1030bf9b5d84SToby Isaac /*@ 1031bf9b5d84SToby Isaac DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest. 1032bf9b5d84SToby Isaac Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be 1033bf9b5d84SToby Isaac some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa. Therefore there are two 1034bf9b5d84SToby Isaac PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates 1035bf9b5d84SToby Isaac pre-adaptation fine cells to post-adaptation coarse cells. 1036bf9b5d84SToby Isaac 1037bf9b5d84SToby Isaac Not collective 1038bf9b5d84SToby Isaac 1039bf9b5d84SToby Isaac Input Parameter: 1040bf9b5d84SToby Isaac dm - the post-adaptation forest 1041bf9b5d84SToby Isaac 1042bf9b5d84SToby Isaac Output Parameter: 10430f17b9e3SToby Isaac preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post- 10440f17b9e3SToby Isaac coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre- 1045bf9b5d84SToby Isaac 1046bf9b5d84SToby Isaac Level: advanced 1047bf9b5d84SToby Isaac 1048bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF() 1049bf9b5d84SToby Isaac @*/ 10500f17b9e3SToby Isaac PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine) 1051bf9b5d84SToby Isaac { 1052bf9b5d84SToby Isaac DM_Forest *forest; 1053bf9b5d84SToby Isaac 1054bf9b5d84SToby Isaac PetscFunctionBegin; 1055bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 10569566063dSJacob Faibussowitsch PetscCall(DMSetUp(dm)); 1057bf9b5d84SToby Isaac forest = (DM_Forest*) dm->data; 1058f885a11aSToby Isaac if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine; 1059f885a11aSToby Isaac if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine; 1060bf9b5d84SToby Isaac PetscFunctionReturn(0); 1061bf9b5d84SToby Isaac } 1062bf9b5d84SToby Isaac 10639be51f97SToby Isaac /*@ 10649be51f97SToby Isaac DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to 10659be51f97SToby Isaac indicate that the diameter of neighboring cells should differ by at most a factor of 2. Subtypes of DMForest may 10669be51f97SToby Isaac only support one particular choice of grading factor. 10679be51f97SToby Isaac 10689be51f97SToby Isaac Logically collective on dm 10699be51f97SToby Isaac 10709be51f97SToby Isaac Input Parameters: 10719be51f97SToby Isaac + dm - the forest 10729be51f97SToby Isaac - grade - the grading factor 10739be51f97SToby Isaac 10749be51f97SToby Isaac Level: advanced 10759be51f97SToby Isaac 10769be51f97SToby Isaac .seealso: DMForestGetGradeFactor() 10779be51f97SToby Isaac @*/ 1078c7eeac06SToby Isaac PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade) 1079c7eeac06SToby Isaac { 1080c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1081c7eeac06SToby Isaac 1082c7eeac06SToby Isaac PetscFunctionBegin; 1083c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 108428b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup"); 1085c7eeac06SToby Isaac forest->gradeFactor = grade; 1086c7eeac06SToby Isaac PetscFunctionReturn(0); 1087c7eeac06SToby Isaac } 1088c7eeac06SToby Isaac 10899be51f97SToby Isaac /*@ 10909be51f97SToby Isaac DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of 10919be51f97SToby Isaac neighboring cells should differ by at most a factor of 2. Subtypes of DMForest may only support one particular 10929be51f97SToby Isaac choice of grading factor. 10939be51f97SToby Isaac 10949be51f97SToby Isaac Not collective 10959be51f97SToby Isaac 10969be51f97SToby Isaac Input Parameter: 10979be51f97SToby Isaac . dm - the forest 10989be51f97SToby Isaac 10999be51f97SToby Isaac Output Parameter: 11009be51f97SToby Isaac . grade - the grading factor 11019be51f97SToby Isaac 11029be51f97SToby Isaac Level: advanced 11039be51f97SToby Isaac 11049be51f97SToby Isaac .seealso: DMForestSetGradeFactor() 11059be51f97SToby Isaac @*/ 1106c7eeac06SToby Isaac PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade) 1107c7eeac06SToby Isaac { 1108c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1109c7eeac06SToby Isaac 1110c7eeac06SToby Isaac PetscFunctionBegin; 1111c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1112c7eeac06SToby Isaac PetscValidIntPointer(grade,2); 1113c7eeac06SToby Isaac *grade = forest->gradeFactor; 1114c7eeac06SToby Isaac PetscFunctionReturn(0); 1115c7eeac06SToby Isaac } 1116c7eeac06SToby Isaac 11179be51f97SToby Isaac /*@ 11189be51f97SToby Isaac DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes 11199be51f97SToby Isaac the cell weight (see DMForestSetCellWeights()) when calculating partitions. The final weight of a cell will be 11209be51f97SToby Isaac (cellWeight) * (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on 11219be51f97SToby Isaac its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the 11229be51f97SToby Isaac computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement. 11239be51f97SToby Isaac 11249be51f97SToby Isaac Logically collective on dm 11259be51f97SToby Isaac 11269be51f97SToby Isaac Input Parameters: 11279be51f97SToby Isaac + dm - the forest 11289be51f97SToby Isaac - weightsFactors - default 1. 11299be51f97SToby Isaac 11309be51f97SToby Isaac Level: advanced 11319be51f97SToby Isaac 11329be51f97SToby Isaac .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights() 11339be51f97SToby Isaac @*/ 1134ef51cf95SToby Isaac PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor) 1135c7eeac06SToby Isaac { 1136c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1137c7eeac06SToby Isaac 1138c7eeac06SToby Isaac PetscFunctionBegin; 1139c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 114028b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup"); 1141c7eeac06SToby Isaac forest->weightsFactor = weightsFactor; 1142c7eeac06SToby Isaac PetscFunctionReturn(0); 1143c7eeac06SToby Isaac } 1144c7eeac06SToby Isaac 11459be51f97SToby Isaac /*@ 11469be51f97SToby Isaac DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see 11479be51f97SToby Isaac DMForestSetCellWeights()) when calculating partitions. The final weight of a cell will be (cellWeight) * 11489be51f97SToby Isaac (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on its level; a 11499be51f97SToby Isaac factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation 11509be51f97SToby Isaac associated with a cell is multiplied by a factor of 2 for each additional level of refinement. 11519be51f97SToby Isaac 11529be51f97SToby Isaac Not collective 11539be51f97SToby Isaac 11549be51f97SToby Isaac Input Parameter: 11559be51f97SToby Isaac . dm - the forest 11569be51f97SToby Isaac 11579be51f97SToby Isaac Output Parameter: 11589be51f97SToby Isaac . weightsFactors - default 1. 11599be51f97SToby Isaac 11609be51f97SToby Isaac Level: advanced 11619be51f97SToby Isaac 11629be51f97SToby Isaac .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights() 11639be51f97SToby Isaac @*/ 1164ef51cf95SToby Isaac PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor) 1165c7eeac06SToby Isaac { 1166c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1167c7eeac06SToby Isaac 1168c7eeac06SToby Isaac PetscFunctionBegin; 1169c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1170c7eeac06SToby Isaac PetscValidRealPointer(weightsFactor,2); 1171c7eeac06SToby Isaac *weightsFactor = forest->weightsFactor; 1172c7eeac06SToby Isaac PetscFunctionReturn(0); 1173c7eeac06SToby Isaac } 1174c7eeac06SToby Isaac 11759be51f97SToby Isaac /*@ 11769be51f97SToby Isaac DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process 11779be51f97SToby Isaac 11789be51f97SToby Isaac Not collective 11799be51f97SToby Isaac 11809be51f97SToby Isaac Input Parameter: 11819be51f97SToby Isaac . dm - the forest 11829be51f97SToby Isaac 11839be51f97SToby Isaac Output Parameters: 11849be51f97SToby Isaac + cStart - the first cell on this process 11859be51f97SToby Isaac - cEnd - one after the final cell on this process 11869be51f97SToby Isaac 11871a244344SSatish Balay Level: intermediate 11889be51f97SToby Isaac 11899be51f97SToby Isaac .seealso: DMForestGetCellSF() 11909be51f97SToby Isaac @*/ 1191c7eeac06SToby Isaac PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd) 1192c7eeac06SToby Isaac { 1193c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1194c7eeac06SToby Isaac 1195c7eeac06SToby Isaac PetscFunctionBegin; 1196c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1197c7eeac06SToby Isaac PetscValidIntPointer(cStart,2); 1198064a246eSJacob Faibussowitsch PetscValidIntPointer(cEnd,3); 1199c7eeac06SToby Isaac if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) { 12009566063dSJacob Faibussowitsch PetscCall(forest->createcellchart(dm,&forest->cStart,&forest->cEnd)); 1201c7eeac06SToby Isaac } 1202c7eeac06SToby Isaac *cStart = forest->cStart; 1203c7eeac06SToby Isaac *cEnd = forest->cEnd; 1204c7eeac06SToby Isaac PetscFunctionReturn(0); 1205c7eeac06SToby Isaac } 1206c7eeac06SToby Isaac 12079be51f97SToby Isaac /*@ 12089be51f97SToby Isaac DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes 12099be51f97SToby Isaac 12109be51f97SToby Isaac Not collective 12119be51f97SToby Isaac 12129be51f97SToby Isaac Input Parameter: 12139be51f97SToby Isaac . dm - the forest 12149be51f97SToby Isaac 12159be51f97SToby Isaac Output Parameter: 12169be51f97SToby Isaac . cellSF - the PetscSF 12179be51f97SToby Isaac 12181a244344SSatish Balay Level: intermediate 12199be51f97SToby Isaac 12209be51f97SToby Isaac .seealso: DMForestGetCellChart() 12219be51f97SToby Isaac @*/ 1222c7eeac06SToby Isaac PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF) 1223c7eeac06SToby Isaac { 1224c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1225c7eeac06SToby Isaac 1226c7eeac06SToby Isaac PetscFunctionBegin; 1227c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1228c7eeac06SToby Isaac PetscValidPointer(cellSF,2); 1229c7eeac06SToby Isaac if ((!forest->cellSF) && forest->createcellsf) { 12309566063dSJacob Faibussowitsch PetscCall(forest->createcellsf(dm,&forest->cellSF)); 1231c7eeac06SToby Isaac } 1232c7eeac06SToby Isaac *cellSF = forest->cellSF; 1233c7eeac06SToby Isaac PetscFunctionReturn(0); 1234c7eeac06SToby Isaac } 1235c7eeac06SToby Isaac 12369be51f97SToby Isaac /*@C 12379be51f97SToby Isaac DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see 12389be51f97SToby Isaac DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination). The 1239cd3c525cSToby Isaac interpretation of the label values is up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, 1240cd3c525cSToby Isaac DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have been reserved as choices that should be accepted by all subtypes. 12419be51f97SToby Isaac 12429be51f97SToby Isaac Logically collective on dm 12439be51f97SToby Isaac 12449be51f97SToby Isaac Input Parameters: 12459be51f97SToby Isaac - dm - the forest 1246a1b0c543SToby Isaac + adaptLabel - the label in the pre-adaptation forest 12479be51f97SToby Isaac 12489be51f97SToby Isaac Level: intermediate 12499be51f97SToby Isaac 12509be51f97SToby Isaac .seealso DMForestGetAdaptivityLabel() 12519be51f97SToby Isaac @*/ 1252a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityLabel(DM dm, DMLabel adaptLabel) 1253c7eeac06SToby Isaac { 1254c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1255c7eeac06SToby Isaac 1256c7eeac06SToby Isaac PetscFunctionBegin; 1257c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 12581fd50544SStefano Zampini if (adaptLabel) PetscValidHeaderSpecific(adaptLabel,DMLABEL_CLASSID,2); 12599566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)adaptLabel)); 12609566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&forest->adaptLabel)); 12611fd50544SStefano Zampini forest->adaptLabel = adaptLabel; 1262c7eeac06SToby Isaac PetscFunctionReturn(0); 1263c7eeac06SToby Isaac } 1264c7eeac06SToby Isaac 12659be51f97SToby Isaac /*@C 12669be51f97SToby Isaac DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that 12679be51f97SToby Isaac holds the adaptation flags (refinement, coarsening, or some combination). The interpretation of the label values is 1268cd3c525cSToby Isaac up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have 1269cd3c525cSToby Isaac been reserved as choices that should be accepted by all subtypes. 12709be51f97SToby Isaac 12719be51f97SToby Isaac Not collective 12729be51f97SToby Isaac 12739be51f97SToby Isaac Input Parameter: 12749be51f97SToby Isaac . dm - the forest 12759be51f97SToby Isaac 12769be51f97SToby Isaac Output Parameter: 12779be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest 12789be51f97SToby Isaac 12799be51f97SToby Isaac Level: intermediate 12809be51f97SToby Isaac 12819be51f97SToby Isaac .seealso DMForestSetAdaptivityLabel() 12829be51f97SToby Isaac @*/ 1283a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityLabel(DM dm, DMLabel *adaptLabel) 1284c7eeac06SToby Isaac { 1285c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1286c7eeac06SToby Isaac 1287c7eeac06SToby Isaac PetscFunctionBegin; 1288c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1289ba936b91SToby Isaac *adaptLabel = forest->adaptLabel; 1290c7eeac06SToby Isaac PetscFunctionReturn(0); 1291c7eeac06SToby Isaac } 1292c7eeac06SToby Isaac 12939be51f97SToby Isaac /*@ 12949be51f97SToby Isaac DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current 12959be51f97SToby Isaac process: weights are used to determine parallel partitioning. Partitions will be created so that each process's 12969be51f97SToby Isaac ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight 12979be51f97SToby Isaac of 1. 12989be51f97SToby Isaac 12999be51f97SToby Isaac Logically collective on dm 13009be51f97SToby Isaac 13019be51f97SToby Isaac Input Parameters: 13029be51f97SToby Isaac + dm - the forest 13039be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1. 13049be51f97SToby Isaac - copyMode - how weights should reference weights 13059be51f97SToby Isaac 13069be51f97SToby Isaac Level: advanced 13079be51f97SToby Isaac 13089be51f97SToby Isaac .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity() 13099be51f97SToby Isaac @*/ 1310c7eeac06SToby Isaac PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode) 1311c7eeac06SToby Isaac { 1312c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1313c7eeac06SToby Isaac PetscInt cStart, cEnd; 1314c7eeac06SToby Isaac 1315c7eeac06SToby Isaac PetscFunctionBegin; 1316c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 13179566063dSJacob Faibussowitsch PetscCall(DMForestGetCellChart(dm,&cStart,&cEnd)); 131863a3b9bcSJacob Faibussowitsch PetscCheck(cEnd >= cStart,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%" PetscInt_FMT ",%" PetscInt_FMT ") is not valid",cStart,cEnd); 1319c7eeac06SToby Isaac if (copyMode == PETSC_COPY_VALUES) { 1320c7eeac06SToby Isaac if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) { 13219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(cEnd-cStart,&forest->cellWeights)); 1322c7eeac06SToby Isaac } 13239566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(forest->cellWeights,weights,cEnd-cStart)); 1324c7eeac06SToby Isaac forest->cellWeightsCopyMode = PETSC_OWN_POINTER; 1325c7eeac06SToby Isaac PetscFunctionReturn(0); 1326c7eeac06SToby Isaac } 1327c7eeac06SToby Isaac if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) { 13289566063dSJacob Faibussowitsch PetscCall(PetscFree(forest->cellWeights)); 1329c7eeac06SToby Isaac } 1330c7eeac06SToby Isaac forest->cellWeights = weights; 1331c7eeac06SToby Isaac forest->cellWeightsCopyMode = copyMode; 1332c7eeac06SToby Isaac PetscFunctionReturn(0); 1333c7eeac06SToby Isaac } 1334c7eeac06SToby Isaac 13359be51f97SToby Isaac /*@ 13369be51f97SToby Isaac DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current 13379be51f97SToby Isaac process: weights are used to determine parallel partitioning. Partitions will be created so that each process's 13389be51f97SToby Isaac ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight 13399be51f97SToby Isaac of 1. 13409be51f97SToby Isaac 13419be51f97SToby Isaac Not collective 13429be51f97SToby Isaac 13439be51f97SToby Isaac Input Parameter: 13449be51f97SToby Isaac . dm - the forest 13459be51f97SToby Isaac 13469be51f97SToby Isaac Output Parameter: 13479be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1. 13489be51f97SToby Isaac 13499be51f97SToby Isaac Level: advanced 13509be51f97SToby Isaac 13519be51f97SToby Isaac .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity() 13529be51f97SToby Isaac @*/ 1353c7eeac06SToby Isaac PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights) 1354c7eeac06SToby Isaac { 1355c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1356c7eeac06SToby Isaac 1357c7eeac06SToby Isaac PetscFunctionBegin; 1358c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1359c7eeac06SToby Isaac PetscValidPointer(weights,2); 1360c7eeac06SToby Isaac *weights = forest->cellWeights; 1361c7eeac06SToby Isaac PetscFunctionReturn(0); 1362c7eeac06SToby Isaac } 1363c7eeac06SToby Isaac 13649be51f97SToby Isaac /*@ 13659be51f97SToby Isaac DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning 13669be51f97SToby Isaac a pre-adaptation forest (see DMForestGetAdaptivityForest()). After partitioning, the ratio of the weight of each 13679be51f97SToby Isaac process's cells to the process's capacity will be roughly equal for all processes. A capacity of 0 indicates that 13689be51f97SToby Isaac the current process should not have any cells after repartitioning. 13699be51f97SToby Isaac 13709be51f97SToby Isaac Logically Collective on dm 13719be51f97SToby Isaac 13729be51f97SToby Isaac Input parameters: 13739be51f97SToby Isaac + dm - the forest 13749be51f97SToby Isaac - capacity - this process's capacity 13759be51f97SToby Isaac 13769be51f97SToby Isaac Level: advanced 13779be51f97SToby Isaac 13789be51f97SToby Isaac .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor() 13799be51f97SToby Isaac @*/ 1380c7eeac06SToby Isaac PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity) 1381c7eeac06SToby Isaac { 1382c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1383c7eeac06SToby Isaac 1384c7eeac06SToby Isaac PetscFunctionBegin; 1385c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 138628b400f6SJacob Faibussowitsch PetscCheck(!dm->setupcalled,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup"); 138763a3b9bcSJacob Faibussowitsch PetscCheck(capacity >= 0.,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %g",(double)capacity); 1388c7eeac06SToby Isaac forest->weightCapacity = capacity; 1389c7eeac06SToby Isaac PetscFunctionReturn(0); 1390c7eeac06SToby Isaac } 1391c7eeac06SToby Isaac 13929be51f97SToby Isaac /*@ 13939be51f97SToby Isaac DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see 13949be51f97SToby Isaac DMForestGetAdaptivityForest()). After partitioning, the ratio of the weight of each process's cells to the 13959be51f97SToby Isaac process's capacity will be roughly equal for all processes. A capacity of 0 indicates that the current process 13969be51f97SToby Isaac should not have any cells after repartitioning. 13979be51f97SToby Isaac 13989be51f97SToby Isaac Not collective 13999be51f97SToby Isaac 14009be51f97SToby Isaac Input parameter: 14019be51f97SToby Isaac . dm - the forest 14029be51f97SToby Isaac 14039be51f97SToby Isaac Output parameter: 14049be51f97SToby Isaac . capacity - this process's capacity 14059be51f97SToby Isaac 14069be51f97SToby Isaac Level: advanced 14079be51f97SToby Isaac 14089be51f97SToby Isaac .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor() 14099be51f97SToby Isaac @*/ 1410c7eeac06SToby Isaac PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity) 1411c7eeac06SToby Isaac { 1412c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1413c7eeac06SToby Isaac 1414c7eeac06SToby Isaac PetscFunctionBegin; 1415c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1416c7eeac06SToby Isaac PetscValidRealPointer(capacity,2); 1417c7eeac06SToby Isaac *capacity = forest->weightCapacity; 1418c7eeac06SToby Isaac PetscFunctionReturn(0); 1419c7eeac06SToby Isaac } 1420c7eeac06SToby Isaac 14210709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm) 1422db4d5e8cSToby Isaac { 142356ba9f64SToby Isaac PetscBool flg, flg1, flg2, flg3, flg4; 1424dd8e54a2SToby Isaac DMForestTopology oldTopo; 1425c7eeac06SToby Isaac char stringBuffer[256]; 1426dd8e54a2SToby Isaac PetscViewer viewer; 1427dd8e54a2SToby Isaac PetscViewerFormat format; 142856ba9f64SToby Isaac PetscInt adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade; 1429c7eeac06SToby Isaac PetscReal weightsFactor; 1430c7eeac06SToby Isaac DMForestAdaptivityStrategy adaptStrategy; 1431db4d5e8cSToby Isaac 1432db4d5e8cSToby Isaac PetscFunctionBegin; 1433db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 14349566063dSJacob Faibussowitsch PetscCall(DMForestGetTopology(dm, &oldTopo)); 1435d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"DMForest Options"); 14369566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,sizeof(stringBuffer),&flg1)); 14379566063dSJacob Faibussowitsch PetscCall(PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2)); 14389566063dSJacob Faibussowitsch PetscCall(PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3)); 14399566063dSJacob Faibussowitsch PetscCall(PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4)); 1440*1dca8a05SBarry 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}"); 144156ba9f64SToby Isaac if (flg1) { 14429566063dSJacob Faibussowitsch PetscCall(DMForestSetTopology(dm,(DMForestTopology)stringBuffer)); 14439566063dSJacob Faibussowitsch PetscCall(DMForestSetBaseDM(dm,NULL)); 14449566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityForest(dm,NULL)); 144556ba9f64SToby Isaac } 144656ba9f64SToby Isaac if (flg2) { 1447dd8e54a2SToby Isaac DM base; 1448dd8e54a2SToby Isaac 14499566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm),&base)); 14509566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer,format)); 14519566063dSJacob Faibussowitsch PetscCall(DMLoad(base,viewer)); 14529566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 14539566063dSJacob Faibussowitsch PetscCall(DMForestSetBaseDM(dm,base)); 14549566063dSJacob Faibussowitsch PetscCall(DMDestroy(&base)); 14559566063dSJacob Faibussowitsch PetscCall(DMForestSetTopology(dm,NULL)); 14569566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityForest(dm,NULL)); 1457dd8e54a2SToby Isaac } 145856ba9f64SToby Isaac if (flg3) { 1459dd8e54a2SToby Isaac DM coarse; 1460dd8e54a2SToby Isaac 14619566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm),&coarse)); 14629566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer,format)); 14639566063dSJacob Faibussowitsch PetscCall(DMLoad(coarse,viewer)); 14649566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 14659566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityForest(dm,coarse)); 14669566063dSJacob Faibussowitsch PetscCall(DMDestroy(&coarse)); 14679566063dSJacob Faibussowitsch PetscCall(DMForestSetTopology(dm,NULL)); 14689566063dSJacob Faibussowitsch PetscCall(DMForestSetBaseDM(dm,NULL)); 1469dd8e54a2SToby Isaac } 147056ba9f64SToby Isaac if (flg4) { 1471dd8e54a2SToby Isaac DM fine; 1472dd8e54a2SToby Isaac 14739566063dSJacob Faibussowitsch PetscCall(DMCreate(PetscObjectComm((PetscObject)dm),&fine)); 14749566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer,format)); 14759566063dSJacob Faibussowitsch PetscCall(DMLoad(fine,viewer)); 14769566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 14779566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityForest(dm,fine)); 14789566063dSJacob Faibussowitsch PetscCall(DMDestroy(&fine)); 14799566063dSJacob Faibussowitsch PetscCall(DMForestSetTopology(dm,NULL)); 14809566063dSJacob Faibussowitsch PetscCall(DMForestSetBaseDM(dm,NULL)); 1481dd8e54a2SToby Isaac } 14829566063dSJacob Faibussowitsch PetscCall(DMForestGetAdjacencyDimension(dm,&adjDim)); 14839566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg,0)); 1484dd8e54a2SToby Isaac if (flg) { 14859566063dSJacob Faibussowitsch PetscCall(DMForestSetAdjacencyDimension(dm,adjDim)); 1486f885a11aSToby Isaac } else { 14879566063dSJacob Faibussowitsch PetscCall(DMForestGetAdjacencyCodimension(dm,&adjCodim)); 14889566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg,1)); 1489dd8e54a2SToby Isaac if (flg) { 14909566063dSJacob Faibussowitsch PetscCall(DMForestSetAdjacencyCodimension(dm,adjCodim)); 1491dd8e54a2SToby Isaac } 1492dd8e54a2SToby Isaac } 14939566063dSJacob Faibussowitsch PetscCall(DMForestGetPartitionOverlap(dm,&overlap)); 14949566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg,0)); 1495dd8e54a2SToby Isaac if (flg) { 14969566063dSJacob Faibussowitsch PetscCall(DMForestSetPartitionOverlap(dm,overlap)); 1497dd8e54a2SToby Isaac } 1498a6121fbdSMatthew G. Knepley #if 0 14999566063dSJacob 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)); 1500a6121fbdSMatthew G. Knepley if (flg) { 15019566063dSJacob Faibussowitsch PetscCall(DMForestSetMinimumRefinement(dm,minRefinement)); 15029566063dSJacob Faibussowitsch PetscCall(DMForestSetInitialRefinement(dm,minRefinement)); 1503a6121fbdSMatthew G. Knepley } 15049566063dSJacob 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)); 1505a6121fbdSMatthew G. Knepley if (flg) { 15069566063dSJacob Faibussowitsch PetscCall(DMForestSetMinimumRefinement(dm,0)); 15079566063dSJacob Faibussowitsch PetscCall(DMForestSetInitialRefinement(dm,initRefinement)); 1508a6121fbdSMatthew G. Knepley } 1509a6121fbdSMatthew G. Knepley #endif 15109566063dSJacob Faibussowitsch PetscCall(DMForestGetMinimumRefinement(dm,&minRefinement)); 15119566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg,0)); 1512dd8e54a2SToby Isaac if (flg) { 15139566063dSJacob Faibussowitsch PetscCall(DMForestSetMinimumRefinement(dm,minRefinement)); 1514db4d5e8cSToby Isaac } 15159566063dSJacob Faibussowitsch PetscCall(DMForestGetInitialRefinement(dm,&initRefinement)); 15169566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg,0)); 151756ba9f64SToby Isaac if (flg) { 15189566063dSJacob Faibussowitsch PetscCall(DMForestSetInitialRefinement(dm,initRefinement)); 151956ba9f64SToby Isaac } 15209566063dSJacob Faibussowitsch PetscCall(DMForestGetMaximumRefinement(dm,&maxRefinement)); 15219566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg,0)); 1522c7eeac06SToby Isaac if (flg) { 15239566063dSJacob Faibussowitsch PetscCall(DMForestSetMaximumRefinement(dm,maxRefinement)); 1524c7eeac06SToby Isaac } 15259566063dSJacob Faibussowitsch PetscCall(DMForestGetAdaptivityStrategy(dm,&adaptStrategy)); 15269566063dSJacob Faibussowitsch PetscCall(PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,sizeof(stringBuffer),&flg)); 1527c7eeac06SToby Isaac if (flg) { 15289566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer)); 1529c7eeac06SToby Isaac } 15309566063dSJacob Faibussowitsch PetscCall(DMForestGetGradeFactor(dm,&grade)); 15319566063dSJacob Faibussowitsch PetscCall(PetscOptionsBoundedInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg,0)); 1532c7eeac06SToby Isaac if (flg) { 15339566063dSJacob Faibussowitsch PetscCall(DMForestSetGradeFactor(dm,grade)); 1534c7eeac06SToby Isaac } 15359566063dSJacob Faibussowitsch PetscCall(DMForestGetCellWeightFactor(dm,&weightsFactor)); 15369566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg)); 1537c7eeac06SToby Isaac if (flg) { 15389566063dSJacob Faibussowitsch PetscCall(DMForestSetCellWeightFactor(dm,weightsFactor)); 1539c7eeac06SToby Isaac } 1540d0609cedSBarry Smith PetscOptionsHeadEnd(); 1541db4d5e8cSToby Isaac PetscFunctionReturn(0); 1542db4d5e8cSToby Isaac } 1543db4d5e8cSToby Isaac 1544276c5506SMatthew G. Knepley PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 1545d8984e3bSMatthew G. Knepley { 1546d8984e3bSMatthew G. Knepley PetscFunctionBegin; 15479566063dSJacob Faibussowitsch if (subdm) PetscCall(DMClone(dm, subdm)); 15489566063dSJacob Faibussowitsch PetscCall(DMCreateSectionSubDM(dm, numFields, fields, is, subdm)); 1549d8984e3bSMatthew G. Knepley PetscFunctionReturn(0); 1550d8984e3bSMatthew G. Knepley } 1551d8984e3bSMatthew G. Knepley 15525421bac9SToby Isaac PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined) 15535421bac9SToby Isaac { 15545421bac9SToby Isaac DMLabel refine; 15555421bac9SToby Isaac DM fineDM; 15565421bac9SToby Isaac 15575421bac9SToby Isaac PetscFunctionBegin; 15589566063dSJacob Faibussowitsch PetscCall(DMGetFineDM(dm,&fineDM)); 15595421bac9SToby Isaac if (fineDM) { 15609566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fineDM)); 15615421bac9SToby Isaac *dmRefined = fineDM; 15625421bac9SToby Isaac PetscFunctionReturn(0); 15635421bac9SToby Isaac } 15649566063dSJacob Faibussowitsch PetscCall(DMForestTemplate(dm,comm,dmRefined)); 15659566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm,"refine",&refine)); 15665421bac9SToby Isaac if (!refine) { 15679566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "refine",&refine)); 15689566063dSJacob Faibussowitsch PetscCall(DMLabelSetDefaultValue(refine,DM_ADAPT_REFINE)); 1569d67d17b1SMatthew G. Knepley } else { 15709566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject) refine)); 1571a1b0c543SToby Isaac } 15729566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityLabel(*dmRefined,refine)); 15739566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&refine)); 15745421bac9SToby Isaac PetscFunctionReturn(0); 15755421bac9SToby Isaac } 15765421bac9SToby Isaac 15775421bac9SToby Isaac PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened) 15785421bac9SToby Isaac { 15795421bac9SToby Isaac DMLabel coarsen; 15805421bac9SToby Isaac DM coarseDM; 15815421bac9SToby Isaac 15825421bac9SToby Isaac PetscFunctionBegin; 15834098eed7SToby Isaac { 15844098eed7SToby Isaac PetscMPIInt mpiComparison; 15854098eed7SToby Isaac MPI_Comm dmcomm = PetscObjectComm((PetscObject)dm); 15864098eed7SToby Isaac 15879566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_compare(comm, dmcomm, &mpiComparison)); 1588*1dca8a05SBarry Smith PetscCheck(mpiComparison == MPI_IDENT || mpiComparison == MPI_CONGRUENT,dmcomm,PETSC_ERR_SUP,"No support for different communicators yet"); 15894098eed7SToby Isaac } 15909566063dSJacob Faibussowitsch PetscCall(DMGetCoarseDM(dm,&coarseDM)); 15915421bac9SToby Isaac if (coarseDM) { 15929566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)coarseDM)); 15935421bac9SToby Isaac *dmCoarsened = coarseDM; 15945421bac9SToby Isaac PetscFunctionReturn(0); 15955421bac9SToby Isaac } 15969566063dSJacob Faibussowitsch PetscCall(DMForestTemplate(dm,comm,dmCoarsened)); 15979566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityPurpose(*dmCoarsened,DM_ADAPT_COARSEN)); 15989566063dSJacob Faibussowitsch PetscCall(DMGetLabel(dm,"coarsen",&coarsen)); 15995421bac9SToby Isaac if (!coarsen) { 16009566063dSJacob Faibussowitsch PetscCall(DMLabelCreate(PETSC_COMM_SELF, "coarsen",&coarsen)); 16019566063dSJacob Faibussowitsch PetscCall(DMLabelSetDefaultValue(coarsen,DM_ADAPT_COARSEN)); 1602a1b0c543SToby Isaac } else { 16039566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject) coarsen)); 16045421bac9SToby Isaac } 16059566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityLabel(*dmCoarsened,coarsen)); 16069566063dSJacob Faibussowitsch PetscCall(DMLabelDestroy(&coarsen)); 16075421bac9SToby Isaac PetscFunctionReturn(0); 16085421bac9SToby Isaac } 16095421bac9SToby Isaac 16109fe9e680SJoe Wallwork PetscErrorCode DMAdaptLabel_Forest(DM dm, PETSC_UNUSED Vec metric, DMLabel label, PETSC_UNUSED DMLabel rgLabel, DM *adaptedDM) 161109350103SToby Isaac { 161209350103SToby Isaac PetscBool success; 161309350103SToby Isaac 161409350103SToby Isaac PetscFunctionBegin; 16159566063dSJacob Faibussowitsch PetscCall(DMForestTemplate(dm,PetscObjectComm((PetscObject)dm),adaptedDM)); 16169566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityLabel(*adaptedDM,label)); 16179566063dSJacob Faibussowitsch PetscCall(DMSetUp(*adaptedDM)); 16189566063dSJacob Faibussowitsch PetscCall(DMForestGetAdaptivitySuccess(*adaptedDM,&success)); 161909350103SToby Isaac if (!success) { 16209566063dSJacob Faibussowitsch PetscCall(DMDestroy(adaptedDM)); 162109350103SToby Isaac *adaptedDM = NULL; 162209350103SToby Isaac } 162309350103SToby Isaac PetscFunctionReturn(0); 162409350103SToby Isaac } 162509350103SToby Isaac 1626d222f98bSToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm) 1627d222f98bSToby Isaac { 1628d222f98bSToby Isaac PetscFunctionBegin; 16299566063dSJacob Faibussowitsch PetscCall(PetscMemzero(dm->ops,sizeof(*(dm->ops)))); 1630d222f98bSToby Isaac 1631d222f98bSToby Isaac dm->ops->clone = DMClone_Forest; 1632d222f98bSToby Isaac dm->ops->setfromoptions = DMSetFromOptions_Forest; 1633d222f98bSToby Isaac dm->ops->destroy = DMDestroy_Forest; 1634d8984e3bSMatthew G. Knepley dm->ops->createsubdm = DMCreateSubDM_Forest; 16355421bac9SToby Isaac dm->ops->refine = DMRefine_Forest; 16365421bac9SToby Isaac dm->ops->coarsen = DMCoarsen_Forest; 1637d222f98bSToby Isaac PetscFunctionReturn(0); 1638d222f98bSToby Isaac } 1639d222f98bSToby Isaac 16409be51f97SToby Isaac /*MC 16419be51f97SToby Isaac 1642bae1f979SBarry Smith DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh. Forests usually have a base DM 1643bae1f979SBarry Smith (see DMForestGetBaseDM()), from which it is refined. The refinement and partitioning of forests is considered 1644bae1f979SBarry Smith immutable after DMSetUp() is called. To adapt a mesh, one should call DMForestTemplate() to create a new mesh that 1645bae1f979SBarry Smith will default to being identical to it, specify how that mesh should differ, and then calling DMSetUp() on the new 1646bae1f979SBarry Smith mesh. 1647bae1f979SBarry Smith 1648bae1f979SBarry Smith To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the 1649bae1f979SBarry Smith previous mesh whose values indicate which cells should be refined (DM_ADAPT_REFINE) or coarsened (DM_ADAPT_COARSEN) 1650bae1f979SBarry Smith and how (subtypes are free to allow additional values for things like anisotropic refinement). The label should be 1651bae1f979SBarry Smith given to the *new* mesh with DMForestSetAdaptivityLabel(). 16529be51f97SToby Isaac 16539be51f97SToby Isaac Level: advanced 16549be51f97SToby Isaac 16559be51f97SToby Isaac .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel() 16569be51f97SToby Isaac M*/ 16579be51f97SToby Isaac 1658db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm) 1659db4d5e8cSToby Isaac { 1660db4d5e8cSToby Isaac DM_Forest *forest; 1661db4d5e8cSToby Isaac 1662db4d5e8cSToby Isaac PetscFunctionBegin; 1663db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 16649566063dSJacob Faibussowitsch PetscCall(PetscNewLog(dm,&forest)); 1665db4d5e8cSToby Isaac dm->dim = 0; 1666db4d5e8cSToby Isaac dm->data = forest; 1667db4d5e8cSToby Isaac forest->refct = 1; 1668db4d5e8cSToby Isaac forest->data = NULL; 1669db4d5e8cSToby Isaac forest->topology = NULL; 1670cd3c525cSToby Isaac forest->adapt = NULL; 1671db4d5e8cSToby Isaac forest->base = NULL; 16726a87ffbfSToby Isaac forest->adaptPurpose = DM_ADAPT_DETERMINE; 1673db4d5e8cSToby Isaac forest->adjDim = PETSC_DEFAULT; 1674db4d5e8cSToby Isaac forest->overlap = PETSC_DEFAULT; 1675db4d5e8cSToby Isaac forest->minRefinement = PETSC_DEFAULT; 1676db4d5e8cSToby Isaac forest->maxRefinement = PETSC_DEFAULT; 167756ba9f64SToby Isaac forest->initRefinement = PETSC_DEFAULT; 1678c7eeac06SToby Isaac forest->cStart = PETSC_DETERMINE; 1679c7eeac06SToby Isaac forest->cEnd = PETSC_DETERMINE; 1680cd3c525cSToby Isaac forest->cellSF = NULL; 1681ebdf65a2SToby Isaac forest->adaptLabel = NULL; 1682db4d5e8cSToby Isaac forest->gradeFactor = 2; 1683db4d5e8cSToby Isaac forest->cellWeights = NULL; 1684db4d5e8cSToby Isaac forest->cellWeightsCopyMode = PETSC_USE_POINTER; 1685db4d5e8cSToby Isaac forest->weightsFactor = 1.; 1686db4d5e8cSToby Isaac forest->weightCapacity = 1.; 16879566063dSJacob Faibussowitsch PetscCall(DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL)); 16889566063dSJacob Faibussowitsch PetscCall(DMInitialize_Forest(dm)); 1689db4d5e8cSToby Isaac PetscFunctionReturn(0); 1690db4d5e8cSToby Isaac } 1691