19be51f97SToby Isaac #include <petsc/private/dmforestimpl.h> /*I "petscdmforest.h" I*/ 29be51f97SToby Isaac #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 3ef19d27cSToby Isaac #include <petscsf.h> 4db4d5e8cSToby Isaac 527d4645fSToby Isaac PetscBool DMForestPackageInitialized = PETSC_FALSE; 627d4645fSToby Isaac 727d4645fSToby Isaac typedef struct _DMForestTypeLink *DMForestTypeLink; 827d4645fSToby Isaac 927d4645fSToby Isaac struct _DMForestTypeLink 1027d4645fSToby Isaac { 1127d4645fSToby Isaac char *name; 1227d4645fSToby Isaac DMForestTypeLink next; 1327d4645fSToby Isaac }; 1427d4645fSToby Isaac 1527d4645fSToby Isaac DMForestTypeLink DMForestTypeList; 1627d4645fSToby Isaac 1727d4645fSToby Isaac #undef __FUNCT__ 1827d4645fSToby Isaac #define __FUNCT__ "DMForestPackageFinalize" 1927d4645fSToby Isaac static PetscErrorCode DMForestPackageFinalize(void) 2027d4645fSToby Isaac { 2127d4645fSToby Isaac DMForestTypeLink oldLink, link = DMForestTypeList; 2227d4645fSToby Isaac PetscErrorCode ierr; 2327d4645fSToby Isaac 2427d4645fSToby Isaac PetscFunctionBegin; 2527d4645fSToby Isaac while (link) { 2627d4645fSToby Isaac oldLink = link; 2727d4645fSToby Isaac ierr = PetscFree(oldLink->name); 2827d4645fSToby Isaac link = oldLink->next; 2927d4645fSToby Isaac ierr = PetscFree(oldLink);CHKERRQ(ierr); 3027d4645fSToby Isaac } 3127d4645fSToby Isaac PetscFunctionReturn(0); 3227d4645fSToby Isaac } 3327d4645fSToby Isaac 3427d4645fSToby Isaac #undef __FUNCT__ 3527d4645fSToby Isaac #define __FUNCT__ "DMForestPackageInitialize" 3627d4645fSToby Isaac static PetscErrorCode DMForestPackageInitialize(void) 3727d4645fSToby Isaac { 3827d4645fSToby Isaac PetscErrorCode ierr; 3927d4645fSToby Isaac 4027d4645fSToby Isaac PetscFunctionBegin; 4127d4645fSToby Isaac if (DMForestPackageInitialized) PetscFunctionReturn(0); 4227d4645fSToby Isaac DMForestPackageInitialized = PETSC_TRUE; 4327d4645fSToby Isaac ierr = DMForestRegisterType(DMFOREST);CHKERRQ(ierr); 4427d4645fSToby Isaac ierr = PetscRegisterFinalize(DMForestPackageFinalize);CHKERRQ(ierr); 4527d4645fSToby Isaac PetscFunctionReturn(0); 4627d4645fSToby Isaac } 4727d4645fSToby Isaac 4827d4645fSToby Isaac #undef __FUNCT__ 4927d4645fSToby Isaac #define __FUNCT__ "DMForestRegisterType" 509be51f97SToby Isaac /*@C 519be51f97SToby Isaac DMForestRegisterType - Registers a DMType as a subtype of DMFOREST (so that DMIsForest() will be correct) 529be51f97SToby Isaac 539be51f97SToby Isaac Not Collective 549be51f97SToby Isaac 559be51f97SToby Isaac Input parameter: 569be51f97SToby Isaac . name - the name of the type 579be51f97SToby Isaac 589be51f97SToby Isaac Level: advanced 599be51f97SToby Isaac 609be51f97SToby Isaac .seealso: DMFOREST, DMIsForest() 619be51f97SToby Isaac @*/ 6227d4645fSToby Isaac PetscErrorCode DMForestRegisterType(DMType name) 6327d4645fSToby Isaac { 6427d4645fSToby Isaac DMForestTypeLink link; 6527d4645fSToby Isaac PetscErrorCode ierr; 6627d4645fSToby Isaac 6727d4645fSToby Isaac PetscFunctionBegin; 6827d4645fSToby Isaac ierr = DMForestPackageInitialize();CHKERRQ(ierr); 6927d4645fSToby Isaac ierr = PetscNew(&link);CHKERRQ(ierr); 7027d4645fSToby Isaac ierr = PetscStrallocpy(name,&link->name);CHKERRQ(ierr); 7127d4645fSToby Isaac link->next = DMForestTypeList; 7227d4645fSToby Isaac DMForestTypeList = link; 7327d4645fSToby Isaac PetscFunctionReturn(0); 7427d4645fSToby Isaac } 7527d4645fSToby Isaac 7627d4645fSToby Isaac #undef __FUNCT__ 7727d4645fSToby Isaac #define __FUNCT__ "DMIsForest" 789be51f97SToby Isaac /*@ 799be51f97SToby Isaac DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes 809be51f97SToby Isaac 819be51f97SToby Isaac Not Collective 829be51f97SToby Isaac 839be51f97SToby Isaac Input parameter: 849be51f97SToby Isaac . dm - the DM object 859be51f97SToby Isaac 869be51f97SToby Isaac Output parameter: 879be51f97SToby Isaac . isForest - whether dm is a subtype of DMFOREST 889be51f97SToby Isaac 899be51f97SToby Isaac Level: intermediate 909be51f97SToby Isaac 919be51f97SToby Isaac .seealso: DMFOREST, DMForestRegisterType() 929be51f97SToby Isaac @*/ 9327d4645fSToby Isaac PetscErrorCode DMIsForest(DM dm, PetscBool *isForest) 9427d4645fSToby Isaac { 9527d4645fSToby Isaac DMForestTypeLink link = DMForestTypeList; 9627d4645fSToby Isaac PetscErrorCode ierr; 9727d4645fSToby Isaac 9827d4645fSToby Isaac PetscFunctionBegin; 9927d4645fSToby Isaac while (link) { 10027d4645fSToby Isaac PetscBool sameType; 10127d4645fSToby Isaac ierr = PetscObjectTypeCompare((PetscObject)dm,link->name,&sameType);CHKERRQ(ierr); 10227d4645fSToby Isaac if (sameType) { 10327d4645fSToby Isaac *isForest = PETSC_TRUE; 10427d4645fSToby Isaac PetscFunctionReturn(0); 10527d4645fSToby Isaac } 10627d4645fSToby Isaac link = link->next; 10727d4645fSToby Isaac } 10827d4645fSToby Isaac *isForest = PETSC_FALSE; 10927d4645fSToby Isaac PetscFunctionReturn(0); 11027d4645fSToby Isaac } 11127d4645fSToby Isaac 112db4d5e8cSToby Isaac #undef __FUNCT__ 113a0452a8eSToby Isaac #define __FUNCT__ "DMForestTemplate" 1149be51f97SToby Isaac /*@ 1159be51f97SToby Isaac DMForestTemplate - Create a new DM that will be adapted from a source DM. The new DM reproduces the configuration 1169be51f97SToby 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 1179be51f97SToby Isaac (by, e.g., refinement or repartitioning). The source DM is also set as the adaptivity source DM of the new DM (see 1189be51f97SToby Isaac DMForestSetAdaptivityForest()). 1199be51f97SToby Isaac 1209be51f97SToby Isaac Collective on dm 1219be51f97SToby Isaac 1229be51f97SToby Isaac Input Parameters: 1239be51f97SToby Isaac + dm - the source DM object 1249be51f97SToby Isaac - comm - the communicator for the new DM (this communicator is currently ignored, but is present so that DMForestTemplate() can be used within DMCoarsen()) 1259be51f97SToby Isaac 1269be51f97SToby Isaac Output Parameter: 1279be51f97SToby Isaac . tdm - the new DM object 1289be51f97SToby Isaac 1299be51f97SToby Isaac Level: intermediate 1309be51f97SToby Isaac 1319be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest() 1329be51f97SToby Isaac @*/ 1339be51f97SToby Isaac PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm) 134a0452a8eSToby Isaac { 135a0452a8eSToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 13620e8089bSToby Isaac DMType type; 137a0452a8eSToby Isaac DM base; 138a0452a8eSToby Isaac DMForestTopology topology; 139a0452a8eSToby Isaac PetscInt dim, overlap, ref, factor; 140a0452a8eSToby Isaac DMForestAdaptivityStrategy strat; 141795844e7SToby Isaac PetscDS ds; 142795844e7SToby Isaac void *ctx; 143a0452a8eSToby Isaac PetscErrorCode ierr; 144a0452a8eSToby Isaac 145a0452a8eSToby Isaac PetscFunctionBegin; 146a0452a8eSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 14720e8089bSToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),tdm);CHKERRQ(ierr); 14820e8089bSToby Isaac ierr = DMGetType(dm,&type);CHKERRQ(ierr); 14920e8089bSToby Isaac ierr = DMSetType(*tdm,type);CHKERRQ(ierr); 150a0452a8eSToby Isaac ierr = DMForestGetBaseDM(dm,&base);CHKERRQ(ierr); 15120e8089bSToby Isaac ierr = DMForestSetBaseDM(*tdm,base);CHKERRQ(ierr); 152a0452a8eSToby Isaac ierr = DMForestGetTopology(dm,&topology);CHKERRQ(ierr); 15320e8089bSToby Isaac ierr = DMForestSetTopology(*tdm,topology);CHKERRQ(ierr); 154a0452a8eSToby Isaac ierr = DMForestGetAdjacencyDimension(dm,&dim);CHKERRQ(ierr); 15520e8089bSToby Isaac ierr = DMForestSetAdjacencyDimension(*tdm,dim);CHKERRQ(ierr); 156a0452a8eSToby Isaac ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr); 15720e8089bSToby Isaac ierr = DMForestSetPartitionOverlap(*tdm,overlap);CHKERRQ(ierr); 158a0452a8eSToby Isaac ierr = DMForestGetMinimumRefinement(dm,&ref);CHKERRQ(ierr); 15920e8089bSToby Isaac ierr = DMForestSetMinimumRefinement(*tdm,ref);CHKERRQ(ierr); 160a0452a8eSToby Isaac ierr = DMForestGetMaximumRefinement(dm,&ref);CHKERRQ(ierr); 16120e8089bSToby Isaac ierr = DMForestSetMaximumRefinement(*tdm,ref);CHKERRQ(ierr); 162a0452a8eSToby Isaac ierr = DMForestGetAdaptivityStrategy(dm,&strat);CHKERRQ(ierr); 16320e8089bSToby Isaac ierr = DMForestSetAdaptivityStrategy(*tdm,strat);CHKERRQ(ierr); 164a0452a8eSToby Isaac ierr = DMForestGetGradeFactor(dm,&factor);CHKERRQ(ierr); 16520e8089bSToby Isaac ierr = DMForestSetGradeFactor(*tdm,factor);CHKERRQ(ierr); 166a0452a8eSToby Isaac if (forest->ftemplate) { 16720e8089bSToby Isaac ierr = (forest->ftemplate) (dm, *tdm);CHKERRQ(ierr); 168a0452a8eSToby Isaac } 16920e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(*tdm,dm);CHKERRQ(ierr); 170795844e7SToby Isaac ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 171795844e7SToby Isaac ierr = DMSetDS(*tdm,ds);CHKERRQ(ierr); 172795844e7SToby Isaac ierr = DMGetApplicationContext(dm,&ctx);CHKERRQ(ierr); 173795844e7SToby Isaac ierr = DMSetApplicationContext(*tdm,&ctx);CHKERRQ(ierr); 174795844e7SToby Isaac if (dm->maxCell) { 175795844e7SToby Isaac const PetscReal *maxCell, *L; 176795844e7SToby Isaac const DMBoundaryType *bd; 177795844e7SToby Isaac 178795844e7SToby Isaac ierr = DMGetPeriodicity(dm,&maxCell,&L,&bd);CHKERRQ(ierr); 179795844e7SToby Isaac ierr = DMSetPeriodicity(*tdm,maxCell,L,bd);CHKERRQ(ierr); 180795844e7SToby Isaac } 181a0452a8eSToby Isaac PetscFunctionReturn(0); 182a0452a8eSToby Isaac } 183a0452a8eSToby Isaac 18401d9d024SToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm); 18501d9d024SToby Isaac 186a0452a8eSToby Isaac #undef __FUNCT__ 187db4d5e8cSToby Isaac #define __FUNCT__ "DMClone_Forest" 188db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm) 189db4d5e8cSToby Isaac { 190db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 191db4d5e8cSToby Isaac const char *type; 192db4d5e8cSToby Isaac PetscErrorCode ierr; 193db4d5e8cSToby Isaac 194db4d5e8cSToby Isaac PetscFunctionBegin; 195db4d5e8cSToby Isaac forest->refct++; 196db4d5e8cSToby Isaac (*newdm)->data = forest; 197db4d5e8cSToby Isaac ierr = PetscObjectGetType((PetscObject) dm, &type);CHKERRQ(ierr); 198db4d5e8cSToby Isaac ierr = PetscObjectChangeTypeName((PetscObject) *newdm, type);CHKERRQ(ierr); 19901d9d024SToby Isaac ierr = DMInitialize_Forest(*newdm);CHKERRQ(ierr); 200db4d5e8cSToby Isaac PetscFunctionReturn(0); 201db4d5e8cSToby Isaac } 202db4d5e8cSToby Isaac 203db4d5e8cSToby Isaac #undef __FUNCT__ 204db4d5e8cSToby Isaac #define __FUNCT__ "DMDestroy_Forest" 205d222f98bSToby Isaac static PetscErrorCode DMDestroy_Forest(DM dm) 206db4d5e8cSToby Isaac { 207db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 208db4d5e8cSToby Isaac PetscErrorCode ierr; 209db4d5e8cSToby Isaac 210db4d5e8cSToby Isaac PetscFunctionBegin; 211db4d5e8cSToby Isaac if (--forest->refct > 0) PetscFunctionReturn(0); 212d222f98bSToby Isaac if (forest->destroy) {ierr = forest->destroy(dm);CHKERRQ(ierr);} 213db4d5e8cSToby Isaac ierr = PetscSFDestroy(&forest->cellSF);CHKERRQ(ierr); 214*bf9b5d84SToby Isaac ierr = PetscSFDestroy(&forest->adaptSFCoarseToFine);CHKERRQ(ierr); 215*bf9b5d84SToby Isaac ierr = PetscSFDestroy(&forest->adaptSFFineToCoarse);CHKERRQ(ierr); 216ebdf65a2SToby Isaac ierr = PetscFree(forest->adaptLabel);CHKERRQ(ierr); 2179a81d013SToby Isaac ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr); 21856ba9f64SToby Isaac ierr = DMDestroy(&forest->base);CHKERRQ(ierr); 219ba936b91SToby Isaac ierr = DMDestroy(&forest->adapt);CHKERRQ(ierr); 22030f902e7SToby Isaac ierr = PetscFree(forest->topology);CHKERRQ(ierr); 22130f902e7SToby Isaac ierr = PetscFree(forest);CHKERRQ(ierr); 222db4d5e8cSToby Isaac PetscFunctionReturn(0); 223db4d5e8cSToby Isaac } 224db4d5e8cSToby Isaac 225db4d5e8cSToby Isaac #undef __FUNCT__ 226dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetTopology" 2279be51f97SToby Isaac /*@C 2289be51f97SToby Isaac DMForestSetTopology - Set the topology of a DMForest during the pre-setup phase. The topology is a string (e.g. 2299be51f97SToby Isaac "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest durint 2309be51f97SToby Isaac DMSetUp(). 2319be51f97SToby Isaac 2329be51f97SToby Isaac Logically collective on dm 2339be51f97SToby Isaac 2349be51f97SToby Isaac Input parameters: 2359be51f97SToby Isaac + dm - the forest 2369be51f97SToby Isaac - topology - the topology of the forest 2379be51f97SToby Isaac 2389be51f97SToby Isaac Level: intermediate 2399be51f97SToby Isaac 2409be51f97SToby Isaac .seealso(): DMForestGetTopology(), DMForestSetBaseDM() 2419be51f97SToby Isaac @*/ 242dd8e54a2SToby Isaac PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology) 243db4d5e8cSToby Isaac { 244db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 245db4d5e8cSToby Isaac PetscErrorCode ierr; 246db4d5e8cSToby Isaac 247db4d5e8cSToby Isaac PetscFunctionBegin; 248db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 249ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the topology after setup"); 250dd8e54a2SToby Isaac ierr = PetscFree(forest->topology);CHKERRQ(ierr); 251dd8e54a2SToby Isaac ierr = PetscStrallocpy((const char *)topology,(char **) &forest->topology);CHKERRQ(ierr); 252db4d5e8cSToby Isaac PetscFunctionReturn(0); 253db4d5e8cSToby Isaac } 254db4d5e8cSToby Isaac 255db4d5e8cSToby Isaac #undef __FUNCT__ 256dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetTopology" 2579be51f97SToby Isaac /*@C 2589be51f97SToby Isaac DMForestGetTopology - Get a string describing the topology of a DMForest. 2599be51f97SToby Isaac 2609be51f97SToby Isaac Not collective 2619be51f97SToby Isaac 2629be51f97SToby Isaac Input parameter: 2639be51f97SToby Isaac . dm - the forest 2649be51f97SToby Isaac 2659be51f97SToby Isaac Output parameter: 2669be51f97SToby Isaac . topology - the topology of the forest (e.g., 'cube', 'shell') 2679be51f97SToby Isaac 2689be51f97SToby Isaac Level: intermediate 2699be51f97SToby Isaac 2709be51f97SToby Isaac .seealso: DMForestSetTopology() 2719be51f97SToby Isaac @*/ 272dd8e54a2SToby Isaac PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology) 273dd8e54a2SToby Isaac { 274dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 275dd8e54a2SToby Isaac 276dd8e54a2SToby Isaac PetscFunctionBegin; 277dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 278dd8e54a2SToby Isaac PetscValidPointer(topology,2); 279dd8e54a2SToby Isaac *topology = forest->topology; 280dd8e54a2SToby Isaac PetscFunctionReturn(0); 281dd8e54a2SToby Isaac } 282dd8e54a2SToby Isaac 283dd8e54a2SToby Isaac #undef __FUNCT__ 284dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetBaseDM" 2859be51f97SToby Isaac /*@ 2869be51f97SToby Isaac DMForestSetBaseDM - During the pre-setup phase, set the DM that defines the base mesh of a DMForest forest. The 2879be51f97SToby Isaac forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its 2889be51f97SToby Isaac base. In general, two forest must share a bse to be comparable, to do things like construct interpolators. 2899be51f97SToby Isaac 2909be51f97SToby Isaac Logically collective on dm 2919be51f97SToby Isaac 2929be51f97SToby Isaac Input Parameters: 2939be51f97SToby Isaac + dm - the forest 2949be51f97SToby Isaac - base - the base DM of the forest 2959be51f97SToby Isaac 2969be51f97SToby Isaac Level: intermediate 2979be51f97SToby Isaac 2989be51f97SToby Isaac .seealso(): DMForestGetBaseDM() 2999be51f97SToby Isaac @*/ 300dd8e54a2SToby Isaac PetscErrorCode DMForestSetBaseDM(DM dm, DM base) 301dd8e54a2SToby Isaac { 302dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 303dd8e54a2SToby Isaac PetscInt dim, dimEmbed; 304dd8e54a2SToby Isaac PetscErrorCode ierr; 305dd8e54a2SToby Isaac 306dd8e54a2SToby Isaac PetscFunctionBegin; 307dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 308ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the base after setup"); 309dd8e54a2SToby Isaac ierr = PetscObjectReference((PetscObject)base);CHKERRQ(ierr); 310dd8e54a2SToby Isaac ierr = DMDestroy(&forest->base);CHKERRQ(ierr); 311dd8e54a2SToby Isaac forest->base = base; 312a0452a8eSToby Isaac if (base) { 313a0452a8eSToby Isaac PetscValidHeaderSpecific(base, DM_CLASSID, 2); 314dd8e54a2SToby Isaac ierr = DMGetDimension(base,&dim);CHKERRQ(ierr); 315dd8e54a2SToby Isaac ierr = DMSetDimension(dm,dim);CHKERRQ(ierr); 316dd8e54a2SToby Isaac ierr = DMGetCoordinateDim(base,&dimEmbed);CHKERRQ(ierr); 317dd8e54a2SToby Isaac ierr = DMSetCoordinateDim(dm,dimEmbed);CHKERRQ(ierr); 318a0452a8eSToby Isaac } 319dd8e54a2SToby Isaac PetscFunctionReturn(0); 320dd8e54a2SToby Isaac } 321dd8e54a2SToby Isaac 322dd8e54a2SToby Isaac #undef __FUNCT__ 323dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetBaseDM" 3249be51f97SToby Isaac /*@ 3259be51f97SToby Isaac DMForestGetBaseDM - Get the base DM of a DMForest forest. The forest will be hierarchically refined from the base, 3269be51f97SToby Isaac and all refinements/coarsenings of the forest will share its base. In general, two forest must share a bse to be 3279be51f97SToby Isaac comparable, to do things like construct interpolators. 3289be51f97SToby Isaac 3299be51f97SToby Isaac Not collective 3309be51f97SToby Isaac 3319be51f97SToby Isaac Input Parameter: 3329be51f97SToby Isaac . dm - the forest 3339be51f97SToby Isaac 3349be51f97SToby Isaac Output Parameter: 3359be51f97SToby Isaac . base - the base DM of the forest 3369be51f97SToby Isaac 3379be51f97SToby Isaac Level: intermediate 3389be51f97SToby Isaac 3399be51f97SToby Isaac .seealso(); DMForestSetBaseDM() 3409be51f97SToby Isaac @*/ 341dd8e54a2SToby Isaac PetscErrorCode DMForestGetBaseDM(DM dm, DM *base) 342dd8e54a2SToby Isaac { 343dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 344dd8e54a2SToby Isaac 345dd8e54a2SToby Isaac PetscFunctionBegin; 346dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 347dd8e54a2SToby Isaac PetscValidPointer(base, 2); 348dd8e54a2SToby Isaac *base = forest->base; 349dd8e54a2SToby Isaac PetscFunctionReturn(0); 350dd8e54a2SToby Isaac } 351dd8e54a2SToby Isaac 352dd8e54a2SToby Isaac #undef __FUNCT__ 353ba936b91SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityForest" 3549be51f97SToby Isaac /*@ 3559be51f97SToby Isaac DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be 3569be51f97SToby Isaac adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) im DMSetUp(). Usually not needed 3579be51f97SToby Isaac by users directly: DMForestTemplate() constructs a new forest to be adapted from an old forest and calls this 3589be51f97SToby Isaac routine. 3599be51f97SToby Isaac 3609be51f97SToby Isaac Logically collective on dm 3619be51f97SToby Isaac 3629be51f97SToby Isaac Input Parameter: 3639be51f97SToby Isaac + dm - the new forest, which will be constructed from adapt 3649be51f97SToby Isaac - adapt - the old forest 3659be51f97SToby Isaac 3669be51f97SToby Isaac Level: intermediate 3679be51f97SToby Isaac 3689be51f97SToby Isaac .seealso: DMForestGetAdaptivityForest(), DMForestSetAdaptivityPurpose() 3699be51f97SToby Isaac @*/ 370ba936b91SToby Isaac PetscErrorCode DMForestSetAdaptivityForest(DM dm,DM adapt) 371dd8e54a2SToby Isaac { 372ba936b91SToby Isaac DM_Forest *forest; 373dd8e54a2SToby Isaac PetscErrorCode ierr; 374dd8e54a2SToby Isaac 375dd8e54a2SToby Isaac PetscFunctionBegin; 376dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 377ba936b91SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 378ba936b91SToby Isaac forest = (DM_Forest *) dm->data; 379ba936b91SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup"); 38026d9498aSToby Isaac switch (forest->adaptPurpose) { 38126d9498aSToby Isaac case DM_FOREST_KEEP: 382ba936b91SToby Isaac ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr); 383ba936b91SToby Isaac ierr = DMDestroy(&(forest->adapt));CHKERRQ(ierr); 384ba936b91SToby Isaac forest->adapt = adapt; 38526d9498aSToby Isaac break; 38626d9498aSToby Isaac case DM_FOREST_REFINE: 38726d9498aSToby Isaac ierr = DMSetCoarseDM(dm,adapt);CHKERRQ(ierr); 38826d9498aSToby Isaac break; 38926d9498aSToby Isaac case DM_FOREST_COARSEN: 39026d9498aSToby Isaac ierr = DMSetFineDM(dm,adapt);CHKERRQ(ierr); 39126d9498aSToby Isaac break; 39226d9498aSToby Isaac default: 39326d9498aSToby Isaac SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose"); 39426d9498aSToby Isaac } 395dd8e54a2SToby Isaac PetscFunctionReturn(0); 396dd8e54a2SToby Isaac } 397dd8e54a2SToby Isaac 398dd8e54a2SToby Isaac #undef __FUNCT__ 399ba936b91SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityForest" 4009be51f97SToby Isaac /*@ 4019be51f97SToby Isaac DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted. 4029be51f97SToby Isaac 4039be51f97SToby Isaac Not collective 4049be51f97SToby Isaac 4059be51f97SToby Isaac Input Parameter: 4069be51f97SToby Isaac . dm - the forest 4079be51f97SToby Isaac 4089be51f97SToby Isaac Output Parameter: 4099be51f97SToby Isaac . adapt - the forest from which dm is/was adapted 4109be51f97SToby Isaac 4119be51f97SToby Isaac Level: intermediate 4129be51f97SToby Isaac 4139be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose() 4149be51f97SToby Isaac @*/ 415ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt) 416dd8e54a2SToby Isaac { 417ba936b91SToby Isaac DM_Forest *forest; 41826d9498aSToby Isaac PetscErrorCode ierr; 419dd8e54a2SToby Isaac 420dd8e54a2SToby Isaac PetscFunctionBegin; 421dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 422ba936b91SToby Isaac forest = (DM_Forest *) dm->data; 42326d9498aSToby Isaac switch (forest->adaptPurpose) { 42426d9498aSToby Isaac case DM_FOREST_KEEP: 425ba936b91SToby Isaac *adapt = forest->adapt; 42626d9498aSToby Isaac break; 42726d9498aSToby Isaac case DM_FOREST_REFINE: 42826d9498aSToby Isaac ierr = DMGetCoarseDM(dm,adapt);CHKERRQ(ierr); 42926d9498aSToby Isaac break; 43026d9498aSToby Isaac case DM_FOREST_COARSEN: 43126d9498aSToby Isaac ierr = DMGetFineDM(dm,adapt);CHKERRQ(ierr); 43226d9498aSToby Isaac break; 43326d9498aSToby Isaac default: 43426d9498aSToby Isaac SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose"); 43526d9498aSToby Isaac } 43626d9498aSToby Isaac PetscFunctionReturn(0); 43726d9498aSToby Isaac } 43826d9498aSToby Isaac 43926d9498aSToby Isaac #undef __FUNCT__ 44026d9498aSToby Isaac #define __FUNCT__ "DMForestSetAdaptivityPurpose" 4419be51f97SToby Isaac /*@ 4429be51f97SToby Isaac DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its 4439be51f97SToby Isaac source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_FOREST_REFINE), coarsening 4449be51f97SToby Isaac (DM_FOREST_COARSEN), or undefined (DM_FOREST_NONE). This only matters for the purposes of reference counting: 4459be51f97SToby Isaac during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse 4469be51f97SToby Isaac relationship (see DMSetFineDM()/DMSetCoarseDM()). If the purpose is not refinement or coarsening, and the user does 4479be51f97SToby Isaac not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can 4489be51f97SToby Isaac cause a memory leak. This method is used by subtypes of DMForest when automatically constructing mesh hierarchies. 4499be51f97SToby Isaac 4509be51f97SToby Isaac Logically collective on dm 4519be51f97SToby Isaac 4529be51f97SToby Isaac Input Parameters: 4539be51f97SToby Isaac + dm - the forest 4549be51f97SToby Isaac - purpose - the adaptivity purpose (DM_FOREST_NONE/DM_FOREST_REFINE/DM_FOREST_COARSEN) 4559be51f97SToby Isaac 4569be51f97SToby Isaac Level: advanced 4579be51f97SToby Isaac 4589be51f97SToby Isaac .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest() 4599be51f97SToby Isaac @*/ 46026d9498aSToby Isaac PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMForestAdaptivityPurpose purpose) 46126d9498aSToby Isaac { 46226d9498aSToby Isaac DM_Forest *forest; 46326d9498aSToby Isaac PetscErrorCode ierr; 46426d9498aSToby Isaac 46526d9498aSToby Isaac PetscFunctionBegin; 46626d9498aSToby Isaac forest = (DM_Forest *) dm->data; 4679be51f97SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup"); 46826d9498aSToby Isaac if (purpose != forest->adaptPurpose) { 46926d9498aSToby Isaac DM adapt; 47026d9498aSToby Isaac 47126d9498aSToby Isaac ierr = DMForestGetAdaptivityForest(dm,&adapt);CHKERRQ(ierr); 47226d9498aSToby Isaac ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr); 47326d9498aSToby Isaac ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr); 47426d9498aSToby Isaac forest->adaptPurpose = purpose; 47526d9498aSToby Isaac ierr = DMForestSetAdaptivityForest(dm,adapt);CHKERRQ(ierr); 47626d9498aSToby Isaac ierr = DMDestroy(&adapt);CHKERRQ(ierr); 47726d9498aSToby Isaac } 478dd8e54a2SToby Isaac PetscFunctionReturn(0); 479dd8e54a2SToby Isaac } 480dd8e54a2SToby Isaac 481dd8e54a2SToby Isaac #undef __FUNCT__ 482dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetAdjacencyDimension" 4839be51f97SToby Isaac /*@ 4849be51f97SToby Isaac DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine 4859be51f97SToby Isaac cell adjacency (for the purposes of partitioning and overlap). 4869be51f97SToby Isaac 4879be51f97SToby Isaac Logically collective on dm 4889be51f97SToby Isaac 4899be51f97SToby Isaac Input Parameters: 4909be51f97SToby Isaac + dm - the forest 4919be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency) 4929be51f97SToby Isaac 4939be51f97SToby Isaac Level: intermediate 4949be51f97SToby Isaac 4959be51f97SToby Isaac .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap() 4969be51f97SToby Isaac @*/ 497dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim) 498dd8e54a2SToby Isaac { 499dd8e54a2SToby Isaac PetscInt dim; 500dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 501dd8e54a2SToby Isaac PetscErrorCode ierr; 502dd8e54a2SToby Isaac 503dd8e54a2SToby Isaac PetscFunctionBegin; 504dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 505ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup"); 506dd8e54a2SToby Isaac if (adjDim < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %d", adjDim); 507dd8e54a2SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 508dd8e54a2SToby Isaac if (adjDim > dim) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %d: %d", dim, adjDim); 509dd8e54a2SToby Isaac forest->adjDim = adjDim; 510dd8e54a2SToby Isaac PetscFunctionReturn(0); 511dd8e54a2SToby Isaac } 512dd8e54a2SToby Isaac 513dd8e54a2SToby Isaac #undef __FUNCT__ 514dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetAdjacencyCodimension" 5159be51f97SToby Isaac /*@ 5169be51f97SToby Isaac DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that, 5179be51f97SToby Isaac e.g., adjacency based on facets can be specified by codimension 1 in all cases) 5189be51f97SToby Isaac 5199be51f97SToby Isaac Logically collective on dm 5209be51f97SToby Isaac 5219be51f97SToby Isaac Input Parameters: 5229be51f97SToby Isaac + dm - the forest 5239be51f97SToby Isaac - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices 5249be51f97SToby Isaac 5259be51f97SToby Isaac Level: intermediate 5269be51f97SToby Isaac 5279be51f97SToby Isaac .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension() 5289be51f97SToby Isaac @*/ 529dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim) 530dd8e54a2SToby Isaac { 531dd8e54a2SToby Isaac PetscInt dim; 532dd8e54a2SToby Isaac PetscErrorCode ierr; 533dd8e54a2SToby Isaac 534dd8e54a2SToby Isaac PetscFunctionBegin; 535dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 536dd8e54a2SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 537dd8e54a2SToby Isaac ierr = DMForestSetAdjacencyDimension(dm,dim-adjCodim);CHKERRQ(ierr); 538dd8e54a2SToby Isaac PetscFunctionReturn(0); 539dd8e54a2SToby Isaac } 540dd8e54a2SToby Isaac 541dd8e54a2SToby Isaac #undef __FUNCT__ 542dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetAdjacencyDimension" 5439be51f97SToby Isaac /*@ 5449be51f97SToby Isaac DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the 5459be51f97SToby Isaac purposes of partitioning and overlap). 5469be51f97SToby Isaac 5479be51f97SToby Isaac Not collective 5489be51f97SToby Isaac 5499be51f97SToby Isaac Input Parameter: 5509be51f97SToby Isaac . dm - the forest 5519be51f97SToby Isaac 5529be51f97SToby Isaac Output Parameter: 5539be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency) 5549be51f97SToby Isaac 5559be51f97SToby Isaac Level: intermediate 5569be51f97SToby Isaac 5579be51f97SToby Isaac .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap() 5589be51f97SToby Isaac @*/ 559dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim) 560dd8e54a2SToby Isaac { 561dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 562dd8e54a2SToby Isaac 563dd8e54a2SToby Isaac PetscFunctionBegin; 564dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 565dd8e54a2SToby Isaac PetscValidIntPointer(adjDim,2); 566dd8e54a2SToby Isaac *adjDim = forest->adjDim; 567dd8e54a2SToby Isaac PetscFunctionReturn(0); 568dd8e54a2SToby Isaac } 569dd8e54a2SToby Isaac 570dd8e54a2SToby Isaac #undef __FUNCT__ 571dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetAdjacencyCodimension" 5729be51f97SToby Isaac /*@ 5739be51f97SToby Isaac DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that, 5749be51f97SToby Isaac e.g., adjacency based on facets can be specified by codimension 1 in all cases) 5759be51f97SToby Isaac 5769be51f97SToby Isaac Not collective 5779be51f97SToby Isaac 5789be51f97SToby Isaac Input Parameter: 5799be51f97SToby Isaac . dm - the forest 5809be51f97SToby Isaac 5819be51f97SToby Isaac Output Parameter: 5829be51f97SToby Isaac . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices 5839be51f97SToby Isaac 5849be51f97SToby Isaac Level: intermediate 5859be51f97SToby Isaac 5869be51f97SToby Isaac .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension() 5879be51f97SToby Isaac @*/ 588dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim) 589dd8e54a2SToby Isaac { 590dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 591dd8e54a2SToby Isaac PetscInt dim; 592dd8e54a2SToby Isaac PetscErrorCode ierr; 593dd8e54a2SToby Isaac 594dd8e54a2SToby Isaac PetscFunctionBegin; 595dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 596dd8e54a2SToby Isaac PetscValidIntPointer(adjCodim,2); 597dd8e54a2SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 598dd8e54a2SToby Isaac *adjCodim = dim - forest->adjDim; 599dd8e54a2SToby Isaac PetscFunctionReturn(0); 600dd8e54a2SToby Isaac } 601dd8e54a2SToby Isaac 602dd8e54a2SToby Isaac #undef __FUNCT__ 603ef51cf95SToby Isaac #define __FUNCT__ "DMForestSetPartitionOverlap" 6049be51f97SToby Isaac /*@ 6059be51f97SToby Isaac DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel 6069be51f97SToby Isaac partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding 6079be51f97SToby Isaac adjacent cells 6089be51f97SToby Isaac 6099be51f97SToby Isaac Logically collective on dm 6109be51f97SToby Isaac 6119be51f97SToby Isaac Input Parameters: 6129be51f97SToby Isaac + dm - the forest 6139be51f97SToby Isaac - overlap - default 0 6149be51f97SToby Isaac 6159be51f97SToby Isaac Level: intermediate 6169be51f97SToby Isaac 6179be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension() 6189be51f97SToby Isaac @*/ 619dd8e54a2SToby Isaac PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap) 620dd8e54a2SToby Isaac { 621dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 622dd8e54a2SToby Isaac 623dd8e54a2SToby Isaac PetscFunctionBegin; 624dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 625ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup"); 626dd8e54a2SToby Isaac if (overlap < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %d", overlap); 627dd8e54a2SToby Isaac forest->overlap = overlap; 628dd8e54a2SToby Isaac PetscFunctionReturn(0); 629dd8e54a2SToby Isaac } 630dd8e54a2SToby Isaac 631dd8e54a2SToby Isaac #undef __FUNCT__ 632dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetPartitionOverlap" 6339be51f97SToby Isaac /*@ 6349be51f97SToby Isaac DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values 6359be51f97SToby Isaac > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells 6369be51f97SToby Isaac 6379be51f97SToby Isaac Not collective 6389be51f97SToby Isaac 6399be51f97SToby Isaac Input Parameter: 6409be51f97SToby Isaac . dm - the forest 6419be51f97SToby Isaac 6429be51f97SToby Isaac Output Parameter: 6439be51f97SToby Isaac . overlap - default 0 6449be51f97SToby Isaac 6459be51f97SToby Isaac Level: intermediate 6469be51f97SToby Isaac 6479be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension() 6489be51f97SToby Isaac @*/ 649dd8e54a2SToby Isaac PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap) 650dd8e54a2SToby Isaac { 651dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 652dd8e54a2SToby Isaac 653dd8e54a2SToby Isaac PetscFunctionBegin; 654dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 655dd8e54a2SToby Isaac PetscValidIntPointer(overlap,2); 656dd8e54a2SToby Isaac *overlap = forest->overlap; 657dd8e54a2SToby Isaac PetscFunctionReturn(0); 658dd8e54a2SToby Isaac } 659dd8e54a2SToby Isaac 660dd8e54a2SToby Isaac #undef __FUNCT__ 661dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetMinimumRefinement" 6629be51f97SToby Isaac /*@ 6639be51f97SToby Isaac DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base 6649be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. If the forest is being created by coarsening a previous forest 6659be51f97SToby Isaac (see DMForestGetAdaptivityForest()) this limits the amount of coarsening. 6669be51f97SToby Isaac 6679be51f97SToby Isaac Logically collective on dm 6689be51f97SToby Isaac 6699be51f97SToby Isaac Input Parameters: 6709be51f97SToby Isaac + dm - the forest 6719be51f97SToby Isaac - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 6729be51f97SToby Isaac 6739be51f97SToby Isaac Level: intermediate 6749be51f97SToby Isaac 6759be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 6769be51f97SToby Isaac @*/ 677dd8e54a2SToby Isaac PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement) 678dd8e54a2SToby Isaac { 679dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 680dd8e54a2SToby Isaac 681dd8e54a2SToby Isaac PetscFunctionBegin; 682dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 683ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup"); 684dd8e54a2SToby Isaac forest->minRefinement = minRefinement; 685dd8e54a2SToby Isaac PetscFunctionReturn(0); 686dd8e54a2SToby Isaac } 687dd8e54a2SToby Isaac 688dd8e54a2SToby Isaac #undef __FUNCT__ 689dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetMinimumRefinement" 6909be51f97SToby Isaac /*@ 6919be51f97SToby Isaac DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see 6929be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. If the forest is being created by coarsening a previous forest (see 6939be51f97SToby Isaac DMForestGetAdaptivityForest()), this limits the amount of coarsening. 6949be51f97SToby Isaac 6959be51f97SToby Isaac Not collective 6969be51f97SToby Isaac 6979be51f97SToby Isaac Input Parameter: 6989be51f97SToby Isaac . dm - the forest 6999be51f97SToby Isaac 7009be51f97SToby Isaac Output Parameter: 7019be51f97SToby Isaac . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7029be51f97SToby Isaac 7039be51f97SToby Isaac Level: intermediate 7049be51f97SToby Isaac 7059be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestGetMaximumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 7069be51f97SToby Isaac @*/ 707dd8e54a2SToby Isaac PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement) 708dd8e54a2SToby Isaac { 709dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 710dd8e54a2SToby Isaac 711dd8e54a2SToby Isaac PetscFunctionBegin; 712dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 713dd8e54a2SToby Isaac PetscValidIntPointer(minRefinement,2); 714dd8e54a2SToby Isaac *minRefinement = forest->minRefinement; 715dd8e54a2SToby Isaac PetscFunctionReturn(0); 716dd8e54a2SToby Isaac } 717dd8e54a2SToby Isaac 718dd8e54a2SToby Isaac #undef __FUNCT__ 71956ba9f64SToby Isaac #define __FUNCT__ "DMForestSetInitialRefinement" 7209be51f97SToby Isaac /*@ 7219be51f97SToby Isaac DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base 7229be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. 7239be51f97SToby Isaac 7249be51f97SToby Isaac Logically collective on dm 7259be51f97SToby Isaac 7269be51f97SToby Isaac Input Parameters: 7279be51f97SToby Isaac + dm - the forest 7289be51f97SToby Isaac - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7299be51f97SToby Isaac 7309be51f97SToby Isaac Level: intermediate 7319be51f97SToby Isaac 7329be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM() 7339be51f97SToby Isaac @*/ 73456ba9f64SToby Isaac PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement) 73556ba9f64SToby Isaac { 73656ba9f64SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 73756ba9f64SToby Isaac 73856ba9f64SToby Isaac PetscFunctionBegin; 73956ba9f64SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 74056ba9f64SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup"); 74156ba9f64SToby Isaac forest->initRefinement = initRefinement; 74256ba9f64SToby Isaac PetscFunctionReturn(0); 74356ba9f64SToby Isaac } 74456ba9f64SToby Isaac 74556ba9f64SToby Isaac #undef __FUNCT__ 74656ba9f64SToby Isaac #define __FUNCT__ "DMForestGetInitialRefinement" 7479be51f97SToby Isaac /*@ 7489be51f97SToby Isaac DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see 7499be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. 7509be51f97SToby Isaac 7519be51f97SToby Isaac Not collective 7529be51f97SToby Isaac 7539be51f97SToby Isaac Input Parameter: 7549be51f97SToby Isaac . dm - the forest 7559be51f97SToby Isaac 7569be51f97SToby Isaac Output Paramater: 7579be51f97SToby Isaac . initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7589be51f97SToby Isaac 7599be51f97SToby Isaac Level: intermediate 7609be51f97SToby Isaac 7619be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM() 7629be51f97SToby Isaac @*/ 76356ba9f64SToby Isaac PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement) 76456ba9f64SToby Isaac { 76556ba9f64SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 76656ba9f64SToby Isaac 76756ba9f64SToby Isaac PetscFunctionBegin; 76856ba9f64SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 76956ba9f64SToby Isaac PetscValidIntPointer(initRefinement,2); 77056ba9f64SToby Isaac *initRefinement = forest->initRefinement; 77156ba9f64SToby Isaac PetscFunctionReturn(0); 77256ba9f64SToby Isaac } 77356ba9f64SToby Isaac 77456ba9f64SToby Isaac #undef __FUNCT__ 775c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetMaximumRefinement" 7769be51f97SToby Isaac /*@ 7779be51f97SToby Isaac DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base 7789be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. If the forest is being created by refining a previous forest 7799be51f97SToby Isaac (see DMForestGetAdaptivityForest()), this limits the amount of refinement. 7809be51f97SToby Isaac 7819be51f97SToby Isaac Logically collective on dm 7829be51f97SToby Isaac 7839be51f97SToby Isaac Input Parameters: 7849be51f97SToby Isaac + dm - the forest 7859be51f97SToby Isaac - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7869be51f97SToby Isaac 7879be51f97SToby Isaac Level: intermediate 7889be51f97SToby Isaac 7899be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityDM() 7909be51f97SToby Isaac @*/ 791c7eeac06SToby Isaac PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement) 792dd8e54a2SToby Isaac { 793dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 794dd8e54a2SToby Isaac 795dd8e54a2SToby Isaac PetscFunctionBegin; 796dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 797ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup"); 798c7eeac06SToby Isaac forest->maxRefinement = maxRefinement; 799dd8e54a2SToby Isaac PetscFunctionReturn(0); 800dd8e54a2SToby Isaac } 801dd8e54a2SToby Isaac 802dd8e54a2SToby Isaac #undef __FUNCT__ 803c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetMaximumRefinement" 8049be51f97SToby Isaac /*@ 8059be51f97SToby Isaac DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see 8069be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. If the forest is being created by refining a previous forest (see 8079be51f97SToby Isaac DMForestGetAdaptivityForest()), this limits the amount of refinement. 8089be51f97SToby Isaac 8099be51f97SToby Isaac Not collective 8109be51f97SToby Isaac 8119be51f97SToby Isaac Input Parameter: 8129be51f97SToby Isaac . dm - the forest 8139be51f97SToby Isaac 8149be51f97SToby Isaac Output Parameter: 8159be51f97SToby Isaac . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 8169be51f97SToby Isaac 8179be51f97SToby Isaac Level: intermediate 8189be51f97SToby Isaac 8199be51f97SToby Isaac .seealso: DMForestSetMaximumRefinement(), DMForestGetMinimumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 8209be51f97SToby Isaac @*/ 821c7eeac06SToby Isaac PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement) 822dd8e54a2SToby Isaac { 823dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 824dd8e54a2SToby Isaac 825dd8e54a2SToby Isaac PetscFunctionBegin; 826dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 827c7eeac06SToby Isaac PetscValidIntPointer(maxRefinement,2); 828c7eeac06SToby Isaac *maxRefinement = forest->maxRefinement; 829dd8e54a2SToby Isaac PetscFunctionReturn(0); 830dd8e54a2SToby Isaac } 831c7eeac06SToby Isaac 832c7eeac06SToby Isaac #undef __FUNCT__ 833c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityStrategy" 8349be51f97SToby Isaac /*@C 8359be51f97SToby Isaac DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes. 8369be51f97SToby Isaac Subtypes of DMForest may define their own strategies. Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree 8379be51f97SToby Isaac for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to 8389be51f97SToby Isaac specify refinement/coarsening. 8399be51f97SToby Isaac 8409be51f97SToby Isaac Logically collective on dm 8419be51f97SToby Isaac 8429be51f97SToby Isaac Input Parameters: 8439be51f97SToby Isaac + dm - the forest 8449be51f97SToby Isaac - adaptStrategy - default DMFORESTADAPTALL 8459be51f97SToby Isaac 8469be51f97SToby Isaac Level: advanced 8479be51f97SToby Isaac 8489be51f97SToby Isaac .seealso: DMForestGetAdaptivityStrategy() 8499be51f97SToby Isaac @*/ 850c7eeac06SToby Isaac PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy) 851c7eeac06SToby Isaac { 852c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 853c7eeac06SToby Isaac PetscErrorCode ierr; 854c7eeac06SToby Isaac 855c7eeac06SToby Isaac PetscFunctionBegin; 856c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 857c7eeac06SToby Isaac ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr); 858a73e2921SToby Isaac ierr = PetscStrallocpy((const char *) adaptStrategy,(char **)&forest->adaptStrategy);CHKERRQ(ierr); 859c7eeac06SToby Isaac PetscFunctionReturn(0); 860c7eeac06SToby Isaac } 861c7eeac06SToby Isaac 862c7eeac06SToby Isaac #undef __FUNCT__ 863c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityStrategy" 8649be51f97SToby Isaac /*@C 8659be51f97SToby Isaac DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes. Subtypes 8669be51f97SToby Isaac of DMForest may define their own strategies. Two default strategies are DMFORESTADAPTALL, which indicates that all 8679be51f97SToby Isaac processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only 8689be51f97SToby Isaac one process needs to specify refinement/coarsening. 8699be51f97SToby Isaac 8709be51f97SToby Isaac Not collective 8719be51f97SToby Isaac 8729be51f97SToby Isaac Input Parameter: 8739be51f97SToby Isaac . dm - the forest 8749be51f97SToby Isaac 8759be51f97SToby Isaac Output Parameter: 8769be51f97SToby Isaac . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL) 8779be51f97SToby Isaac 8789be51f97SToby Isaac Level: advanced 8799be51f97SToby Isaac 8809be51f97SToby Isaac .seealso: DMForestSetAdaptivityStrategy() 8819be51f97SToby Isaac @*/ 882c7eeac06SToby Isaac PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy) 883c7eeac06SToby Isaac { 884c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 885c7eeac06SToby Isaac 886c7eeac06SToby Isaac PetscFunctionBegin; 887c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 888c7eeac06SToby Isaac PetscValidPointer(adaptStrategy,2); 889c7eeac06SToby Isaac *adaptStrategy = forest->adaptStrategy; 890c7eeac06SToby Isaac PetscFunctionReturn(0); 891c7eeac06SToby Isaac } 892c7eeac06SToby Isaac 893c7eeac06SToby Isaac #undef __FUNCT__ 894*bf9b5d84SToby Isaac #define __FUNCT__ "DMForestSetComputeAdaptivitySF" 895*bf9b5d84SToby Isaac /*@ 896*bf9b5d84SToby Isaac DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed 897*bf9b5d84SToby 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(). 898*bf9b5d84SToby Isaac 899*bf9b5d84SToby Isaac Logically collective on dm 900*bf9b5d84SToby Isaac 901*bf9b5d84SToby Isaac Input Parameters: 902*bf9b5d84SToby Isaac + dm - the post-adaptation forest 903*bf9b5d84SToby Isaac - computeSF - default PETSC_TRUE 904*bf9b5d84SToby Isaac 905*bf9b5d84SToby Isaac Level: advanced 906*bf9b5d84SToby Isaac 907*bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF() 908*bf9b5d84SToby Isaac @*/ 909*bf9b5d84SToby Isaac PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF) 910*bf9b5d84SToby Isaac { 911*bf9b5d84SToby Isaac DM_Forest *forest; 912*bf9b5d84SToby Isaac 913*bf9b5d84SToby Isaac PetscFunctionBegin; 914*bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 915*bf9b5d84SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called"); 916*bf9b5d84SToby Isaac forest = (DM_Forest *) dm->data; 917*bf9b5d84SToby Isaac forest->computeAdaptSF = computeSF; 918*bf9b5d84SToby Isaac PetscFunctionReturn(0); 919*bf9b5d84SToby Isaac } 920*bf9b5d84SToby Isaac 921*bf9b5d84SToby Isaac #undef __FUNCT__ 922*bf9b5d84SToby Isaac #define __FUNCT__ "DMForestGetComputeAdaptivitySF" 923*bf9b5d84SToby Isaac /*@ 924*bf9b5d84SToby Isaac DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the 925*bf9b5d84SToby Isaac pre-adaptation forest to the post-adaptiation forest. After DMSetUp() is called, these transfer PetscSFs can be 926*bf9b5d84SToby Isaac accessed with DMForestGetAdaptivitySF(). 927*bf9b5d84SToby Isaac 928*bf9b5d84SToby Isaac Not collective 929*bf9b5d84SToby Isaac 930*bf9b5d84SToby Isaac Input Parameter: 931*bf9b5d84SToby Isaac . dm - the post-adaptation forest 932*bf9b5d84SToby Isaac 933*bf9b5d84SToby Isaac Output Parameter: 934*bf9b5d84SToby Isaac . computeSF - default PETSC_TRUE 935*bf9b5d84SToby Isaac 936*bf9b5d84SToby Isaac Level: advanced 937*bf9b5d84SToby Isaac 938*bf9b5d84SToby Isaac .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF() 939*bf9b5d84SToby Isaac @*/ 940*bf9b5d84SToby Isaac PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF) 941*bf9b5d84SToby Isaac { 942*bf9b5d84SToby Isaac DM_Forest *forest; 943*bf9b5d84SToby Isaac 944*bf9b5d84SToby Isaac PetscFunctionBegin; 945*bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 946*bf9b5d84SToby Isaac forest = (DM_Forest *) dm->data; 947*bf9b5d84SToby Isaac *computeSF = forest->computeAdaptSF; 948*bf9b5d84SToby Isaac PetscFunctionReturn(0); 949*bf9b5d84SToby Isaac } 950*bf9b5d84SToby Isaac 951*bf9b5d84SToby Isaac #undef __FUNCT__ 952*bf9b5d84SToby Isaac #define __FUNCT__ "DMForestGetAdaptivitySF" 953*bf9b5d84SToby Isaac /*@ 954*bf9b5d84SToby Isaac DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest. 955*bf9b5d84SToby Isaac Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be 956*bf9b5d84SToby Isaac some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa. Therefore there are two 957*bf9b5d84SToby Isaac PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates 958*bf9b5d84SToby Isaac pre-adaptation fine cells to post-adaptation coarse cells. 959*bf9b5d84SToby Isaac 960*bf9b5d84SToby Isaac Not collective 961*bf9b5d84SToby Isaac 962*bf9b5d84SToby Isaac Input Parameter: 963*bf9b5d84SToby Isaac dm - the post-adaptation forest 964*bf9b5d84SToby Isaac 965*bf9b5d84SToby Isaac Output Parameter: 966*bf9b5d84SToby Isaac coarseToFineSF - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post- 967*bf9b5d84SToby Isaac fineToCoarseSF - pre-adaptation fine cells to post-adaptation coarse cells: BCase goes from post- to pre- 968*bf9b5d84SToby Isaac 969*bf9b5d84SToby Isaac Level: advanced 970*bf9b5d84SToby Isaac 971*bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF() 972*bf9b5d84SToby Isaac @*/ 973*bf9b5d84SToby Isaac PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *coarseToFineSF, PetscSF *fineToCoarseSF) 974*bf9b5d84SToby Isaac { 975*bf9b5d84SToby Isaac DM_Forest *forest; 976*bf9b5d84SToby Isaac PetscErrorCode ierr; 977*bf9b5d84SToby Isaac 978*bf9b5d84SToby Isaac PetscFunctionBegin; 979*bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 980*bf9b5d84SToby Isaac ierr = DMSetUp(dm);CHKERRQ(ierr); 981*bf9b5d84SToby Isaac forest = (DM_Forest *) dm->data; 982*bf9b5d84SToby Isaac if (coarseToFineSF) { 983*bf9b5d84SToby Isaac *coarseToFineSF = forest->adaptSFCoarseToFine; 984*bf9b5d84SToby Isaac } 985*bf9b5d84SToby Isaac if (fineToCoarseSF) { 986*bf9b5d84SToby Isaac *fineToCoarseSF = forest->adaptSFFineToCoarse; 987*bf9b5d84SToby Isaac } 988*bf9b5d84SToby Isaac PetscFunctionReturn(0); 989*bf9b5d84SToby Isaac } 990*bf9b5d84SToby Isaac 991*bf9b5d84SToby Isaac #undef __FUNCT__ 992c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetGradeFactor" 9939be51f97SToby Isaac /*@ 9949be51f97SToby Isaac DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to 9959be51f97SToby Isaac indicate that the diameter of neighboring cells should differ by at most a factor of 2. Subtypes of DMForest may 9969be51f97SToby Isaac only support one particular choice of grading factor. 9979be51f97SToby Isaac 9989be51f97SToby Isaac Logically collective on dm 9999be51f97SToby Isaac 10009be51f97SToby Isaac Input Parameters: 10019be51f97SToby Isaac + dm - the forest 10029be51f97SToby Isaac - grade - the grading factor 10039be51f97SToby Isaac 10049be51f97SToby Isaac Level: advanced 10059be51f97SToby Isaac 10069be51f97SToby Isaac .seealso: DMForestGetGradeFactor() 10079be51f97SToby Isaac @*/ 1008c7eeac06SToby Isaac PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade) 1009c7eeac06SToby Isaac { 1010c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1011c7eeac06SToby Isaac 1012c7eeac06SToby Isaac PetscFunctionBegin; 1013c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1014ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup"); 1015c7eeac06SToby Isaac forest->gradeFactor = grade; 1016c7eeac06SToby Isaac PetscFunctionReturn(0); 1017c7eeac06SToby Isaac } 1018c7eeac06SToby Isaac 1019c7eeac06SToby Isaac #undef __FUNCT__ 1020c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetGradeFactor" 10219be51f97SToby Isaac /*@ 10229be51f97SToby Isaac DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of 10239be51f97SToby Isaac neighboring cells should differ by at most a factor of 2. Subtypes of DMForest may only support one particular 10249be51f97SToby Isaac choice of grading factor. 10259be51f97SToby Isaac 10269be51f97SToby Isaac Not collective 10279be51f97SToby Isaac 10289be51f97SToby Isaac Input Parameter: 10299be51f97SToby Isaac . dm - the forest 10309be51f97SToby Isaac 10319be51f97SToby Isaac Output Parameter: 10329be51f97SToby Isaac . grade - the grading factor 10339be51f97SToby Isaac 10349be51f97SToby Isaac Level: advanced 10359be51f97SToby Isaac 10369be51f97SToby Isaac .seealso: DMForestSetGradeFactor() 10379be51f97SToby Isaac @*/ 1038c7eeac06SToby Isaac PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade) 1039c7eeac06SToby Isaac { 1040c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1041c7eeac06SToby Isaac 1042c7eeac06SToby Isaac PetscFunctionBegin; 1043c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1044c7eeac06SToby Isaac PetscValidIntPointer(grade,2); 1045c7eeac06SToby Isaac *grade = forest->gradeFactor; 1046c7eeac06SToby Isaac PetscFunctionReturn(0); 1047c7eeac06SToby Isaac } 1048c7eeac06SToby Isaac 1049c7eeac06SToby Isaac #undef __FUNCT__ 1050ef51cf95SToby Isaac #define __FUNCT__ "DMForestSetCellWeightFactor" 10519be51f97SToby Isaac /*@ 10529be51f97SToby Isaac DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes 10539be51f97SToby Isaac the cell weight (see DMForestSetCellWeights()) when calculating partitions. The final weight of a cell will be 10549be51f97SToby Isaac (cellWeight) * (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on 10559be51f97SToby Isaac its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the 10569be51f97SToby Isaac computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement. 10579be51f97SToby Isaac 10589be51f97SToby Isaac Logically collective on dm 10599be51f97SToby Isaac 10609be51f97SToby Isaac Input Parameters: 10619be51f97SToby Isaac + dm - the forest 10629be51f97SToby Isaac - weightsFactors - default 1. 10639be51f97SToby Isaac 10649be51f97SToby Isaac Level: advanced 10659be51f97SToby Isaac 10669be51f97SToby Isaac .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights() 10679be51f97SToby Isaac @*/ 1068ef51cf95SToby Isaac PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor) 1069c7eeac06SToby Isaac { 1070c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1071c7eeac06SToby Isaac 1072c7eeac06SToby Isaac PetscFunctionBegin; 1073c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1074ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup"); 1075c7eeac06SToby Isaac forest->weightsFactor = weightsFactor; 1076c7eeac06SToby Isaac PetscFunctionReturn(0); 1077c7eeac06SToby Isaac } 1078c7eeac06SToby Isaac 1079c7eeac06SToby Isaac #undef __FUNCT__ 1080ef51cf95SToby Isaac #define __FUNCT__ "DMForestGetCellWeightFactor" 10819be51f97SToby Isaac /*@ 10829be51f97SToby Isaac DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see 10839be51f97SToby Isaac DMForestSetCellWeights()) when calculating partitions. The final weight of a cell will be (cellWeight) * 10849be51f97SToby Isaac (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on its level; a 10859be51f97SToby Isaac factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation 10869be51f97SToby Isaac associated with a cell is multiplied by a factor of 2 for each additional level of refinement. 10879be51f97SToby Isaac 10889be51f97SToby Isaac Not collective 10899be51f97SToby Isaac 10909be51f97SToby Isaac Input Parameter: 10919be51f97SToby Isaac . dm - the forest 10929be51f97SToby Isaac 10939be51f97SToby Isaac Output Parameter: 10949be51f97SToby Isaac . weightsFactors - default 1. 10959be51f97SToby Isaac 10969be51f97SToby Isaac Level: advanced 10979be51f97SToby Isaac 10989be51f97SToby Isaac .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights() 10999be51f97SToby Isaac @*/ 1100ef51cf95SToby Isaac PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor) 1101c7eeac06SToby Isaac { 1102c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1103c7eeac06SToby Isaac 1104c7eeac06SToby Isaac PetscFunctionBegin; 1105c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1106c7eeac06SToby Isaac PetscValidRealPointer(weightsFactor,2); 1107c7eeac06SToby Isaac *weightsFactor = forest->weightsFactor; 1108c7eeac06SToby Isaac PetscFunctionReturn(0); 1109c7eeac06SToby Isaac } 1110c7eeac06SToby Isaac 1111c7eeac06SToby Isaac #undef __FUNCT__ 1112c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellChart" 11139be51f97SToby Isaac /*@ 11149be51f97SToby Isaac DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process 11159be51f97SToby Isaac 11169be51f97SToby Isaac Not collective 11179be51f97SToby Isaac 11189be51f97SToby Isaac Input Parameter: 11199be51f97SToby Isaac . dm - the forest 11209be51f97SToby Isaac 11219be51f97SToby Isaac Output Parameters: 11229be51f97SToby Isaac + cStart - the first cell on this process 11239be51f97SToby Isaac - cEnd - one after the final cell on this process 11249be51f97SToby Isaac 11259be51f97SToby Isaac level: intermediate 11269be51f97SToby Isaac 11279be51f97SToby Isaac .seealso: DMForestGetCellSF() 11289be51f97SToby Isaac @*/ 1129c7eeac06SToby Isaac PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd) 1130c7eeac06SToby Isaac { 1131c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1132c7eeac06SToby Isaac PetscErrorCode ierr; 1133c7eeac06SToby Isaac 1134c7eeac06SToby Isaac PetscFunctionBegin; 1135c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1136c7eeac06SToby Isaac PetscValidIntPointer(cStart,2); 1137c7eeac06SToby Isaac PetscValidIntPointer(cEnd,2); 1138c7eeac06SToby Isaac if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) { 1139c7eeac06SToby Isaac ierr = forest->createcellchart(dm,&forest->cStart,&forest->cEnd);CHKERRQ(ierr); 1140c7eeac06SToby Isaac } 1141c7eeac06SToby Isaac *cStart = forest->cStart; 1142c7eeac06SToby Isaac *cEnd = forest->cEnd; 1143c7eeac06SToby Isaac PetscFunctionReturn(0); 1144c7eeac06SToby Isaac } 1145c7eeac06SToby Isaac 1146c7eeac06SToby Isaac #undef __FUNCT__ 1147c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellSF" 11489be51f97SToby Isaac /*@ 11499be51f97SToby Isaac DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes 11509be51f97SToby Isaac 11519be51f97SToby Isaac Not collective 11529be51f97SToby Isaac 11539be51f97SToby Isaac Input Parameter: 11549be51f97SToby Isaac . dm - the forest 11559be51f97SToby Isaac 11569be51f97SToby Isaac Output Parameter: 11579be51f97SToby Isaac . cellSF - the PetscSF 11589be51f97SToby Isaac 11599be51f97SToby Isaac level: intermediate 11609be51f97SToby Isaac 11619be51f97SToby Isaac .seealso: DMForestGetCellChart() 11629be51f97SToby Isaac @*/ 1163c7eeac06SToby Isaac PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF) 1164c7eeac06SToby Isaac { 1165c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1166c7eeac06SToby Isaac PetscErrorCode ierr; 1167c7eeac06SToby Isaac 1168c7eeac06SToby Isaac PetscFunctionBegin; 1169c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1170c7eeac06SToby Isaac PetscValidPointer(cellSF,2); 1171c7eeac06SToby Isaac if ((!forest->cellSF) && forest->createcellsf) { 1172c7eeac06SToby Isaac ierr = forest->createcellsf(dm,&forest->cellSF);CHKERRQ(ierr); 1173c7eeac06SToby Isaac } 1174c7eeac06SToby Isaac *cellSF = forest->cellSF; 1175c7eeac06SToby Isaac PetscFunctionReturn(0); 1176c7eeac06SToby Isaac } 1177c7eeac06SToby Isaac 1178c7eeac06SToby Isaac #undef __FUNCT__ 1179ebdf65a2SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityLabel" 11809be51f97SToby Isaac /*@C 11819be51f97SToby Isaac DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see 11829be51f97SToby Isaac DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination). The 11839be51f97SToby Isaac interpretation of the label values is up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and 11849be51f97SToby Isaac DM_FOREST_COARSEN have been reserved as choices that should be accepted by all subtypes. 11859be51f97SToby Isaac 11869be51f97SToby Isaac Logically collective on dm 11879be51f97SToby Isaac 11889be51f97SToby Isaac Input Parameters: 11899be51f97SToby Isaac - dm - the forest 11909be51f97SToby Isaac + adaptLabel - the name of the label in the pre-adaptation forest 11919be51f97SToby Isaac 11929be51f97SToby Isaac Level: intermediate 11939be51f97SToby Isaac 11949be51f97SToby Isaac .seealso DMForestGetAdaptivityLabel() 11959be51f97SToby Isaac @*/ 1196ebdf65a2SToby Isaac PetscErrorCode DMForestSetAdaptivityLabel(DM dm, const char * adaptLabel) 1197c7eeac06SToby Isaac { 1198c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1199c7eeac06SToby Isaac PetscErrorCode ierr; 1200c7eeac06SToby Isaac 1201c7eeac06SToby Isaac PetscFunctionBegin; 1202c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1203ebdf65a2SToby Isaac ierr = PetscFree(forest->adaptLabel);CHKERRQ(ierr); 1204ebdf65a2SToby Isaac ierr = PetscStrallocpy(adaptLabel,&forest->adaptLabel);CHKERRQ(ierr); 1205c7eeac06SToby Isaac PetscFunctionReturn(0); 1206c7eeac06SToby Isaac } 1207c7eeac06SToby Isaac 1208c7eeac06SToby Isaac #undef __FUNCT__ 1209ebdf65a2SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityLabel" 12109be51f97SToby Isaac /*@C 12119be51f97SToby Isaac DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that 12129be51f97SToby Isaac holds the adaptation flags (refinement, coarsening, or some combination). The interpretation of the label values is 12139be51f97SToby Isaac up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and DM_FOREST_COARSEN have been reserved as 12149be51f97SToby Isaac choices that should be accepted by all subtypes. 12159be51f97SToby Isaac 12169be51f97SToby Isaac Not collective 12179be51f97SToby Isaac 12189be51f97SToby Isaac Input Parameter: 12199be51f97SToby Isaac . dm - the forest 12209be51f97SToby Isaac 12219be51f97SToby Isaac Output Parameter: 12229be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest 12239be51f97SToby Isaac 12249be51f97SToby Isaac Level: intermediate 12259be51f97SToby Isaac 12269be51f97SToby Isaac .seealso DMForestSetAdaptivityLabel() 12279be51f97SToby Isaac @*/ 1228ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityLabel(DM dm, const char ** adaptLabel) 1229c7eeac06SToby Isaac { 1230c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1231c7eeac06SToby Isaac 1232c7eeac06SToby Isaac PetscFunctionBegin; 1233c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1234ba936b91SToby Isaac *adaptLabel = forest->adaptLabel; 1235c7eeac06SToby Isaac PetscFunctionReturn(0); 1236c7eeac06SToby Isaac } 1237c7eeac06SToby Isaac 1238c7eeac06SToby Isaac #undef __FUNCT__ 1239c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetCellWeights" 12409be51f97SToby Isaac /*@ 12419be51f97SToby Isaac DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current 12429be51f97SToby Isaac process: weights are used to determine parallel partitioning. Partitions will be created so that each process's 12439be51f97SToby Isaac ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight 12449be51f97SToby Isaac of 1. 12459be51f97SToby Isaac 12469be51f97SToby Isaac Logically collective on dm 12479be51f97SToby Isaac 12489be51f97SToby Isaac Input Parameters: 12499be51f97SToby Isaac + dm - the forest 12509be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1. 12519be51f97SToby Isaac - copyMode - how weights should reference weights 12529be51f97SToby Isaac 12539be51f97SToby Isaac Level: advanced 12549be51f97SToby Isaac 12559be51f97SToby Isaac .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity() 12569be51f97SToby Isaac @*/ 1257c7eeac06SToby Isaac PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode) 1258c7eeac06SToby Isaac { 1259c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1260c7eeac06SToby Isaac PetscInt cStart, cEnd; 1261c7eeac06SToby Isaac PetscErrorCode ierr; 1262c7eeac06SToby Isaac 1263c7eeac06SToby Isaac PetscFunctionBegin; 1264c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1265c7eeac06SToby Isaac ierr = DMForestGetCellChart(dm,&cStart,&cEnd);CHKERRQ(ierr); 1266c7eeac06SToby Isaac if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd); 1267c7eeac06SToby Isaac if (copyMode == PETSC_COPY_VALUES) { 1268c7eeac06SToby Isaac if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) { 1269c7eeac06SToby Isaac ierr = PetscMalloc1(cEnd-cStart,&forest->cellWeights);CHKERRQ(ierr); 1270c7eeac06SToby Isaac } 1271c7eeac06SToby Isaac ierr = PetscMemcpy(forest->cellWeights,weights,(cEnd-cStart)*sizeof(*weights));CHKERRQ(ierr); 1272c7eeac06SToby Isaac forest->cellWeightsCopyMode = PETSC_OWN_POINTER; 1273c7eeac06SToby Isaac PetscFunctionReturn(0); 1274c7eeac06SToby Isaac } 1275c7eeac06SToby Isaac if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) { 1276c7eeac06SToby Isaac ierr = PetscFree(forest->cellWeights);CHKERRQ(ierr); 1277c7eeac06SToby Isaac } 1278c7eeac06SToby Isaac forest->cellWeights = weights; 1279c7eeac06SToby Isaac forest->cellWeightsCopyMode = copyMode; 1280c7eeac06SToby Isaac PetscFunctionReturn(0); 1281c7eeac06SToby Isaac } 1282c7eeac06SToby Isaac 1283c7eeac06SToby Isaac #undef __FUNCT__ 1284c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellWeights" 12859be51f97SToby Isaac /*@ 12869be51f97SToby Isaac DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current 12879be51f97SToby Isaac process: weights are used to determine parallel partitioning. Partitions will be created so that each process's 12889be51f97SToby Isaac ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight 12899be51f97SToby Isaac of 1. 12909be51f97SToby Isaac 12919be51f97SToby Isaac Not collective 12929be51f97SToby Isaac 12939be51f97SToby Isaac Input Parameter: 12949be51f97SToby Isaac . dm - the forest 12959be51f97SToby Isaac 12969be51f97SToby Isaac Output Parameter: 12979be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1. 12989be51f97SToby Isaac 12999be51f97SToby Isaac Level: advanced 13009be51f97SToby Isaac 13019be51f97SToby Isaac .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity() 13029be51f97SToby Isaac @*/ 1303c7eeac06SToby Isaac PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights) 1304c7eeac06SToby Isaac { 1305c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1306c7eeac06SToby Isaac 1307c7eeac06SToby Isaac PetscFunctionBegin; 1308c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1309c7eeac06SToby Isaac PetscValidPointer(weights,2); 1310c7eeac06SToby Isaac *weights = forest->cellWeights; 1311c7eeac06SToby Isaac PetscFunctionReturn(0); 1312c7eeac06SToby Isaac } 1313c7eeac06SToby Isaac 1314c7eeac06SToby Isaac #undef __FUNCT__ 1315c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetWeightCapacity" 13169be51f97SToby Isaac /*@ 13179be51f97SToby Isaac DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning 13189be51f97SToby Isaac a pre-adaptation forest (see DMForestGetAdaptivityForest()). After partitioning, the ratio of the weight of each 13199be51f97SToby Isaac process's cells to the process's capacity will be roughly equal for all processes. A capacity of 0 indicates that 13209be51f97SToby Isaac the current process should not have any cells after repartitioning. 13219be51f97SToby Isaac 13229be51f97SToby Isaac Logically Collective on dm 13239be51f97SToby Isaac 13249be51f97SToby Isaac Input parameters: 13259be51f97SToby Isaac + dm - the forest 13269be51f97SToby Isaac - capacity - this process's capacity 13279be51f97SToby Isaac 13289be51f97SToby Isaac Level: advanced 13299be51f97SToby Isaac 13309be51f97SToby Isaac .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor() 13319be51f97SToby Isaac @*/ 1332c7eeac06SToby Isaac PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity) 1333c7eeac06SToby Isaac { 1334c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1335c7eeac06SToby Isaac 1336c7eeac06SToby Isaac PetscFunctionBegin; 1337c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1338ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup"); 1339c7eeac06SToby Isaac if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity); 1340c7eeac06SToby Isaac forest->weightCapacity = capacity; 1341c7eeac06SToby Isaac PetscFunctionReturn(0); 1342c7eeac06SToby Isaac } 1343c7eeac06SToby Isaac 1344c7eeac06SToby Isaac #undef __FUNCT__ 1345c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetWeightCapacity" 13469be51f97SToby Isaac /*@ 13479be51f97SToby Isaac DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see 13489be51f97SToby Isaac DMForestGetAdaptivityForest()). After partitioning, the ratio of the weight of each process's cells to the 13499be51f97SToby Isaac process's capacity will be roughly equal for all processes. A capacity of 0 indicates that the current process 13509be51f97SToby Isaac should not have any cells after repartitioning. 13519be51f97SToby Isaac 13529be51f97SToby Isaac Not collective 13539be51f97SToby Isaac 13549be51f97SToby Isaac Input parameter: 13559be51f97SToby Isaac . dm - the forest 13569be51f97SToby Isaac 13579be51f97SToby Isaac Output parameter: 13589be51f97SToby Isaac . capacity - this process's capacity 13599be51f97SToby Isaac 13609be51f97SToby Isaac Level: advanced 13619be51f97SToby Isaac 13629be51f97SToby Isaac .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor() 13639be51f97SToby Isaac @*/ 1364c7eeac06SToby Isaac PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity) 1365c7eeac06SToby Isaac { 1366c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1367c7eeac06SToby Isaac 1368c7eeac06SToby Isaac PetscFunctionBegin; 1369c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1370c7eeac06SToby Isaac PetscValidRealPointer(capacity,2); 1371c7eeac06SToby Isaac *capacity = forest->weightCapacity; 1372c7eeac06SToby Isaac PetscFunctionReturn(0); 1373c7eeac06SToby Isaac } 1374c7eeac06SToby Isaac 1375dd8e54a2SToby Isaac #undef __FUNCT__ 1376db4d5e8cSToby Isaac #define __FUNCT__ "DMSetFromOptions_Forest" 13770709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm) 1378db4d5e8cSToby Isaac { 1379db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 138056ba9f64SToby Isaac PetscBool flg, flg1, flg2, flg3, flg4; 1381dd8e54a2SToby Isaac DMForestTopology oldTopo; 1382c7eeac06SToby Isaac char stringBuffer[256]; 1383dd8e54a2SToby Isaac PetscViewer viewer; 1384dd8e54a2SToby Isaac PetscViewerFormat format; 138556ba9f64SToby Isaac PetscInt adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade; 1386c7eeac06SToby Isaac PetscReal weightsFactor; 1387c7eeac06SToby Isaac DMForestAdaptivityStrategy adaptStrategy; 1388db4d5e8cSToby Isaac PetscErrorCode ierr; 1389db4d5e8cSToby Isaac 1390db4d5e8cSToby Isaac PetscFunctionBegin; 1391db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 139258762b62SToby Isaac forest->setfromoptionscalled = PETSC_TRUE; 1393dd8e54a2SToby Isaac ierr = DMForestGetTopology(dm, &oldTopo);CHKERRQ(ierr); 1394a3eda1eaSToby Isaac ierr = PetscOptionsHead(PetscOptionsObject,"DMForest Options");CHKERRQ(ierr); 139556ba9f64SToby Isaac ierr = PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);CHKERRQ(ierr); 139656ba9f64SToby Isaac ierr = PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);CHKERRQ(ierr); 139756ba9f64SToby Isaac ierr = PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);CHKERRQ(ierr); 139856ba9f64SToby Isaac ierr = PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);CHKERRQ(ierr); 139956ba9f64SToby Isaac if ((PetscInt) flg1 + (PetscInt) flg2 + (PetscInt) flg3 + (PetscInt) flg4 > 1) { 140056ba9f64SToby Isaac SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Specify only one of -dm_forest_{topology,base_dm,coarse_forest,fine_forest}"); 1401dd8e54a2SToby Isaac } 140256ba9f64SToby Isaac if (flg1) { 140356ba9f64SToby Isaac ierr = DMForestSetTopology(dm,(DMForestTopology)stringBuffer);CHKERRQ(ierr); 140456ba9f64SToby Isaac ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr); 140520e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr); 140656ba9f64SToby Isaac } 140756ba9f64SToby Isaac if (flg2) { 1408dd8e54a2SToby Isaac DM base; 1409dd8e54a2SToby Isaac 1410dd8e54a2SToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),&base);CHKERRQ(ierr); 1411dd8e54a2SToby Isaac ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 1412dd8e54a2SToby Isaac ierr = DMLoad(base,viewer);CHKERRQ(ierr); 1413dd8e54a2SToby Isaac ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1414dd8e54a2SToby Isaac ierr = DMForestSetBaseDM(dm,base);CHKERRQ(ierr); 1415dd8e54a2SToby Isaac ierr = DMDestroy(&base);CHKERRQ(ierr); 141656ba9f64SToby Isaac ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr); 141720e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr); 1418dd8e54a2SToby Isaac } 141956ba9f64SToby Isaac if (flg3) { 1420dd8e54a2SToby Isaac DM coarse; 1421dd8e54a2SToby Isaac 1422dd8e54a2SToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),&coarse);CHKERRQ(ierr); 1423dd8e54a2SToby Isaac ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 1424dd8e54a2SToby Isaac ierr = DMLoad(coarse,viewer);CHKERRQ(ierr); 1425dd8e54a2SToby Isaac ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 142620e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,coarse);CHKERRQ(ierr); 1427dd8e54a2SToby Isaac ierr = DMDestroy(&coarse);CHKERRQ(ierr); 142856ba9f64SToby Isaac ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr); 142956ba9f64SToby Isaac ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr); 1430dd8e54a2SToby Isaac } 143156ba9f64SToby Isaac if (flg4) { 1432dd8e54a2SToby Isaac DM fine; 1433dd8e54a2SToby Isaac 1434dd8e54a2SToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),&fine);CHKERRQ(ierr); 1435dd8e54a2SToby Isaac ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 1436dd8e54a2SToby Isaac ierr = DMLoad(fine,viewer);CHKERRQ(ierr); 1437dd8e54a2SToby Isaac ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 143820e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,fine);CHKERRQ(ierr); 1439dd8e54a2SToby Isaac ierr = DMDestroy(&fine);CHKERRQ(ierr); 144056ba9f64SToby Isaac ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr); 144156ba9f64SToby Isaac ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr); 1442dd8e54a2SToby Isaac } 1443dd8e54a2SToby Isaac ierr = DMForestGetAdjacencyDimension(dm,&adjDim);CHKERRQ(ierr); 1444dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg);CHKERRQ(ierr); 1445dd8e54a2SToby Isaac if (flg) { 1446dd8e54a2SToby Isaac ierr = DMForestSetAdjacencyDimension(dm,adjDim);CHKERRQ(ierr); 1447dd8e54a2SToby Isaac } 1448dd8e54a2SToby Isaac else { 1449dd8e54a2SToby Isaac ierr = DMForestGetAdjacencyCodimension(dm,&adjCodim);CHKERRQ(ierr); 1450dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg);CHKERRQ(ierr); 1451dd8e54a2SToby Isaac if (flg) { 1452dd8e54a2SToby Isaac ierr = DMForestSetAdjacencyCodimension(dm,adjCodim);CHKERRQ(ierr); 1453dd8e54a2SToby Isaac } 1454dd8e54a2SToby Isaac } 1455dd8e54a2SToby Isaac ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr); 1456dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg);CHKERRQ(ierr); 1457dd8e54a2SToby Isaac if (flg) { 1458dd8e54a2SToby Isaac ierr = DMForestSetPartitionOverlap(dm,overlap);CHKERRQ(ierr); 1459dd8e54a2SToby Isaac } 1460a6121fbdSMatthew G. Knepley #if 0 1461a6121fbdSMatthew G. Knepley ierr = PetscOptionsInt("-dm_refine","equivalent to -dm_forest_set_minimum_refinement and -dm_forest_set_initial_refinement with the same value",NULL,minRefinement,&minRefinement,&flg);CHKERRQ(ierr); 1462a6121fbdSMatthew G. Knepley if (flg) { 1463a6121fbdSMatthew G. Knepley ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr); 1464a6121fbdSMatthew G. Knepley ierr = DMForestSetInitialRefinement(dm,minRefinement);CHKERRQ(ierr); 1465a6121fbdSMatthew G. Knepley } 1466a6121fbdSMatthew G. Knepley ierr = PetscOptionsInt("-dm_refine_hierarchy","equivalent to -dm_forest_set_minimum_refinement 0 and -dm_forest_set_initial_refinement",NULL,initRefinement,&initRefinement,&flg);CHKERRQ(ierr); 1467a6121fbdSMatthew G. Knepley if (flg) { 1468a6121fbdSMatthew G. Knepley ierr = DMForestSetMinimumRefinement(dm,0);CHKERRQ(ierr); 1469a6121fbdSMatthew G. Knepley ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr); 1470a6121fbdSMatthew G. Knepley } 1471a6121fbdSMatthew G. Knepley #endif 1472dd8e54a2SToby Isaac ierr = DMForestGetMinimumRefinement(dm,&minRefinement);CHKERRQ(ierr); 1473dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg);CHKERRQ(ierr); 1474dd8e54a2SToby Isaac if (flg) { 1475dd8e54a2SToby Isaac ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr); 1476db4d5e8cSToby Isaac } 147756ba9f64SToby Isaac ierr = DMForestGetInitialRefinement(dm,&initRefinement);CHKERRQ(ierr); 147856ba9f64SToby Isaac ierr = PetscOptionsInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg);CHKERRQ(ierr); 147956ba9f64SToby Isaac if (flg) { 148056ba9f64SToby Isaac ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr); 148156ba9f64SToby Isaac } 1482c7eeac06SToby Isaac ierr = DMForestGetMaximumRefinement(dm,&maxRefinement);CHKERRQ(ierr); 1483c7eeac06SToby Isaac ierr = PetscOptionsInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg);CHKERRQ(ierr); 1484c7eeac06SToby Isaac if (flg) { 1485c7eeac06SToby Isaac ierr = DMForestSetMaximumRefinement(dm,maxRefinement);CHKERRQ(ierr); 1486c7eeac06SToby Isaac } 1487c7eeac06SToby Isaac ierr = DMForestGetAdaptivityStrategy(dm,&adaptStrategy);CHKERRQ(ierr); 1488c7eeac06SToby Isaac ierr = PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,256,&flg);CHKERRQ(ierr); 1489c7eeac06SToby Isaac if (flg) { 1490c7eeac06SToby Isaac ierr = DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);CHKERRQ(ierr); 1491c7eeac06SToby Isaac } 1492c7eeac06SToby Isaac ierr = DMForestGetGradeFactor(dm,&grade);CHKERRQ(ierr); 1493c7eeac06SToby Isaac ierr = PetscOptionsInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg);CHKERRQ(ierr); 1494c7eeac06SToby Isaac if (flg) { 1495c7eeac06SToby Isaac ierr = DMForestSetGradeFactor(dm,grade);CHKERRQ(ierr); 1496c7eeac06SToby Isaac } 1497c7eeac06SToby Isaac ierr = DMForestGetCellWeightFactor(dm,&weightsFactor);CHKERRQ(ierr); 1498c7eeac06SToby Isaac ierr = PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);CHKERRQ(ierr); 1499c7eeac06SToby Isaac if (flg) { 1500c7eeac06SToby Isaac ierr = DMForestSetCellWeightFactor(dm,weightsFactor);CHKERRQ(ierr); 1501c7eeac06SToby Isaac } 1502db4d5e8cSToby Isaac ierr = PetscOptionsTail();CHKERRQ(ierr); 1503db4d5e8cSToby Isaac PetscFunctionReturn(0); 1504db4d5e8cSToby Isaac } 1505db4d5e8cSToby Isaac 1506db4d5e8cSToby Isaac #undef __FUNCT__ 1507d8984e3bSMatthew G. Knepley #define __FUNCT__ "DMCreateSubDM_Forest" 1508d8984e3bSMatthew G. Knepley PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm) 1509d8984e3bSMatthew G. Knepley { 1510d8984e3bSMatthew G. Knepley PetscErrorCode ierr; 1511d8984e3bSMatthew G. Knepley 1512d8984e3bSMatthew G. Knepley PetscFunctionBegin; 1513d8984e3bSMatthew G. Knepley if (subdm) {ierr = DMClone(dm, subdm);CHKERRQ(ierr);} 1514d8984e3bSMatthew G. Knepley ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr); 1515d8984e3bSMatthew G. Knepley PetscFunctionReturn(0); 1516d8984e3bSMatthew G. Knepley } 1517d8984e3bSMatthew G. Knepley 1518d8984e3bSMatthew G. Knepley #undef __FUNCT__ 15195421bac9SToby Isaac #define __FUNCT__ "DMRefine_Forest" 15205421bac9SToby Isaac PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined) 15215421bac9SToby Isaac { 15225421bac9SToby Isaac DMLabel refine; 15235421bac9SToby Isaac DM fineDM; 15245421bac9SToby Isaac PetscErrorCode ierr; 15255421bac9SToby Isaac 15265421bac9SToby Isaac PetscFunctionBegin; 15275421bac9SToby Isaac ierr = DMGetFineDM(dm,&fineDM);CHKERRQ(ierr); 15285421bac9SToby Isaac if (fineDM) { 15295421bac9SToby Isaac ierr = PetscObjectReference((PetscObject)fineDM);CHKERRQ(ierr); 15305421bac9SToby Isaac *dmRefined = fineDM; 15315421bac9SToby Isaac PetscFunctionReturn(0); 15325421bac9SToby Isaac } 15335421bac9SToby Isaac ierr = DMForestTemplate(dm,comm,dmRefined);CHKERRQ(ierr); 15345421bac9SToby Isaac ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr); 15355421bac9SToby Isaac if (!refine) { 15365421bac9SToby Isaac ierr = DMCreateLabel(dm,"refine");CHKERRQ(ierr); 15375421bac9SToby Isaac ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr); 15385421bac9SToby Isaac ierr = DMLabelSetDefaultValue(refine,DM_FOREST_REFINE);CHKERRQ(ierr); 15395421bac9SToby Isaac } 15405421bac9SToby Isaac ierr = DMForestSetAdaptivityLabel(*dmRefined,"refine");CHKERRQ(ierr); 15415421bac9SToby Isaac PetscFunctionReturn(0); 15425421bac9SToby Isaac } 15435421bac9SToby Isaac 15445421bac9SToby Isaac #undef __FUNCT__ 15455421bac9SToby Isaac #define __FUNCT__ "DMCoarsen_Forest" 15465421bac9SToby Isaac PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened) 15475421bac9SToby Isaac { 15485421bac9SToby Isaac DMLabel coarsen; 15495421bac9SToby Isaac DM coarseDM; 15505421bac9SToby Isaac PetscErrorCode ierr; 15515421bac9SToby Isaac 15525421bac9SToby Isaac PetscFunctionBegin; 15535421bac9SToby Isaac ierr = DMGetCoarseDM(dm,&coarseDM);CHKERRQ(ierr); 15545421bac9SToby Isaac if (coarseDM) { 15555421bac9SToby Isaac ierr = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr); 15565421bac9SToby Isaac *dmCoarsened = coarseDM; 15575421bac9SToby Isaac PetscFunctionReturn(0); 15585421bac9SToby Isaac } 15595421bac9SToby Isaac ierr = DMForestTemplate(dm,comm,dmCoarsened);CHKERRQ(ierr); 15605421bac9SToby Isaac ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr); 15615421bac9SToby Isaac if (!coarsen) { 15625421bac9SToby Isaac ierr = DMCreateLabel(dm,"coarsen");CHKERRQ(ierr); 15635421bac9SToby Isaac ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr); 15645421bac9SToby Isaac ierr = DMLabelSetDefaultValue(coarsen,DM_FOREST_COARSEN);CHKERRQ(ierr); 15655421bac9SToby Isaac } 15665421bac9SToby Isaac ierr = DMForestSetAdaptivityLabel(*dmCoarsened,"coarsen");CHKERRQ(ierr); 15675421bac9SToby Isaac PetscFunctionReturn(0); 15685421bac9SToby Isaac } 15695421bac9SToby Isaac 15705421bac9SToby Isaac #undef __FUNCT__ 1571d222f98bSToby Isaac #define __FUNCT__ "DMInitialize_Forest" 1572d222f98bSToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm) 1573d222f98bSToby Isaac { 1574d222f98bSToby Isaac PetscErrorCode ierr; 1575d222f98bSToby Isaac 1576d222f98bSToby Isaac PetscFunctionBegin; 1577d222f98bSToby Isaac ierr = PetscMemzero(dm->ops,sizeof(*(dm->ops)));CHKERRQ(ierr); 1578d222f98bSToby Isaac 1579d222f98bSToby Isaac dm->ops->clone = DMClone_Forest; 1580d222f98bSToby Isaac dm->ops->setfromoptions = DMSetFromOptions_Forest; 1581d222f98bSToby Isaac dm->ops->destroy = DMDestroy_Forest; 1582d8984e3bSMatthew G. Knepley dm->ops->createsubdm = DMCreateSubDM_Forest; 15835421bac9SToby Isaac dm->ops->refine = DMRefine_Forest; 15845421bac9SToby Isaac dm->ops->coarsen = DMCoarsen_Forest; 1585d222f98bSToby Isaac PetscFunctionReturn(0); 1586d222f98bSToby Isaac } 1587d222f98bSToby Isaac 15889be51f97SToby Isaac /*MC 15899be51f97SToby Isaac DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh. Forests usually have a base DM (see DMForestGetBaseDM()), from which it is refined. The refinement and partitioning of forests is considered immutable after DMSetUp() is called. To adapt a mesh, one should call DMForestTemplate() to create a new mesh that will default to being identical to it, specify how that mesh should differ, and then calling DMSetUp() on the new mesh. 15909be51f97SToby Isaac 15919be51f97SToby Isaac To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the previous mesh whose values indicate which cells should be refined (DM_FOREST_REFINE) or coarsened (DM_FOREST_COARSEN) and how (subtypes are free to allow additional values for things like anisotropic refinement). The name of the label should be given to the *new* mesh with DMForestSetAdaptivityLabel(). 15929be51f97SToby Isaac 15939be51f97SToby Isaac Level: advanced 15949be51f97SToby Isaac 15959be51f97SToby Isaac .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel() 15969be51f97SToby Isaac M*/ 15979be51f97SToby Isaac 1598d222f98bSToby Isaac #undef __FUNCT__ 1599db4d5e8cSToby Isaac #define __FUNCT__ "DMCreate_Forest" 1600db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm) 1601db4d5e8cSToby Isaac { 1602db4d5e8cSToby Isaac DM_Forest *forest; 1603db4d5e8cSToby Isaac PetscErrorCode ierr; 1604db4d5e8cSToby Isaac 1605db4d5e8cSToby Isaac PetscFunctionBegin; 1606db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1607db4d5e8cSToby Isaac ierr = PetscNewLog(dm,&forest);CHKERRQ(ierr); 1608db4d5e8cSToby Isaac dm->dim = 0; 1609db4d5e8cSToby Isaac dm->data = forest; 1610db4d5e8cSToby Isaac forest->refct = 1; 1611db4d5e8cSToby Isaac forest->data = NULL; 161258762b62SToby Isaac forest->setfromoptionscalled = PETSC_FALSE; 1613db4d5e8cSToby Isaac forest->topology = NULL; 1614db4d5e8cSToby Isaac forest->base = NULL; 1615db4d5e8cSToby Isaac forest->adjDim = PETSC_DEFAULT; 1616db4d5e8cSToby Isaac forest->overlap = PETSC_DEFAULT; 1617db4d5e8cSToby Isaac forest->minRefinement = PETSC_DEFAULT; 1618db4d5e8cSToby Isaac forest->maxRefinement = PETSC_DEFAULT; 161956ba9f64SToby Isaac forest->initRefinement = PETSC_DEFAULT; 1620c7eeac06SToby Isaac forest->cStart = PETSC_DETERMINE; 1621c7eeac06SToby Isaac forest->cEnd = PETSC_DETERMINE; 1622db4d5e8cSToby Isaac forest->cellSF = 0; 1623ebdf65a2SToby Isaac forest->adaptLabel = NULL; 1624db4d5e8cSToby Isaac forest->gradeFactor = 2; 1625db4d5e8cSToby Isaac forest->cellWeights = NULL; 1626db4d5e8cSToby Isaac forest->cellWeightsCopyMode = PETSC_USE_POINTER; 1627db4d5e8cSToby Isaac forest->weightsFactor = 1.; 1628db4d5e8cSToby Isaac forest->weightCapacity = 1.; 1629a73e2921SToby Isaac ierr = DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);CHKERRQ(ierr); 1630d222f98bSToby Isaac ierr = DMInitialize_Forest(dm);CHKERRQ(ierr); 1631db4d5e8cSToby Isaac PetscFunctionReturn(0); 1632db4d5e8cSToby Isaac } 1633db4d5e8cSToby Isaac 1634