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; 27f416af30SBarry Smith ierr = PetscFree(oldLink->name);CHKERRQ(ierr); 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); 2140f17b9e3SToby Isaac ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr); 2150f17b9e3SToby Isaac ierr = PetscSFDestroy(&forest->coarseToPreFine);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__ 48256c0450aSToby Isaac #define __FUNCT__ "DMForestGetAdaptivityPurpose" 48356c0450aSToby Isaac /*@ 48456c0450aSToby Isaac DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with 48556c0450aSToby Isaac DMForestSetAdaptivityForest()) for the purpose of refinement (DM_FOREST_REFINE), coarsening (DM_FOREST_COARSEN), or 48656c0450aSToby Isaac undefined (DM_FOREST_NONE). This only matters for the purposes of reference counting: during DMDestroy(), cyclic 48756c0450aSToby Isaac references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see 48856c0450aSToby Isaac DMSetFineDM()/DMSetCoarseDM()). If the purpose is not refinement or coarsening, and the user does not maintain a 48956c0450aSToby Isaac reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory 49056c0450aSToby Isaac leak. This method is used by subtypes of DMForest when automatically constructing mesh hierarchies. 49156c0450aSToby Isaac 49256c0450aSToby Isaac Not collective 49356c0450aSToby Isaac 49456c0450aSToby Isaac Input Parameter: 49556c0450aSToby Isaac . dm - the forest 49656c0450aSToby Isaac 49756c0450aSToby Isaac Output Parameter: 49856c0450aSToby Isaac . purpose - the adaptivity purpose (DM_FOREST_NONE/DM_FOREST_REFINE/DM_FOREST_COARSEN) 49956c0450aSToby Isaac 50056c0450aSToby Isaac Level: advanced 50156c0450aSToby Isaac 50256c0450aSToby Isaac .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest() 50356c0450aSToby Isaac @*/ 50456c0450aSToby Isaac PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMForestAdaptivityPurpose *purpose) 50556c0450aSToby Isaac { 50656c0450aSToby Isaac DM_Forest *forest; 50756c0450aSToby Isaac 50856c0450aSToby Isaac PetscFunctionBegin; 50956c0450aSToby Isaac forest = (DM_Forest *) dm->data; 51056c0450aSToby Isaac *purpose = forest->adaptPurpose; 51156c0450aSToby Isaac PetscFunctionReturn(0); 51256c0450aSToby Isaac } 51356c0450aSToby Isaac 51456c0450aSToby Isaac #undef __FUNCT__ 515dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetAdjacencyDimension" 5169be51f97SToby Isaac /*@ 5179be51f97SToby Isaac DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine 5189be51f97SToby Isaac cell adjacency (for the purposes of partitioning and overlap). 5199be51f97SToby Isaac 5209be51f97SToby Isaac Logically collective on dm 5219be51f97SToby Isaac 5229be51f97SToby Isaac Input Parameters: 5239be51f97SToby Isaac + dm - the forest 5249be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency) 5259be51f97SToby Isaac 5269be51f97SToby Isaac Level: intermediate 5279be51f97SToby Isaac 5289be51f97SToby Isaac .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap() 5299be51f97SToby Isaac @*/ 530dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim) 531dd8e54a2SToby Isaac { 532dd8e54a2SToby Isaac PetscInt dim; 533dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 534dd8e54a2SToby Isaac PetscErrorCode ierr; 535dd8e54a2SToby Isaac 536dd8e54a2SToby Isaac PetscFunctionBegin; 537dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 538ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup"); 539dd8e54a2SToby Isaac if (adjDim < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %d", adjDim); 540dd8e54a2SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 541dd8e54a2SToby Isaac if (adjDim > dim) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %d: %d", dim, adjDim); 542dd8e54a2SToby Isaac forest->adjDim = adjDim; 543dd8e54a2SToby Isaac PetscFunctionReturn(0); 544dd8e54a2SToby Isaac } 545dd8e54a2SToby Isaac 546dd8e54a2SToby Isaac #undef __FUNCT__ 547dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetAdjacencyCodimension" 5489be51f97SToby Isaac /*@ 5499be51f97SToby Isaac DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that, 5509be51f97SToby Isaac e.g., adjacency based on facets can be specified by codimension 1 in all cases) 5519be51f97SToby Isaac 5529be51f97SToby Isaac Logically collective on dm 5539be51f97SToby Isaac 5549be51f97SToby Isaac Input Parameters: 5559be51f97SToby Isaac + dm - the forest 5569be51f97SToby Isaac - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices 5579be51f97SToby Isaac 5589be51f97SToby Isaac Level: intermediate 5599be51f97SToby Isaac 5609be51f97SToby Isaac .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension() 5619be51f97SToby Isaac @*/ 562dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim) 563dd8e54a2SToby Isaac { 564dd8e54a2SToby Isaac PetscInt dim; 565dd8e54a2SToby Isaac PetscErrorCode ierr; 566dd8e54a2SToby Isaac 567dd8e54a2SToby Isaac PetscFunctionBegin; 568dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 569dd8e54a2SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 570dd8e54a2SToby Isaac ierr = DMForestSetAdjacencyDimension(dm,dim-adjCodim);CHKERRQ(ierr); 571dd8e54a2SToby Isaac PetscFunctionReturn(0); 572dd8e54a2SToby Isaac } 573dd8e54a2SToby Isaac 574dd8e54a2SToby Isaac #undef __FUNCT__ 575dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetAdjacencyDimension" 5769be51f97SToby Isaac /*@ 5779be51f97SToby Isaac DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the 5789be51f97SToby Isaac purposes of partitioning and overlap). 5799be51f97SToby Isaac 5809be51f97SToby Isaac Not collective 5819be51f97SToby Isaac 5829be51f97SToby Isaac Input Parameter: 5839be51f97SToby Isaac . dm - the forest 5849be51f97SToby Isaac 5859be51f97SToby Isaac Output Parameter: 5869be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency) 5879be51f97SToby Isaac 5889be51f97SToby Isaac Level: intermediate 5899be51f97SToby Isaac 5909be51f97SToby Isaac .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap() 5919be51f97SToby Isaac @*/ 592dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim) 593dd8e54a2SToby Isaac { 594dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 595dd8e54a2SToby Isaac 596dd8e54a2SToby Isaac PetscFunctionBegin; 597dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 598dd8e54a2SToby Isaac PetscValidIntPointer(adjDim,2); 599dd8e54a2SToby Isaac *adjDim = forest->adjDim; 600dd8e54a2SToby Isaac PetscFunctionReturn(0); 601dd8e54a2SToby Isaac } 602dd8e54a2SToby Isaac 603dd8e54a2SToby Isaac #undef __FUNCT__ 604dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetAdjacencyCodimension" 6059be51f97SToby Isaac /*@ 6069be51f97SToby Isaac DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that, 6079be51f97SToby Isaac e.g., adjacency based on facets can be specified by codimension 1 in all cases) 6089be51f97SToby Isaac 6099be51f97SToby Isaac Not collective 6109be51f97SToby Isaac 6119be51f97SToby Isaac Input Parameter: 6129be51f97SToby Isaac . dm - the forest 6139be51f97SToby Isaac 6149be51f97SToby Isaac Output Parameter: 6159be51f97SToby Isaac . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices 6169be51f97SToby Isaac 6179be51f97SToby Isaac Level: intermediate 6189be51f97SToby Isaac 6199be51f97SToby Isaac .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension() 6209be51f97SToby Isaac @*/ 621dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim) 622dd8e54a2SToby Isaac { 623dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 624dd8e54a2SToby Isaac PetscInt dim; 625dd8e54a2SToby Isaac PetscErrorCode ierr; 626dd8e54a2SToby Isaac 627dd8e54a2SToby Isaac PetscFunctionBegin; 628dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 629dd8e54a2SToby Isaac PetscValidIntPointer(adjCodim,2); 630dd8e54a2SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 631dd8e54a2SToby Isaac *adjCodim = dim - forest->adjDim; 632dd8e54a2SToby Isaac PetscFunctionReturn(0); 633dd8e54a2SToby Isaac } 634dd8e54a2SToby Isaac 635dd8e54a2SToby Isaac #undef __FUNCT__ 636ef51cf95SToby Isaac #define __FUNCT__ "DMForestSetPartitionOverlap" 6379be51f97SToby Isaac /*@ 6389be51f97SToby Isaac DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel 6399be51f97SToby Isaac partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding 6409be51f97SToby Isaac adjacent cells 6419be51f97SToby Isaac 6429be51f97SToby Isaac Logically collective on dm 6439be51f97SToby Isaac 6449be51f97SToby Isaac Input Parameters: 6459be51f97SToby Isaac + dm - the forest 6469be51f97SToby Isaac - overlap - default 0 6479be51f97SToby Isaac 6489be51f97SToby Isaac Level: intermediate 6499be51f97SToby Isaac 6509be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension() 6519be51f97SToby Isaac @*/ 652dd8e54a2SToby Isaac PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap) 653dd8e54a2SToby Isaac { 654dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 655dd8e54a2SToby Isaac 656dd8e54a2SToby Isaac PetscFunctionBegin; 657dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 658ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup"); 659dd8e54a2SToby Isaac if (overlap < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %d", overlap); 660dd8e54a2SToby Isaac forest->overlap = overlap; 661dd8e54a2SToby Isaac PetscFunctionReturn(0); 662dd8e54a2SToby Isaac } 663dd8e54a2SToby Isaac 664dd8e54a2SToby Isaac #undef __FUNCT__ 665dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetPartitionOverlap" 6669be51f97SToby Isaac /*@ 6679be51f97SToby Isaac DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values 6689be51f97SToby Isaac > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells 6699be51f97SToby Isaac 6709be51f97SToby Isaac Not collective 6719be51f97SToby Isaac 6729be51f97SToby Isaac Input Parameter: 6739be51f97SToby Isaac . dm - the forest 6749be51f97SToby Isaac 6759be51f97SToby Isaac Output Parameter: 6769be51f97SToby Isaac . overlap - default 0 6779be51f97SToby Isaac 6789be51f97SToby Isaac Level: intermediate 6799be51f97SToby Isaac 6809be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension() 6819be51f97SToby Isaac @*/ 682dd8e54a2SToby Isaac PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap) 683dd8e54a2SToby Isaac { 684dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 685dd8e54a2SToby Isaac 686dd8e54a2SToby Isaac PetscFunctionBegin; 687dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 688dd8e54a2SToby Isaac PetscValidIntPointer(overlap,2); 689dd8e54a2SToby Isaac *overlap = forest->overlap; 690dd8e54a2SToby Isaac PetscFunctionReturn(0); 691dd8e54a2SToby Isaac } 692dd8e54a2SToby Isaac 693dd8e54a2SToby Isaac #undef __FUNCT__ 694dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetMinimumRefinement" 6959be51f97SToby Isaac /*@ 6969be51f97SToby Isaac DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base 6979be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. If the forest is being created by coarsening a previous forest 6989be51f97SToby Isaac (see DMForestGetAdaptivityForest()) this limits the amount of coarsening. 6999be51f97SToby Isaac 7009be51f97SToby Isaac Logically collective on dm 7019be51f97SToby Isaac 7029be51f97SToby Isaac Input Parameters: 7039be51f97SToby Isaac + dm - the forest 7049be51f97SToby Isaac - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7059be51f97SToby Isaac 7069be51f97SToby Isaac Level: intermediate 7079be51f97SToby Isaac 7089be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 7099be51f97SToby Isaac @*/ 710dd8e54a2SToby Isaac PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement) 711dd8e54a2SToby Isaac { 712dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 713dd8e54a2SToby Isaac 714dd8e54a2SToby Isaac PetscFunctionBegin; 715dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 716ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup"); 717dd8e54a2SToby Isaac forest->minRefinement = minRefinement; 718dd8e54a2SToby Isaac PetscFunctionReturn(0); 719dd8e54a2SToby Isaac } 720dd8e54a2SToby Isaac 721dd8e54a2SToby Isaac #undef __FUNCT__ 722dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetMinimumRefinement" 7239be51f97SToby Isaac /*@ 7249be51f97SToby Isaac DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see 7259be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. If the forest is being created by coarsening a previous forest (see 7269be51f97SToby Isaac DMForestGetAdaptivityForest()), this limits the amount of coarsening. 7279be51f97SToby Isaac 7289be51f97SToby Isaac Not collective 7299be51f97SToby Isaac 7309be51f97SToby Isaac Input Parameter: 7319be51f97SToby Isaac . dm - the forest 7329be51f97SToby Isaac 7339be51f97SToby Isaac Output Parameter: 7349be51f97SToby Isaac . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7359be51f97SToby Isaac 7369be51f97SToby Isaac Level: intermediate 7379be51f97SToby Isaac 7389be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestGetMaximumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 7399be51f97SToby Isaac @*/ 740dd8e54a2SToby Isaac PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement) 741dd8e54a2SToby Isaac { 742dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 743dd8e54a2SToby Isaac 744dd8e54a2SToby Isaac PetscFunctionBegin; 745dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 746dd8e54a2SToby Isaac PetscValidIntPointer(minRefinement,2); 747dd8e54a2SToby Isaac *minRefinement = forest->minRefinement; 748dd8e54a2SToby Isaac PetscFunctionReturn(0); 749dd8e54a2SToby Isaac } 750dd8e54a2SToby Isaac 751dd8e54a2SToby Isaac #undef __FUNCT__ 75256ba9f64SToby Isaac #define __FUNCT__ "DMForestSetInitialRefinement" 7539be51f97SToby Isaac /*@ 7549be51f97SToby Isaac DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base 7559be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. 7569be51f97SToby Isaac 7579be51f97SToby Isaac Logically collective on dm 7589be51f97SToby Isaac 7599be51f97SToby Isaac Input Parameters: 7609be51f97SToby Isaac + dm - the forest 7619be51f97SToby Isaac - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7629be51f97SToby Isaac 7639be51f97SToby Isaac Level: intermediate 7649be51f97SToby Isaac 7659be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM() 7669be51f97SToby Isaac @*/ 76756ba9f64SToby Isaac PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement) 76856ba9f64SToby Isaac { 76956ba9f64SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 77056ba9f64SToby Isaac 77156ba9f64SToby Isaac PetscFunctionBegin; 77256ba9f64SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 77356ba9f64SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup"); 77456ba9f64SToby Isaac forest->initRefinement = initRefinement; 77556ba9f64SToby Isaac PetscFunctionReturn(0); 77656ba9f64SToby Isaac } 77756ba9f64SToby Isaac 77856ba9f64SToby Isaac #undef __FUNCT__ 77956ba9f64SToby Isaac #define __FUNCT__ "DMForestGetInitialRefinement" 7809be51f97SToby Isaac /*@ 7819be51f97SToby Isaac DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see 7829be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. 7839be51f97SToby Isaac 7849be51f97SToby Isaac Not collective 7859be51f97SToby Isaac 7869be51f97SToby Isaac Input Parameter: 7879be51f97SToby Isaac . dm - the forest 7889be51f97SToby Isaac 7899be51f97SToby Isaac Output Paramater: 7909be51f97SToby Isaac . initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7919be51f97SToby Isaac 7929be51f97SToby Isaac Level: intermediate 7939be51f97SToby Isaac 7949be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM() 7959be51f97SToby Isaac @*/ 79656ba9f64SToby Isaac PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement) 79756ba9f64SToby Isaac { 79856ba9f64SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 79956ba9f64SToby Isaac 80056ba9f64SToby Isaac PetscFunctionBegin; 80156ba9f64SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 80256ba9f64SToby Isaac PetscValidIntPointer(initRefinement,2); 80356ba9f64SToby Isaac *initRefinement = forest->initRefinement; 80456ba9f64SToby Isaac PetscFunctionReturn(0); 80556ba9f64SToby Isaac } 80656ba9f64SToby Isaac 80756ba9f64SToby Isaac #undef __FUNCT__ 808c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetMaximumRefinement" 8099be51f97SToby Isaac /*@ 8109be51f97SToby Isaac DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base 8119be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. If the forest is being created by refining a previous forest 8129be51f97SToby Isaac (see DMForestGetAdaptivityForest()), this limits the amount of refinement. 8139be51f97SToby Isaac 8149be51f97SToby Isaac Logically collective on dm 8159be51f97SToby Isaac 8169be51f97SToby Isaac Input Parameters: 8179be51f97SToby Isaac + dm - the forest 8189be51f97SToby Isaac - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 8199be51f97SToby Isaac 8209be51f97SToby Isaac Level: intermediate 8219be51f97SToby Isaac 8229be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityDM() 8239be51f97SToby Isaac @*/ 824c7eeac06SToby Isaac PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement) 825dd8e54a2SToby Isaac { 826dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 827dd8e54a2SToby Isaac 828dd8e54a2SToby Isaac PetscFunctionBegin; 829dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 830ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup"); 831c7eeac06SToby Isaac forest->maxRefinement = maxRefinement; 832dd8e54a2SToby Isaac PetscFunctionReturn(0); 833dd8e54a2SToby Isaac } 834dd8e54a2SToby Isaac 835dd8e54a2SToby Isaac #undef __FUNCT__ 836c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetMaximumRefinement" 8379be51f97SToby Isaac /*@ 8389be51f97SToby Isaac DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see 8399be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. If the forest is being created by refining a previous forest (see 8409be51f97SToby Isaac DMForestGetAdaptivityForest()), this limits the amount of refinement. 8419be51f97SToby Isaac 8429be51f97SToby Isaac Not collective 8439be51f97SToby Isaac 8449be51f97SToby Isaac Input Parameter: 8459be51f97SToby Isaac . dm - the forest 8469be51f97SToby Isaac 8479be51f97SToby Isaac Output Parameter: 8489be51f97SToby Isaac . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 8499be51f97SToby Isaac 8509be51f97SToby Isaac Level: intermediate 8519be51f97SToby Isaac 8529be51f97SToby Isaac .seealso: DMForestSetMaximumRefinement(), DMForestGetMinimumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 8539be51f97SToby Isaac @*/ 854c7eeac06SToby Isaac PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement) 855dd8e54a2SToby Isaac { 856dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 857dd8e54a2SToby Isaac 858dd8e54a2SToby Isaac PetscFunctionBegin; 859dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 860c7eeac06SToby Isaac PetscValidIntPointer(maxRefinement,2); 861c7eeac06SToby Isaac *maxRefinement = forest->maxRefinement; 862dd8e54a2SToby Isaac PetscFunctionReturn(0); 863dd8e54a2SToby Isaac } 864c7eeac06SToby Isaac 865c7eeac06SToby Isaac #undef __FUNCT__ 866c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityStrategy" 8679be51f97SToby Isaac /*@C 8689be51f97SToby Isaac DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes. 8699be51f97SToby Isaac Subtypes of DMForest may define their own strategies. Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree 8709be51f97SToby Isaac for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to 8719be51f97SToby Isaac specify refinement/coarsening. 8729be51f97SToby Isaac 8739be51f97SToby Isaac Logically collective on dm 8749be51f97SToby Isaac 8759be51f97SToby Isaac Input Parameters: 8769be51f97SToby Isaac + dm - the forest 8779be51f97SToby Isaac - adaptStrategy - default DMFORESTADAPTALL 8789be51f97SToby Isaac 8799be51f97SToby Isaac Level: advanced 8809be51f97SToby Isaac 8819be51f97SToby Isaac .seealso: DMForestGetAdaptivityStrategy() 8829be51f97SToby Isaac @*/ 883c7eeac06SToby Isaac PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy) 884c7eeac06SToby Isaac { 885c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 886c7eeac06SToby Isaac PetscErrorCode ierr; 887c7eeac06SToby Isaac 888c7eeac06SToby Isaac PetscFunctionBegin; 889c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 890c7eeac06SToby Isaac ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr); 891a73e2921SToby Isaac ierr = PetscStrallocpy((const char *) adaptStrategy,(char **)&forest->adaptStrategy);CHKERRQ(ierr); 892c7eeac06SToby Isaac PetscFunctionReturn(0); 893c7eeac06SToby Isaac } 894c7eeac06SToby Isaac 895c7eeac06SToby Isaac #undef __FUNCT__ 896c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityStrategy" 8979be51f97SToby Isaac /*@C 8989be51f97SToby Isaac DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes. Subtypes 8999be51f97SToby Isaac of DMForest may define their own strategies. Two default strategies are DMFORESTADAPTALL, which indicates that all 9009be51f97SToby Isaac processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only 9019be51f97SToby Isaac one process needs to specify refinement/coarsening. 9029be51f97SToby Isaac 9039be51f97SToby Isaac Not collective 9049be51f97SToby Isaac 9059be51f97SToby Isaac Input Parameter: 9069be51f97SToby Isaac . dm - the forest 9079be51f97SToby Isaac 9089be51f97SToby Isaac Output Parameter: 9099be51f97SToby Isaac . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL) 9109be51f97SToby Isaac 9119be51f97SToby Isaac Level: advanced 9129be51f97SToby Isaac 9139be51f97SToby Isaac .seealso: DMForestSetAdaptivityStrategy() 9149be51f97SToby Isaac @*/ 915c7eeac06SToby Isaac PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy) 916c7eeac06SToby Isaac { 917c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 918c7eeac06SToby Isaac 919c7eeac06SToby Isaac PetscFunctionBegin; 920c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 921c7eeac06SToby Isaac PetscValidPointer(adaptStrategy,2); 922c7eeac06SToby Isaac *adaptStrategy = forest->adaptStrategy; 923c7eeac06SToby Isaac PetscFunctionReturn(0); 924c7eeac06SToby Isaac } 925c7eeac06SToby Isaac 926c7eeac06SToby Isaac #undef __FUNCT__ 927bf9b5d84SToby Isaac #define __FUNCT__ "DMForestSetComputeAdaptivitySF" 928bf9b5d84SToby Isaac /*@ 929bf9b5d84SToby Isaac DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed 930bf9b5d84SToby 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(). 931bf9b5d84SToby Isaac 932bf9b5d84SToby Isaac Logically collective on dm 933bf9b5d84SToby Isaac 934bf9b5d84SToby Isaac Input Parameters: 935bf9b5d84SToby Isaac + dm - the post-adaptation forest 936bf9b5d84SToby Isaac - computeSF - default PETSC_TRUE 937bf9b5d84SToby Isaac 938bf9b5d84SToby Isaac Level: advanced 939bf9b5d84SToby Isaac 940bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF() 941bf9b5d84SToby Isaac @*/ 942bf9b5d84SToby Isaac PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF) 943bf9b5d84SToby Isaac { 944bf9b5d84SToby Isaac DM_Forest *forest; 945bf9b5d84SToby Isaac 946bf9b5d84SToby Isaac PetscFunctionBegin; 947bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 948bf9b5d84SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called"); 949bf9b5d84SToby Isaac forest = (DM_Forest *) dm->data; 950bf9b5d84SToby Isaac forest->computeAdaptSF = computeSF; 951bf9b5d84SToby Isaac PetscFunctionReturn(0); 952bf9b5d84SToby Isaac } 953bf9b5d84SToby Isaac 954bf9b5d84SToby Isaac #undef __FUNCT__ 955bf9b5d84SToby Isaac #define __FUNCT__ "DMForestGetComputeAdaptivitySF" 956bf9b5d84SToby Isaac /*@ 957bf9b5d84SToby Isaac DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the 958bf9b5d84SToby Isaac pre-adaptation forest to the post-adaptiation forest. After DMSetUp() is called, these transfer PetscSFs can be 959bf9b5d84SToby Isaac accessed with DMForestGetAdaptivitySF(). 960bf9b5d84SToby Isaac 961bf9b5d84SToby Isaac Not collective 962bf9b5d84SToby Isaac 963bf9b5d84SToby Isaac Input Parameter: 964bf9b5d84SToby Isaac . dm - the post-adaptation forest 965bf9b5d84SToby Isaac 966bf9b5d84SToby Isaac Output Parameter: 967bf9b5d84SToby Isaac . computeSF - default PETSC_TRUE 968bf9b5d84SToby Isaac 969bf9b5d84SToby Isaac Level: advanced 970bf9b5d84SToby Isaac 971bf9b5d84SToby Isaac .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF() 972bf9b5d84SToby Isaac @*/ 973bf9b5d84SToby Isaac PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF) 974bf9b5d84SToby Isaac { 975bf9b5d84SToby Isaac DM_Forest *forest; 976bf9b5d84SToby Isaac 977bf9b5d84SToby Isaac PetscFunctionBegin; 978bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 979bf9b5d84SToby Isaac forest = (DM_Forest *) dm->data; 980bf9b5d84SToby Isaac *computeSF = forest->computeAdaptSF; 981bf9b5d84SToby Isaac PetscFunctionReturn(0); 982bf9b5d84SToby Isaac } 983bf9b5d84SToby Isaac 984bf9b5d84SToby Isaac #undef __FUNCT__ 985bf9b5d84SToby Isaac #define __FUNCT__ "DMForestGetAdaptivitySF" 986bf9b5d84SToby Isaac /*@ 987bf9b5d84SToby Isaac DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest. 988bf9b5d84SToby Isaac Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be 989bf9b5d84SToby Isaac some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa. Therefore there are two 990bf9b5d84SToby Isaac PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates 991bf9b5d84SToby Isaac pre-adaptation fine cells to post-adaptation coarse cells. 992bf9b5d84SToby Isaac 993bf9b5d84SToby Isaac Not collective 994bf9b5d84SToby Isaac 995bf9b5d84SToby Isaac Input Parameter: 996bf9b5d84SToby Isaac dm - the post-adaptation forest 997bf9b5d84SToby Isaac 998bf9b5d84SToby Isaac Output Parameter: 9990f17b9e3SToby Isaac preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post- 10000f17b9e3SToby Isaac coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre- 1001bf9b5d84SToby Isaac 1002bf9b5d84SToby Isaac Level: advanced 1003bf9b5d84SToby Isaac 1004bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF() 1005bf9b5d84SToby Isaac @*/ 10060f17b9e3SToby Isaac PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine) 1007bf9b5d84SToby Isaac { 1008bf9b5d84SToby Isaac DM_Forest *forest; 1009bf9b5d84SToby Isaac PetscErrorCode ierr; 1010bf9b5d84SToby Isaac 1011bf9b5d84SToby Isaac PetscFunctionBegin; 1012bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1013bf9b5d84SToby Isaac ierr = DMSetUp(dm);CHKERRQ(ierr); 1014bf9b5d84SToby Isaac forest = (DM_Forest *) dm->data; 10150f17b9e3SToby Isaac if (preCoarseToFine) { 10160f17b9e3SToby Isaac *preCoarseToFine = forest->preCoarseToFine; 1017bf9b5d84SToby Isaac } 10180f17b9e3SToby Isaac if (coarseToPreFine) { 10190f17b9e3SToby Isaac *coarseToPreFine = forest->coarseToPreFine; 1020bf9b5d84SToby Isaac } 1021bf9b5d84SToby Isaac PetscFunctionReturn(0); 1022bf9b5d84SToby Isaac } 1023bf9b5d84SToby Isaac 1024bf9b5d84SToby Isaac #undef __FUNCT__ 1025c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetGradeFactor" 10269be51f97SToby Isaac /*@ 10279be51f97SToby Isaac DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to 10289be51f97SToby Isaac indicate that the diameter of neighboring cells should differ by at most a factor of 2. Subtypes of DMForest may 10299be51f97SToby Isaac only support one particular choice of grading factor. 10309be51f97SToby Isaac 10319be51f97SToby Isaac Logically collective on dm 10329be51f97SToby Isaac 10339be51f97SToby Isaac Input Parameters: 10349be51f97SToby Isaac + dm - the forest 10359be51f97SToby Isaac - grade - the grading factor 10369be51f97SToby Isaac 10379be51f97SToby Isaac Level: advanced 10389be51f97SToby Isaac 10399be51f97SToby Isaac .seealso: DMForestGetGradeFactor() 10409be51f97SToby Isaac @*/ 1041c7eeac06SToby Isaac PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade) 1042c7eeac06SToby Isaac { 1043c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1044c7eeac06SToby Isaac 1045c7eeac06SToby Isaac PetscFunctionBegin; 1046c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1047ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup"); 1048c7eeac06SToby Isaac forest->gradeFactor = grade; 1049c7eeac06SToby Isaac PetscFunctionReturn(0); 1050c7eeac06SToby Isaac } 1051c7eeac06SToby Isaac 1052c7eeac06SToby Isaac #undef __FUNCT__ 1053c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetGradeFactor" 10549be51f97SToby Isaac /*@ 10559be51f97SToby Isaac DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of 10569be51f97SToby Isaac neighboring cells should differ by at most a factor of 2. Subtypes of DMForest may only support one particular 10579be51f97SToby Isaac choice of grading factor. 10589be51f97SToby Isaac 10599be51f97SToby Isaac Not collective 10609be51f97SToby Isaac 10619be51f97SToby Isaac Input Parameter: 10629be51f97SToby Isaac . dm - the forest 10639be51f97SToby Isaac 10649be51f97SToby Isaac Output Parameter: 10659be51f97SToby Isaac . grade - the grading factor 10669be51f97SToby Isaac 10679be51f97SToby Isaac Level: advanced 10689be51f97SToby Isaac 10699be51f97SToby Isaac .seealso: DMForestSetGradeFactor() 10709be51f97SToby Isaac @*/ 1071c7eeac06SToby Isaac PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade) 1072c7eeac06SToby Isaac { 1073c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1074c7eeac06SToby Isaac 1075c7eeac06SToby Isaac PetscFunctionBegin; 1076c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1077c7eeac06SToby Isaac PetscValidIntPointer(grade,2); 1078c7eeac06SToby Isaac *grade = forest->gradeFactor; 1079c7eeac06SToby Isaac PetscFunctionReturn(0); 1080c7eeac06SToby Isaac } 1081c7eeac06SToby Isaac 1082c7eeac06SToby Isaac #undef __FUNCT__ 1083ef51cf95SToby Isaac #define __FUNCT__ "DMForestSetCellWeightFactor" 10849be51f97SToby Isaac /*@ 10859be51f97SToby Isaac DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes 10869be51f97SToby Isaac the cell weight (see DMForestSetCellWeights()) when calculating partitions. The final weight of a cell will be 10879be51f97SToby Isaac (cellWeight) * (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on 10889be51f97SToby Isaac its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the 10899be51f97SToby Isaac computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement. 10909be51f97SToby Isaac 10919be51f97SToby Isaac Logically collective on dm 10929be51f97SToby Isaac 10939be51f97SToby Isaac Input Parameters: 10949be51f97SToby Isaac + dm - the forest 10959be51f97SToby Isaac - weightsFactors - default 1. 10969be51f97SToby Isaac 10979be51f97SToby Isaac Level: advanced 10989be51f97SToby Isaac 10999be51f97SToby Isaac .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights() 11009be51f97SToby Isaac @*/ 1101ef51cf95SToby Isaac PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor) 1102c7eeac06SToby Isaac { 1103c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1104c7eeac06SToby Isaac 1105c7eeac06SToby Isaac PetscFunctionBegin; 1106c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1107ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup"); 1108c7eeac06SToby Isaac forest->weightsFactor = weightsFactor; 1109c7eeac06SToby Isaac PetscFunctionReturn(0); 1110c7eeac06SToby Isaac } 1111c7eeac06SToby Isaac 1112c7eeac06SToby Isaac #undef __FUNCT__ 1113ef51cf95SToby Isaac #define __FUNCT__ "DMForestGetCellWeightFactor" 11149be51f97SToby Isaac /*@ 11159be51f97SToby Isaac DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see 11169be51f97SToby Isaac DMForestSetCellWeights()) when calculating partitions. The final weight of a cell will be (cellWeight) * 11179be51f97SToby Isaac (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on its level; a 11189be51f97SToby Isaac factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation 11199be51f97SToby Isaac associated with a cell is multiplied by a factor of 2 for each additional level of refinement. 11209be51f97SToby Isaac 11219be51f97SToby Isaac Not collective 11229be51f97SToby Isaac 11239be51f97SToby Isaac Input Parameter: 11249be51f97SToby Isaac . dm - the forest 11259be51f97SToby Isaac 11269be51f97SToby Isaac Output Parameter: 11279be51f97SToby Isaac . weightsFactors - default 1. 11289be51f97SToby Isaac 11299be51f97SToby Isaac Level: advanced 11309be51f97SToby Isaac 11319be51f97SToby Isaac .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights() 11329be51f97SToby Isaac @*/ 1133ef51cf95SToby Isaac PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor) 1134c7eeac06SToby Isaac { 1135c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1136c7eeac06SToby Isaac 1137c7eeac06SToby Isaac PetscFunctionBegin; 1138c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1139c7eeac06SToby Isaac PetscValidRealPointer(weightsFactor,2); 1140c7eeac06SToby Isaac *weightsFactor = forest->weightsFactor; 1141c7eeac06SToby Isaac PetscFunctionReturn(0); 1142c7eeac06SToby Isaac } 1143c7eeac06SToby Isaac 1144c7eeac06SToby Isaac #undef __FUNCT__ 1145c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellChart" 11469be51f97SToby Isaac /*@ 11479be51f97SToby Isaac DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process 11489be51f97SToby Isaac 11499be51f97SToby Isaac Not collective 11509be51f97SToby Isaac 11519be51f97SToby Isaac Input Parameter: 11529be51f97SToby Isaac . dm - the forest 11539be51f97SToby Isaac 11549be51f97SToby Isaac Output Parameters: 11559be51f97SToby Isaac + cStart - the first cell on this process 11569be51f97SToby Isaac - cEnd - one after the final cell on this process 11579be51f97SToby Isaac 1158*1a244344SSatish Balay Level: intermediate 11599be51f97SToby Isaac 11609be51f97SToby Isaac .seealso: DMForestGetCellSF() 11619be51f97SToby Isaac @*/ 1162c7eeac06SToby Isaac PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd) 1163c7eeac06SToby Isaac { 1164c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1165c7eeac06SToby Isaac PetscErrorCode ierr; 1166c7eeac06SToby Isaac 1167c7eeac06SToby Isaac PetscFunctionBegin; 1168c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1169c7eeac06SToby Isaac PetscValidIntPointer(cStart,2); 1170c7eeac06SToby Isaac PetscValidIntPointer(cEnd,2); 1171c7eeac06SToby Isaac if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) { 1172c7eeac06SToby Isaac ierr = forest->createcellchart(dm,&forest->cStart,&forest->cEnd);CHKERRQ(ierr); 1173c7eeac06SToby Isaac } 1174c7eeac06SToby Isaac *cStart = forest->cStart; 1175c7eeac06SToby Isaac *cEnd = forest->cEnd; 1176c7eeac06SToby Isaac PetscFunctionReturn(0); 1177c7eeac06SToby Isaac } 1178c7eeac06SToby Isaac 1179c7eeac06SToby Isaac #undef __FUNCT__ 1180c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellSF" 11819be51f97SToby Isaac /*@ 11829be51f97SToby Isaac DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes 11839be51f97SToby Isaac 11849be51f97SToby Isaac Not collective 11859be51f97SToby Isaac 11869be51f97SToby Isaac Input Parameter: 11879be51f97SToby Isaac . dm - the forest 11889be51f97SToby Isaac 11899be51f97SToby Isaac Output Parameter: 11909be51f97SToby Isaac . cellSF - the PetscSF 11919be51f97SToby Isaac 1192*1a244344SSatish Balay Level: intermediate 11939be51f97SToby Isaac 11949be51f97SToby Isaac .seealso: DMForestGetCellChart() 11959be51f97SToby Isaac @*/ 1196c7eeac06SToby Isaac PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF) 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); 1203c7eeac06SToby Isaac PetscValidPointer(cellSF,2); 1204c7eeac06SToby Isaac if ((!forest->cellSF) && forest->createcellsf) { 1205c7eeac06SToby Isaac ierr = forest->createcellsf(dm,&forest->cellSF);CHKERRQ(ierr); 1206c7eeac06SToby Isaac } 1207c7eeac06SToby Isaac *cellSF = forest->cellSF; 1208c7eeac06SToby Isaac PetscFunctionReturn(0); 1209c7eeac06SToby Isaac } 1210c7eeac06SToby Isaac 1211c7eeac06SToby Isaac #undef __FUNCT__ 1212ebdf65a2SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityLabel" 12139be51f97SToby Isaac /*@C 12149be51f97SToby Isaac DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see 12159be51f97SToby Isaac DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination). The 12169be51f97SToby Isaac interpretation of the label values is up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and 12179be51f97SToby Isaac DM_FOREST_COARSEN have been reserved as choices that should be accepted by all subtypes. 12189be51f97SToby Isaac 12199be51f97SToby Isaac Logically collective on dm 12209be51f97SToby Isaac 12219be51f97SToby Isaac Input Parameters: 12229be51f97SToby Isaac - dm - the forest 12239be51f97SToby Isaac + adaptLabel - the name of the label in the pre-adaptation forest 12249be51f97SToby Isaac 12259be51f97SToby Isaac Level: intermediate 12269be51f97SToby Isaac 12279be51f97SToby Isaac .seealso DMForestGetAdaptivityLabel() 12289be51f97SToby Isaac @*/ 1229ebdf65a2SToby Isaac PetscErrorCode DMForestSetAdaptivityLabel(DM dm, const char * adaptLabel) 1230c7eeac06SToby Isaac { 1231c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1232c7eeac06SToby Isaac PetscErrorCode ierr; 1233c7eeac06SToby Isaac 1234c7eeac06SToby Isaac PetscFunctionBegin; 1235c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1236ebdf65a2SToby Isaac ierr = PetscFree(forest->adaptLabel);CHKERRQ(ierr); 1237ebdf65a2SToby Isaac ierr = PetscStrallocpy(adaptLabel,&forest->adaptLabel);CHKERRQ(ierr); 1238c7eeac06SToby Isaac PetscFunctionReturn(0); 1239c7eeac06SToby Isaac } 1240c7eeac06SToby Isaac 1241c7eeac06SToby Isaac #undef __FUNCT__ 1242ebdf65a2SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityLabel" 12439be51f97SToby Isaac /*@C 12449be51f97SToby Isaac DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that 12459be51f97SToby Isaac holds the adaptation flags (refinement, coarsening, or some combination). The interpretation of the label values is 12469be51f97SToby Isaac up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and DM_FOREST_COARSEN have been reserved as 12479be51f97SToby Isaac choices that should be accepted by all subtypes. 12489be51f97SToby Isaac 12499be51f97SToby Isaac Not collective 12509be51f97SToby Isaac 12519be51f97SToby Isaac Input Parameter: 12529be51f97SToby Isaac . dm - the forest 12539be51f97SToby Isaac 12549be51f97SToby Isaac Output Parameter: 12559be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest 12569be51f97SToby Isaac 12579be51f97SToby Isaac Level: intermediate 12589be51f97SToby Isaac 12599be51f97SToby Isaac .seealso DMForestSetAdaptivityLabel() 12609be51f97SToby Isaac @*/ 1261ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityLabel(DM dm, const char ** adaptLabel) 1262c7eeac06SToby Isaac { 1263c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1264c7eeac06SToby Isaac 1265c7eeac06SToby Isaac PetscFunctionBegin; 1266c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1267ba936b91SToby Isaac *adaptLabel = forest->adaptLabel; 1268c7eeac06SToby Isaac PetscFunctionReturn(0); 1269c7eeac06SToby Isaac } 1270c7eeac06SToby Isaac 1271c7eeac06SToby Isaac #undef __FUNCT__ 1272c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetCellWeights" 12739be51f97SToby Isaac /*@ 12749be51f97SToby Isaac DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current 12759be51f97SToby Isaac process: weights are used to determine parallel partitioning. Partitions will be created so that each process's 12769be51f97SToby Isaac ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight 12779be51f97SToby Isaac of 1. 12789be51f97SToby Isaac 12799be51f97SToby Isaac Logically collective on dm 12809be51f97SToby Isaac 12819be51f97SToby Isaac Input Parameters: 12829be51f97SToby Isaac + dm - the forest 12839be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1. 12849be51f97SToby Isaac - copyMode - how weights should reference weights 12859be51f97SToby Isaac 12869be51f97SToby Isaac Level: advanced 12879be51f97SToby Isaac 12889be51f97SToby Isaac .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity() 12899be51f97SToby Isaac @*/ 1290c7eeac06SToby Isaac PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode) 1291c7eeac06SToby Isaac { 1292c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1293c7eeac06SToby Isaac PetscInt cStart, cEnd; 1294c7eeac06SToby Isaac PetscErrorCode ierr; 1295c7eeac06SToby Isaac 1296c7eeac06SToby Isaac PetscFunctionBegin; 1297c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1298c7eeac06SToby Isaac ierr = DMForestGetCellChart(dm,&cStart,&cEnd);CHKERRQ(ierr); 1299c7eeac06SToby Isaac if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd); 1300c7eeac06SToby Isaac if (copyMode == PETSC_COPY_VALUES) { 1301c7eeac06SToby Isaac if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) { 1302c7eeac06SToby Isaac ierr = PetscMalloc1(cEnd-cStart,&forest->cellWeights);CHKERRQ(ierr); 1303c7eeac06SToby Isaac } 1304c7eeac06SToby Isaac ierr = PetscMemcpy(forest->cellWeights,weights,(cEnd-cStart)*sizeof(*weights));CHKERRQ(ierr); 1305c7eeac06SToby Isaac forest->cellWeightsCopyMode = PETSC_OWN_POINTER; 1306c7eeac06SToby Isaac PetscFunctionReturn(0); 1307c7eeac06SToby Isaac } 1308c7eeac06SToby Isaac if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) { 1309c7eeac06SToby Isaac ierr = PetscFree(forest->cellWeights);CHKERRQ(ierr); 1310c7eeac06SToby Isaac } 1311c7eeac06SToby Isaac forest->cellWeights = weights; 1312c7eeac06SToby Isaac forest->cellWeightsCopyMode = copyMode; 1313c7eeac06SToby Isaac PetscFunctionReturn(0); 1314c7eeac06SToby Isaac } 1315c7eeac06SToby Isaac 1316c7eeac06SToby Isaac #undef __FUNCT__ 1317c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellWeights" 13189be51f97SToby Isaac /*@ 13199be51f97SToby Isaac DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current 13209be51f97SToby Isaac process: weights are used to determine parallel partitioning. Partitions will be created so that each process's 13219be51f97SToby Isaac ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight 13229be51f97SToby Isaac of 1. 13239be51f97SToby Isaac 13249be51f97SToby Isaac Not collective 13259be51f97SToby Isaac 13269be51f97SToby Isaac Input Parameter: 13279be51f97SToby Isaac . dm - the forest 13289be51f97SToby Isaac 13299be51f97SToby Isaac Output Parameter: 13309be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1. 13319be51f97SToby Isaac 13329be51f97SToby Isaac Level: advanced 13339be51f97SToby Isaac 13349be51f97SToby Isaac .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity() 13359be51f97SToby Isaac @*/ 1336c7eeac06SToby Isaac PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights) 1337c7eeac06SToby Isaac { 1338c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1339c7eeac06SToby Isaac 1340c7eeac06SToby Isaac PetscFunctionBegin; 1341c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1342c7eeac06SToby Isaac PetscValidPointer(weights,2); 1343c7eeac06SToby Isaac *weights = forest->cellWeights; 1344c7eeac06SToby Isaac PetscFunctionReturn(0); 1345c7eeac06SToby Isaac } 1346c7eeac06SToby Isaac 1347c7eeac06SToby Isaac #undef __FUNCT__ 1348c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetWeightCapacity" 13499be51f97SToby Isaac /*@ 13509be51f97SToby Isaac DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning 13519be51f97SToby Isaac a pre-adaptation forest (see DMForestGetAdaptivityForest()). After partitioning, the ratio of the weight of each 13529be51f97SToby Isaac process's cells to the process's capacity will be roughly equal for all processes. A capacity of 0 indicates that 13539be51f97SToby Isaac the current process should not have any cells after repartitioning. 13549be51f97SToby Isaac 13559be51f97SToby Isaac Logically Collective on dm 13569be51f97SToby Isaac 13579be51f97SToby Isaac Input parameters: 13589be51f97SToby Isaac + dm - the forest 13599be51f97SToby Isaac - capacity - this process's capacity 13609be51f97SToby Isaac 13619be51f97SToby Isaac Level: advanced 13629be51f97SToby Isaac 13639be51f97SToby Isaac .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor() 13649be51f97SToby Isaac @*/ 1365c7eeac06SToby Isaac PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity) 1366c7eeac06SToby Isaac { 1367c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1368c7eeac06SToby Isaac 1369c7eeac06SToby Isaac PetscFunctionBegin; 1370c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1371ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup"); 1372c7eeac06SToby Isaac if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity); 1373c7eeac06SToby Isaac forest->weightCapacity = capacity; 1374c7eeac06SToby Isaac PetscFunctionReturn(0); 1375c7eeac06SToby Isaac } 1376c7eeac06SToby Isaac 1377c7eeac06SToby Isaac #undef __FUNCT__ 1378c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetWeightCapacity" 13799be51f97SToby Isaac /*@ 13809be51f97SToby Isaac DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see 13819be51f97SToby Isaac DMForestGetAdaptivityForest()). After partitioning, the ratio of the weight of each process's cells to the 13829be51f97SToby Isaac process's capacity will be roughly equal for all processes. A capacity of 0 indicates that the current process 13839be51f97SToby Isaac should not have any cells after repartitioning. 13849be51f97SToby Isaac 13859be51f97SToby Isaac Not collective 13869be51f97SToby Isaac 13879be51f97SToby Isaac Input parameter: 13889be51f97SToby Isaac . dm - the forest 13899be51f97SToby Isaac 13909be51f97SToby Isaac Output parameter: 13919be51f97SToby Isaac . capacity - this process's capacity 13929be51f97SToby Isaac 13939be51f97SToby Isaac Level: advanced 13949be51f97SToby Isaac 13959be51f97SToby Isaac .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor() 13969be51f97SToby Isaac @*/ 1397c7eeac06SToby Isaac PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity) 1398c7eeac06SToby Isaac { 1399c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 1400c7eeac06SToby Isaac 1401c7eeac06SToby Isaac PetscFunctionBegin; 1402c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1403c7eeac06SToby Isaac PetscValidRealPointer(capacity,2); 1404c7eeac06SToby Isaac *capacity = forest->weightCapacity; 1405c7eeac06SToby Isaac PetscFunctionReturn(0); 1406c7eeac06SToby Isaac } 1407c7eeac06SToby Isaac 1408dd8e54a2SToby Isaac #undef __FUNCT__ 1409db4d5e8cSToby Isaac #define __FUNCT__ "DMSetFromOptions_Forest" 14100709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm) 1411db4d5e8cSToby Isaac { 1412db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest *) dm->data; 141356ba9f64SToby Isaac PetscBool flg, flg1, flg2, flg3, flg4; 1414dd8e54a2SToby Isaac DMForestTopology oldTopo; 1415c7eeac06SToby Isaac char stringBuffer[256]; 1416dd8e54a2SToby Isaac PetscViewer viewer; 1417dd8e54a2SToby Isaac PetscViewerFormat format; 141856ba9f64SToby Isaac PetscInt adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade; 1419c7eeac06SToby Isaac PetscReal weightsFactor; 1420c7eeac06SToby Isaac DMForestAdaptivityStrategy adaptStrategy; 1421db4d5e8cSToby Isaac PetscErrorCode ierr; 1422db4d5e8cSToby Isaac 1423db4d5e8cSToby Isaac PetscFunctionBegin; 1424db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 142558762b62SToby Isaac forest->setfromoptionscalled = PETSC_TRUE; 1426dd8e54a2SToby Isaac ierr = DMForestGetTopology(dm, &oldTopo);CHKERRQ(ierr); 1427a3eda1eaSToby Isaac ierr = PetscOptionsHead(PetscOptionsObject,"DMForest Options");CHKERRQ(ierr); 142856ba9f64SToby Isaac ierr = PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);CHKERRQ(ierr); 142956ba9f64SToby Isaac ierr = PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);CHKERRQ(ierr); 143056ba9f64SToby Isaac ierr = PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);CHKERRQ(ierr); 143156ba9f64SToby Isaac ierr = PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);CHKERRQ(ierr); 143256ba9f64SToby Isaac if ((PetscInt) flg1 + (PetscInt) flg2 + (PetscInt) flg3 + (PetscInt) flg4 > 1) { 143356ba9f64SToby Isaac SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Specify only one of -dm_forest_{topology,base_dm,coarse_forest,fine_forest}"); 1434dd8e54a2SToby Isaac } 143556ba9f64SToby Isaac if (flg1) { 143656ba9f64SToby Isaac ierr = DMForestSetTopology(dm,(DMForestTopology)stringBuffer);CHKERRQ(ierr); 143756ba9f64SToby Isaac ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr); 143820e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr); 143956ba9f64SToby Isaac } 144056ba9f64SToby Isaac if (flg2) { 1441dd8e54a2SToby Isaac DM base; 1442dd8e54a2SToby Isaac 1443dd8e54a2SToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),&base);CHKERRQ(ierr); 1444dd8e54a2SToby Isaac ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 1445dd8e54a2SToby Isaac ierr = DMLoad(base,viewer);CHKERRQ(ierr); 1446dd8e54a2SToby Isaac ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1447dd8e54a2SToby Isaac ierr = DMForestSetBaseDM(dm,base);CHKERRQ(ierr); 1448dd8e54a2SToby Isaac ierr = DMDestroy(&base);CHKERRQ(ierr); 144956ba9f64SToby Isaac ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr); 145020e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr); 1451dd8e54a2SToby Isaac } 145256ba9f64SToby Isaac if (flg3) { 1453dd8e54a2SToby Isaac DM coarse; 1454dd8e54a2SToby Isaac 1455dd8e54a2SToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),&coarse);CHKERRQ(ierr); 1456dd8e54a2SToby Isaac ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 1457dd8e54a2SToby Isaac ierr = DMLoad(coarse,viewer);CHKERRQ(ierr); 1458dd8e54a2SToby Isaac ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 145920e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,coarse);CHKERRQ(ierr); 1460dd8e54a2SToby Isaac ierr = DMDestroy(&coarse);CHKERRQ(ierr); 146156ba9f64SToby Isaac ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr); 146256ba9f64SToby Isaac ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr); 1463dd8e54a2SToby Isaac } 146456ba9f64SToby Isaac if (flg4) { 1465dd8e54a2SToby Isaac DM fine; 1466dd8e54a2SToby Isaac 1467dd8e54a2SToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),&fine);CHKERRQ(ierr); 1468dd8e54a2SToby Isaac ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 1469dd8e54a2SToby Isaac ierr = DMLoad(fine,viewer);CHKERRQ(ierr); 1470dd8e54a2SToby Isaac ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 147120e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,fine);CHKERRQ(ierr); 1472dd8e54a2SToby Isaac ierr = DMDestroy(&fine);CHKERRQ(ierr); 147356ba9f64SToby Isaac ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr); 147456ba9f64SToby Isaac ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr); 1475dd8e54a2SToby Isaac } 1476dd8e54a2SToby Isaac ierr = DMForestGetAdjacencyDimension(dm,&adjDim);CHKERRQ(ierr); 1477dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg);CHKERRQ(ierr); 1478dd8e54a2SToby Isaac if (flg) { 1479dd8e54a2SToby Isaac ierr = DMForestSetAdjacencyDimension(dm,adjDim);CHKERRQ(ierr); 1480dd8e54a2SToby Isaac } 1481dd8e54a2SToby Isaac else { 1482dd8e54a2SToby Isaac ierr = DMForestGetAdjacencyCodimension(dm,&adjCodim);CHKERRQ(ierr); 1483dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg);CHKERRQ(ierr); 1484dd8e54a2SToby Isaac if (flg) { 1485dd8e54a2SToby Isaac ierr = DMForestSetAdjacencyCodimension(dm,adjCodim);CHKERRQ(ierr); 1486dd8e54a2SToby Isaac } 1487dd8e54a2SToby Isaac } 1488dd8e54a2SToby Isaac ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr); 1489dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg);CHKERRQ(ierr); 1490dd8e54a2SToby Isaac if (flg) { 1491dd8e54a2SToby Isaac ierr = DMForestSetPartitionOverlap(dm,overlap);CHKERRQ(ierr); 1492dd8e54a2SToby Isaac } 1493a6121fbdSMatthew G. Knepley #if 0 1494a6121fbdSMatthew 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); 1495a6121fbdSMatthew G. Knepley if (flg) { 1496a6121fbdSMatthew G. Knepley ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr); 1497a6121fbdSMatthew G. Knepley ierr = DMForestSetInitialRefinement(dm,minRefinement);CHKERRQ(ierr); 1498a6121fbdSMatthew G. Knepley } 1499a6121fbdSMatthew 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); 1500a6121fbdSMatthew G. Knepley if (flg) { 1501a6121fbdSMatthew G. Knepley ierr = DMForestSetMinimumRefinement(dm,0);CHKERRQ(ierr); 1502a6121fbdSMatthew G. Knepley ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr); 1503a6121fbdSMatthew G. Knepley } 1504a6121fbdSMatthew G. Knepley #endif 1505dd8e54a2SToby Isaac ierr = DMForestGetMinimumRefinement(dm,&minRefinement);CHKERRQ(ierr); 1506dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg);CHKERRQ(ierr); 1507dd8e54a2SToby Isaac if (flg) { 1508dd8e54a2SToby Isaac ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr); 1509db4d5e8cSToby Isaac } 151056ba9f64SToby Isaac ierr = DMForestGetInitialRefinement(dm,&initRefinement);CHKERRQ(ierr); 151156ba9f64SToby Isaac ierr = PetscOptionsInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg);CHKERRQ(ierr); 151256ba9f64SToby Isaac if (flg) { 151356ba9f64SToby Isaac ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr); 151456ba9f64SToby Isaac } 1515c7eeac06SToby Isaac ierr = DMForestGetMaximumRefinement(dm,&maxRefinement);CHKERRQ(ierr); 1516c7eeac06SToby Isaac ierr = PetscOptionsInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg);CHKERRQ(ierr); 1517c7eeac06SToby Isaac if (flg) { 1518c7eeac06SToby Isaac ierr = DMForestSetMaximumRefinement(dm,maxRefinement);CHKERRQ(ierr); 1519c7eeac06SToby Isaac } 1520c7eeac06SToby Isaac ierr = DMForestGetAdaptivityStrategy(dm,&adaptStrategy);CHKERRQ(ierr); 1521c7eeac06SToby Isaac ierr = PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,256,&flg);CHKERRQ(ierr); 1522c7eeac06SToby Isaac if (flg) { 1523c7eeac06SToby Isaac ierr = DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);CHKERRQ(ierr); 1524c7eeac06SToby Isaac } 1525c7eeac06SToby Isaac ierr = DMForestGetGradeFactor(dm,&grade);CHKERRQ(ierr); 1526c7eeac06SToby Isaac ierr = PetscOptionsInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg);CHKERRQ(ierr); 1527c7eeac06SToby Isaac if (flg) { 1528c7eeac06SToby Isaac ierr = DMForestSetGradeFactor(dm,grade);CHKERRQ(ierr); 1529c7eeac06SToby Isaac } 1530c7eeac06SToby Isaac ierr = DMForestGetCellWeightFactor(dm,&weightsFactor);CHKERRQ(ierr); 1531c7eeac06SToby Isaac ierr = PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);CHKERRQ(ierr); 1532c7eeac06SToby Isaac if (flg) { 1533c7eeac06SToby Isaac ierr = DMForestSetCellWeightFactor(dm,weightsFactor);CHKERRQ(ierr); 1534c7eeac06SToby Isaac } 1535db4d5e8cSToby Isaac ierr = PetscOptionsTail();CHKERRQ(ierr); 1536db4d5e8cSToby Isaac PetscFunctionReturn(0); 1537db4d5e8cSToby Isaac } 1538db4d5e8cSToby Isaac 1539db4d5e8cSToby Isaac #undef __FUNCT__ 1540d8984e3bSMatthew G. Knepley #define __FUNCT__ "DMCreateSubDM_Forest" 1541d8984e3bSMatthew G. Knepley PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm) 1542d8984e3bSMatthew G. Knepley { 1543d8984e3bSMatthew G. Knepley PetscErrorCode ierr; 1544d8984e3bSMatthew G. Knepley 1545d8984e3bSMatthew G. Knepley PetscFunctionBegin; 1546d8984e3bSMatthew G. Knepley if (subdm) {ierr = DMClone(dm, subdm);CHKERRQ(ierr);} 1547d8984e3bSMatthew G. Knepley ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr); 1548d8984e3bSMatthew G. Knepley PetscFunctionReturn(0); 1549d8984e3bSMatthew G. Knepley } 1550d8984e3bSMatthew G. Knepley 1551d8984e3bSMatthew G. Knepley #undef __FUNCT__ 15525421bac9SToby Isaac #define __FUNCT__ "DMRefine_Forest" 15535421bac9SToby Isaac PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined) 15545421bac9SToby Isaac { 15555421bac9SToby Isaac DMLabel refine; 15565421bac9SToby Isaac DM fineDM; 15575421bac9SToby Isaac PetscErrorCode ierr; 15585421bac9SToby Isaac 15595421bac9SToby Isaac PetscFunctionBegin; 15605421bac9SToby Isaac ierr = DMGetFineDM(dm,&fineDM);CHKERRQ(ierr); 15615421bac9SToby Isaac if (fineDM) { 15625421bac9SToby Isaac ierr = PetscObjectReference((PetscObject)fineDM);CHKERRQ(ierr); 15635421bac9SToby Isaac *dmRefined = fineDM; 15645421bac9SToby Isaac PetscFunctionReturn(0); 15655421bac9SToby Isaac } 15665421bac9SToby Isaac ierr = DMForestTemplate(dm,comm,dmRefined);CHKERRQ(ierr); 15675421bac9SToby Isaac ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr); 15685421bac9SToby Isaac if (!refine) { 15695421bac9SToby Isaac ierr = DMCreateLabel(dm,"refine");CHKERRQ(ierr); 15705421bac9SToby Isaac ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr); 15715421bac9SToby Isaac ierr = DMLabelSetDefaultValue(refine,DM_FOREST_REFINE);CHKERRQ(ierr); 15725421bac9SToby Isaac } 15735421bac9SToby Isaac ierr = DMForestSetAdaptivityLabel(*dmRefined,"refine");CHKERRQ(ierr); 15745421bac9SToby Isaac PetscFunctionReturn(0); 15755421bac9SToby Isaac } 15765421bac9SToby Isaac 15775421bac9SToby Isaac #undef __FUNCT__ 15785421bac9SToby Isaac #define __FUNCT__ "DMCoarsen_Forest" 15795421bac9SToby Isaac PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened) 15805421bac9SToby Isaac { 15815421bac9SToby Isaac DMLabel coarsen; 15825421bac9SToby Isaac DM coarseDM; 15835421bac9SToby Isaac PetscErrorCode ierr; 15845421bac9SToby Isaac 15855421bac9SToby Isaac PetscFunctionBegin; 15864098eed7SToby Isaac { 15874098eed7SToby Isaac PetscMPIInt mpiComparison; 15884098eed7SToby Isaac MPI_Comm dmcomm = PetscObjectComm((PetscObject)dm); 15894098eed7SToby Isaac 15904098eed7SToby Isaac ierr = MPI_Comm_compare(comm, dmcomm, &mpiComparison);CHKERRQ(ierr); 15914098eed7SToby Isaac if (mpiComparison != MPI_IDENT && mpiComparison != MPI_CONGRUENT) { 15924098eed7SToby Isaac SETERRQ(dmcomm,PETSC_ERR_SUP,"No support for different communicators yet"); 15934098eed7SToby Isaac } 15944098eed7SToby Isaac } 15955421bac9SToby Isaac ierr = DMGetCoarseDM(dm,&coarseDM);CHKERRQ(ierr); 15965421bac9SToby Isaac if (coarseDM) { 15975421bac9SToby Isaac ierr = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr); 15985421bac9SToby Isaac *dmCoarsened = coarseDM; 15995421bac9SToby Isaac PetscFunctionReturn(0); 16005421bac9SToby Isaac } 16015421bac9SToby Isaac ierr = DMForestTemplate(dm,comm,dmCoarsened);CHKERRQ(ierr); 16024098eed7SToby Isaac ierr = DMForestSetAdaptivityPurpose(coarseDM,DM_FOREST_COARSEN);CHKERRQ(ierr); 16035421bac9SToby Isaac ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr); 16045421bac9SToby Isaac if (!coarsen) { 16055421bac9SToby Isaac ierr = DMCreateLabel(dm,"coarsen");CHKERRQ(ierr); 16065421bac9SToby Isaac ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr); 16075421bac9SToby Isaac ierr = DMLabelSetDefaultValue(coarsen,DM_FOREST_COARSEN);CHKERRQ(ierr); 16085421bac9SToby Isaac } 16095421bac9SToby Isaac ierr = DMForestSetAdaptivityLabel(*dmCoarsened,"coarsen");CHKERRQ(ierr); 16105421bac9SToby Isaac PetscFunctionReturn(0); 16115421bac9SToby Isaac } 16125421bac9SToby Isaac 16135421bac9SToby Isaac #undef __FUNCT__ 1614d222f98bSToby Isaac #define __FUNCT__ "DMInitialize_Forest" 1615d222f98bSToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm) 1616d222f98bSToby Isaac { 1617d222f98bSToby Isaac PetscErrorCode ierr; 1618d222f98bSToby Isaac 1619d222f98bSToby Isaac PetscFunctionBegin; 1620d222f98bSToby Isaac ierr = PetscMemzero(dm->ops,sizeof(*(dm->ops)));CHKERRQ(ierr); 1621d222f98bSToby Isaac 1622d222f98bSToby Isaac dm->ops->clone = DMClone_Forest; 1623d222f98bSToby Isaac dm->ops->setfromoptions = DMSetFromOptions_Forest; 1624d222f98bSToby Isaac dm->ops->destroy = DMDestroy_Forest; 1625d8984e3bSMatthew G. Knepley dm->ops->createsubdm = DMCreateSubDM_Forest; 16265421bac9SToby Isaac dm->ops->refine = DMRefine_Forest; 16275421bac9SToby Isaac dm->ops->coarsen = DMCoarsen_Forest; 1628d222f98bSToby Isaac PetscFunctionReturn(0); 1629d222f98bSToby Isaac } 1630d222f98bSToby Isaac 16319be51f97SToby Isaac /*MC 16329be51f97SToby 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. 16339be51f97SToby Isaac 16349be51f97SToby 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(). 16359be51f97SToby Isaac 16369be51f97SToby Isaac Level: advanced 16379be51f97SToby Isaac 16389be51f97SToby Isaac .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel() 16399be51f97SToby Isaac M*/ 16409be51f97SToby Isaac 1641d222f98bSToby Isaac #undef __FUNCT__ 1642db4d5e8cSToby Isaac #define __FUNCT__ "DMCreate_Forest" 1643db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm) 1644db4d5e8cSToby Isaac { 1645db4d5e8cSToby Isaac DM_Forest *forest; 1646db4d5e8cSToby Isaac PetscErrorCode ierr; 1647db4d5e8cSToby Isaac 1648db4d5e8cSToby Isaac PetscFunctionBegin; 1649db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1650db4d5e8cSToby Isaac ierr = PetscNewLog(dm,&forest);CHKERRQ(ierr); 1651db4d5e8cSToby Isaac dm->dim = 0; 1652db4d5e8cSToby Isaac dm->data = forest; 1653db4d5e8cSToby Isaac forest->refct = 1; 1654db4d5e8cSToby Isaac forest->data = NULL; 165558762b62SToby Isaac forest->setfromoptionscalled = PETSC_FALSE; 1656db4d5e8cSToby Isaac forest->topology = NULL; 1657db4d5e8cSToby Isaac forest->base = NULL; 1658db4d5e8cSToby Isaac forest->adjDim = PETSC_DEFAULT; 1659db4d5e8cSToby Isaac forest->overlap = PETSC_DEFAULT; 1660db4d5e8cSToby Isaac forest->minRefinement = PETSC_DEFAULT; 1661db4d5e8cSToby Isaac forest->maxRefinement = PETSC_DEFAULT; 166256ba9f64SToby Isaac forest->initRefinement = PETSC_DEFAULT; 1663c7eeac06SToby Isaac forest->cStart = PETSC_DETERMINE; 1664c7eeac06SToby Isaac forest->cEnd = PETSC_DETERMINE; 1665db4d5e8cSToby Isaac forest->cellSF = 0; 1666ebdf65a2SToby Isaac forest->adaptLabel = NULL; 1667db4d5e8cSToby Isaac forest->gradeFactor = 2; 1668db4d5e8cSToby Isaac forest->cellWeights = NULL; 1669db4d5e8cSToby Isaac forest->cellWeightsCopyMode = PETSC_USE_POINTER; 1670db4d5e8cSToby Isaac forest->weightsFactor = 1.; 1671db4d5e8cSToby Isaac forest->weightCapacity = 1.; 1672a73e2921SToby Isaac ierr = DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);CHKERRQ(ierr); 1673d222f98bSToby Isaac ierr = DMInitialize_Forest(dm);CHKERRQ(ierr); 1674db4d5e8cSToby Isaac PetscFunctionReturn(0); 1675db4d5e8cSToby Isaac } 1676db4d5e8cSToby Isaac 1677