19be51f97SToby Isaac #include <petsc/private/dmforestimpl.h> /*I "petscdmforest.h" I*/ 29be51f97SToby Isaac #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 3a1b0c543SToby Isaac #include <petsc/private/dmlabelimpl.h> /*I "petscdmlabel.h" I*/ 4ef19d27cSToby Isaac #include <petscsf.h> 5db4d5e8cSToby Isaac 627d4645fSToby Isaac PetscBool DMForestPackageInitialized = PETSC_FALSE; 727d4645fSToby Isaac 827d4645fSToby Isaac typedef struct _DMForestTypeLink*DMForestTypeLink; 927d4645fSToby Isaac 1027d4645fSToby Isaac struct _DMForestTypeLink 1127d4645fSToby Isaac { 1227d4645fSToby Isaac char *name; 1327d4645fSToby Isaac DMForestTypeLink next; 1427d4645fSToby Isaac }; 1527d4645fSToby Isaac 1627d4645fSToby Isaac DMForestTypeLink DMForestTypeList; 1727d4645fSToby Isaac 1827d4645fSToby Isaac #undef __FUNCT__ 1927d4645fSToby Isaac #define __FUNCT__ "DMForestPackageFinalize" 2027d4645fSToby Isaac static PetscErrorCode DMForestPackageFinalize(void) 2127d4645fSToby Isaac { 2227d4645fSToby Isaac DMForestTypeLink oldLink, link = DMForestTypeList; 2327d4645fSToby Isaac PetscErrorCode ierr; 2427d4645fSToby Isaac 2527d4645fSToby Isaac PetscFunctionBegin; 2627d4645fSToby Isaac while (link) { 2727d4645fSToby Isaac oldLink = link; 28f416af30SBarry Smith ierr = PetscFree(oldLink->name);CHKERRQ(ierr); 2927d4645fSToby Isaac link = oldLink->next; 3027d4645fSToby Isaac ierr = PetscFree(oldLink);CHKERRQ(ierr); 3127d4645fSToby Isaac } 3227d4645fSToby Isaac PetscFunctionReturn(0); 3327d4645fSToby Isaac } 3427d4645fSToby Isaac 3527d4645fSToby Isaac #undef __FUNCT__ 3627d4645fSToby Isaac #define __FUNCT__ "DMForestPackageInitialize" 3727d4645fSToby Isaac static PetscErrorCode DMForestPackageInitialize(void) 3827d4645fSToby Isaac { 3927d4645fSToby Isaac PetscErrorCode ierr; 4027d4645fSToby Isaac 4127d4645fSToby Isaac PetscFunctionBegin; 4227d4645fSToby Isaac if (DMForestPackageInitialized) PetscFunctionReturn(0); 4327d4645fSToby Isaac DMForestPackageInitialized = PETSC_TRUE; 44f885a11aSToby Isaac 4527d4645fSToby Isaac ierr = DMForestRegisterType(DMFOREST);CHKERRQ(ierr); 4627d4645fSToby Isaac ierr = PetscRegisterFinalize(DMForestPackageFinalize);CHKERRQ(ierr); 4727d4645fSToby Isaac PetscFunctionReturn(0); 4827d4645fSToby Isaac } 4927d4645fSToby Isaac 5027d4645fSToby Isaac #undef __FUNCT__ 5127d4645fSToby Isaac #define __FUNCT__ "DMForestRegisterType" 529be51f97SToby Isaac /*@C 539be51f97SToby Isaac DMForestRegisterType - Registers a DMType as a subtype of DMFOREST (so that DMIsForest() will be correct) 549be51f97SToby Isaac 559be51f97SToby Isaac Not Collective 569be51f97SToby Isaac 579be51f97SToby Isaac Input parameter: 589be51f97SToby Isaac . name - the name of the type 599be51f97SToby Isaac 609be51f97SToby Isaac Level: advanced 619be51f97SToby Isaac 629be51f97SToby Isaac .seealso: DMFOREST, DMIsForest() 639be51f97SToby Isaac @*/ 6427d4645fSToby Isaac PetscErrorCode DMForestRegisterType(DMType name) 6527d4645fSToby Isaac { 6627d4645fSToby Isaac DMForestTypeLink link; 6727d4645fSToby Isaac PetscErrorCode ierr; 6827d4645fSToby Isaac 6927d4645fSToby Isaac PetscFunctionBegin; 7027d4645fSToby Isaac ierr = DMForestPackageInitialize();CHKERRQ(ierr); 7127d4645fSToby Isaac ierr = PetscNew(&link);CHKERRQ(ierr); 7227d4645fSToby Isaac ierr = PetscStrallocpy(name,&link->name);CHKERRQ(ierr); 7327d4645fSToby Isaac link->next = DMForestTypeList; 7427d4645fSToby Isaac DMForestTypeList = link; 7527d4645fSToby Isaac PetscFunctionReturn(0); 7627d4645fSToby Isaac } 7727d4645fSToby Isaac 7827d4645fSToby Isaac #undef __FUNCT__ 7927d4645fSToby Isaac #define __FUNCT__ "DMIsForest" 809be51f97SToby Isaac /*@ 819be51f97SToby Isaac DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes 829be51f97SToby Isaac 839be51f97SToby Isaac Not Collective 849be51f97SToby Isaac 859be51f97SToby Isaac Input parameter: 869be51f97SToby Isaac . dm - the DM object 879be51f97SToby Isaac 889be51f97SToby Isaac Output parameter: 899be51f97SToby Isaac . isForest - whether dm is a subtype of DMFOREST 909be51f97SToby Isaac 919be51f97SToby Isaac Level: intermediate 929be51f97SToby Isaac 939be51f97SToby Isaac .seealso: DMFOREST, DMForestRegisterType() 949be51f97SToby Isaac @*/ 9527d4645fSToby Isaac PetscErrorCode DMIsForest(DM dm, PetscBool *isForest) 9627d4645fSToby Isaac { 9727d4645fSToby Isaac DMForestTypeLink link = DMForestTypeList; 9827d4645fSToby Isaac PetscErrorCode ierr; 9927d4645fSToby Isaac 10027d4645fSToby Isaac PetscFunctionBegin; 10127d4645fSToby Isaac while (link) { 10227d4645fSToby Isaac PetscBool sameType; 10327d4645fSToby Isaac ierr = PetscObjectTypeCompare((PetscObject)dm,link->name,&sameType);CHKERRQ(ierr); 10427d4645fSToby Isaac if (sameType) { 10527d4645fSToby Isaac *isForest = PETSC_TRUE; 10627d4645fSToby Isaac PetscFunctionReturn(0); 10727d4645fSToby Isaac } 10827d4645fSToby Isaac link = link->next; 10927d4645fSToby Isaac } 11027d4645fSToby Isaac *isForest = PETSC_FALSE; 11127d4645fSToby Isaac PetscFunctionReturn(0); 11227d4645fSToby Isaac } 11327d4645fSToby Isaac 114db4d5e8cSToby Isaac #undef __FUNCT__ 115a0452a8eSToby Isaac #define __FUNCT__ "DMForestTemplate" 1169be51f97SToby Isaac /*@ 1179be51f97SToby Isaac DMForestTemplate - Create a new DM that will be adapted from a source DM. The new DM reproduces the configuration 1189be51f97SToby 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 1199be51f97SToby Isaac (by, e.g., refinement or repartitioning). The source DM is also set as the adaptivity source DM of the new DM (see 1209be51f97SToby Isaac DMForestSetAdaptivityForest()). 1219be51f97SToby Isaac 1229be51f97SToby Isaac Collective on dm 1239be51f97SToby Isaac 1249be51f97SToby Isaac Input Parameters: 1259be51f97SToby Isaac + dm - the source DM object 1269be51f97SToby Isaac - comm - the communicator for the new DM (this communicator is currently ignored, but is present so that DMForestTemplate() can be used within DMCoarsen()) 1279be51f97SToby Isaac 1289be51f97SToby Isaac Output Parameter: 1299be51f97SToby Isaac . tdm - the new DM object 1309be51f97SToby Isaac 1319be51f97SToby Isaac Level: intermediate 1329be51f97SToby Isaac 1339be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest() 1349be51f97SToby Isaac @*/ 1359be51f97SToby Isaac PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm) 136a0452a8eSToby Isaac { 137a0452a8eSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 13820e8089bSToby Isaac DMType type; 139a0452a8eSToby Isaac DM base; 140a0452a8eSToby Isaac DMForestTopology topology; 141a0452a8eSToby Isaac PetscInt dim, overlap, ref, factor; 142a0452a8eSToby Isaac DMForestAdaptivityStrategy strat; 143795844e7SToby Isaac PetscDS ds; 144795844e7SToby Isaac void *ctx; 14549fc9a2fSToby Isaac PetscErrorCode (*map)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void*); 1463e58adeeSToby Isaac void *mapCtx; 147a0452a8eSToby Isaac PetscErrorCode ierr; 148a0452a8eSToby Isaac 149a0452a8eSToby Isaac PetscFunctionBegin; 150a0452a8eSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 15120e8089bSToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),tdm);CHKERRQ(ierr); 15220e8089bSToby Isaac ierr = DMGetType(dm,&type);CHKERRQ(ierr); 15320e8089bSToby Isaac ierr = DMSetType(*tdm,type);CHKERRQ(ierr); 154a0452a8eSToby Isaac ierr = DMForestGetBaseDM(dm,&base);CHKERRQ(ierr); 15520e8089bSToby Isaac ierr = DMForestSetBaseDM(*tdm,base);CHKERRQ(ierr); 156a0452a8eSToby Isaac ierr = DMForestGetTopology(dm,&topology);CHKERRQ(ierr); 15720e8089bSToby Isaac ierr = DMForestSetTopology(*tdm,topology);CHKERRQ(ierr); 158a0452a8eSToby Isaac ierr = DMForestGetAdjacencyDimension(dm,&dim);CHKERRQ(ierr); 15920e8089bSToby Isaac ierr = DMForestSetAdjacencyDimension(*tdm,dim);CHKERRQ(ierr); 160a0452a8eSToby Isaac ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr); 16120e8089bSToby Isaac ierr = DMForestSetPartitionOverlap(*tdm,overlap);CHKERRQ(ierr); 162a0452a8eSToby Isaac ierr = DMForestGetMinimumRefinement(dm,&ref);CHKERRQ(ierr); 16320e8089bSToby Isaac ierr = DMForestSetMinimumRefinement(*tdm,ref);CHKERRQ(ierr); 164a0452a8eSToby Isaac ierr = DMForestGetMaximumRefinement(dm,&ref);CHKERRQ(ierr); 16520e8089bSToby Isaac ierr = DMForestSetMaximumRefinement(*tdm,ref);CHKERRQ(ierr); 166a0452a8eSToby Isaac ierr = DMForestGetAdaptivityStrategy(dm,&strat);CHKERRQ(ierr); 16720e8089bSToby Isaac ierr = DMForestSetAdaptivityStrategy(*tdm,strat);CHKERRQ(ierr); 168a0452a8eSToby Isaac ierr = DMForestGetGradeFactor(dm,&factor);CHKERRQ(ierr); 16920e8089bSToby Isaac ierr = DMForestSetGradeFactor(*tdm,factor);CHKERRQ(ierr); 1703e58adeeSToby Isaac ierr = DMForestGetBaseCoordinateMapping(dm,&map,&mapCtx);CHKERRQ(ierr); 1715e8c540aSToby Isaac ierr = DMForestSetBaseCoordinateMapping(*tdm,map,mapCtx);CHKERRQ(ierr); 172a0452a8eSToby Isaac if (forest->ftemplate) { 17320e8089bSToby Isaac ierr = (forest->ftemplate)(dm, *tdm);CHKERRQ(ierr); 174a0452a8eSToby Isaac } 17520e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(*tdm,dm);CHKERRQ(ierr); 176795844e7SToby Isaac ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 177795844e7SToby Isaac ierr = DMSetDS(*tdm,ds);CHKERRQ(ierr); 178795844e7SToby Isaac ierr = DMGetApplicationContext(dm,&ctx);CHKERRQ(ierr); 179795844e7SToby Isaac ierr = DMSetApplicationContext(*tdm,&ctx);CHKERRQ(ierr); 180795844e7SToby Isaac if (dm->maxCell) { 181795844e7SToby Isaac const PetscReal *maxCell, *L; 182795844e7SToby Isaac const DMBoundaryType *bd; 183795844e7SToby Isaac 184795844e7SToby Isaac ierr = DMGetPeriodicity(dm,&maxCell,&L,&bd);CHKERRQ(ierr); 185795844e7SToby Isaac ierr = DMSetPeriodicity(*tdm,maxCell,L,bd);CHKERRQ(ierr); 186795844e7SToby Isaac } 187bff67a9bSToby Isaac ierr = DMCopyBoundary(dm,*tdm);CHKERRQ(ierr); 188a0452a8eSToby Isaac PetscFunctionReturn(0); 189a0452a8eSToby Isaac } 190a0452a8eSToby Isaac 19101d9d024SToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm); 19201d9d024SToby Isaac 193a0452a8eSToby Isaac #undef __FUNCT__ 194db4d5e8cSToby Isaac #define __FUNCT__ "DMClone_Forest" 195db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm) 196db4d5e8cSToby Isaac { 197db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 198db4d5e8cSToby Isaac const char *type; 199db4d5e8cSToby Isaac PetscErrorCode ierr; 200db4d5e8cSToby Isaac 201db4d5e8cSToby Isaac PetscFunctionBegin; 202db4d5e8cSToby Isaac forest->refct++; 203db4d5e8cSToby Isaac (*newdm)->data = forest; 204db4d5e8cSToby Isaac ierr = PetscObjectGetType((PetscObject) dm, &type);CHKERRQ(ierr); 205db4d5e8cSToby Isaac ierr = PetscObjectChangeTypeName((PetscObject) *newdm, type);CHKERRQ(ierr); 20601d9d024SToby Isaac ierr = DMInitialize_Forest(*newdm);CHKERRQ(ierr); 207db4d5e8cSToby Isaac PetscFunctionReturn(0); 208db4d5e8cSToby Isaac } 209db4d5e8cSToby Isaac 210db4d5e8cSToby Isaac #undef __FUNCT__ 211db4d5e8cSToby Isaac #define __FUNCT__ "DMDestroy_Forest" 212d222f98bSToby Isaac static PetscErrorCode DMDestroy_Forest(DM dm) 213db4d5e8cSToby Isaac { 214db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 215db4d5e8cSToby Isaac PetscErrorCode ierr; 216db4d5e8cSToby Isaac 217db4d5e8cSToby Isaac PetscFunctionBegin; 218db4d5e8cSToby Isaac if (--forest->refct > 0) PetscFunctionReturn(0); 219d222f98bSToby Isaac if (forest->destroy) {ierr = forest->destroy(dm);CHKERRQ(ierr);} 220db4d5e8cSToby Isaac ierr = PetscSFDestroy(&forest->cellSF);CHKERRQ(ierr); 2210f17b9e3SToby Isaac ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr); 2220f17b9e3SToby Isaac ierr = PetscSFDestroy(&forest->coarseToPreFine);CHKERRQ(ierr); 223a1b0c543SToby Isaac ierr = DMLabelDestroy(&forest->adaptLabel);CHKERRQ(ierr); 2249a81d013SToby Isaac ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr); 22556ba9f64SToby Isaac ierr = DMDestroy(&forest->base);CHKERRQ(ierr); 226ba936b91SToby Isaac ierr = DMDestroy(&forest->adapt);CHKERRQ(ierr); 22730f902e7SToby Isaac ierr = PetscFree(forest->topology);CHKERRQ(ierr); 22830f902e7SToby Isaac ierr = PetscFree(forest);CHKERRQ(ierr); 229a1b0c543SToby Isaac ierr = PetscObjectComposeFunction((PetscObject)dm,"DMAdaptLabel_C",NULL);CHKERRQ(ierr); 230db4d5e8cSToby Isaac PetscFunctionReturn(0); 231db4d5e8cSToby Isaac } 232db4d5e8cSToby Isaac 233db4d5e8cSToby Isaac #undef __FUNCT__ 234dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetTopology" 2359be51f97SToby Isaac /*@C 2369be51f97SToby Isaac DMForestSetTopology - Set the topology of a DMForest during the pre-setup phase. The topology is a string (e.g. 2379be51f97SToby Isaac "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest durint 2389be51f97SToby Isaac DMSetUp(). 2399be51f97SToby Isaac 2409be51f97SToby Isaac Logically collective on dm 2419be51f97SToby Isaac 2429be51f97SToby Isaac Input parameters: 2439be51f97SToby Isaac + dm - the forest 2449be51f97SToby Isaac - topology - the topology of the forest 2459be51f97SToby Isaac 2469be51f97SToby Isaac Level: intermediate 2479be51f97SToby Isaac 2489be51f97SToby Isaac .seealso(): DMForestGetTopology(), DMForestSetBaseDM() 2499be51f97SToby Isaac @*/ 250dd8e54a2SToby Isaac PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology) 251db4d5e8cSToby Isaac { 252db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 253db4d5e8cSToby Isaac PetscErrorCode ierr; 254db4d5e8cSToby Isaac 255db4d5e8cSToby Isaac PetscFunctionBegin; 256db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 257ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the topology after setup"); 258dd8e54a2SToby Isaac ierr = PetscFree(forest->topology);CHKERRQ(ierr); 259dd8e54a2SToby Isaac ierr = PetscStrallocpy((const char*)topology,(char**) &forest->topology);CHKERRQ(ierr); 260db4d5e8cSToby Isaac PetscFunctionReturn(0); 261db4d5e8cSToby Isaac } 262db4d5e8cSToby Isaac 263db4d5e8cSToby Isaac #undef __FUNCT__ 264dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetTopology" 2659be51f97SToby Isaac /*@C 2669be51f97SToby Isaac DMForestGetTopology - Get a string describing the topology of a DMForest. 2679be51f97SToby Isaac 2689be51f97SToby Isaac Not collective 2699be51f97SToby Isaac 2709be51f97SToby Isaac Input parameter: 2719be51f97SToby Isaac . dm - the forest 2729be51f97SToby Isaac 2739be51f97SToby Isaac Output parameter: 2749be51f97SToby Isaac . topology - the topology of the forest (e.g., 'cube', 'shell') 2759be51f97SToby Isaac 2769be51f97SToby Isaac Level: intermediate 2779be51f97SToby Isaac 2789be51f97SToby Isaac .seealso: DMForestSetTopology() 2799be51f97SToby Isaac @*/ 280dd8e54a2SToby Isaac PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology) 281dd8e54a2SToby Isaac { 282dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 283dd8e54a2SToby Isaac 284dd8e54a2SToby Isaac PetscFunctionBegin; 285dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 286dd8e54a2SToby Isaac PetscValidPointer(topology,2); 287dd8e54a2SToby Isaac *topology = forest->topology; 288dd8e54a2SToby Isaac PetscFunctionReturn(0); 289dd8e54a2SToby Isaac } 290dd8e54a2SToby Isaac 291dd8e54a2SToby Isaac #undef __FUNCT__ 292dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetBaseDM" 2939be51f97SToby Isaac /*@ 2949be51f97SToby Isaac DMForestSetBaseDM - During the pre-setup phase, set the DM that defines the base mesh of a DMForest forest. The 2959be51f97SToby Isaac forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its 2969be51f97SToby Isaac base. In general, two forest must share a bse to be comparable, to do things like construct interpolators. 2979be51f97SToby Isaac 2989be51f97SToby Isaac Logically collective on dm 2999be51f97SToby Isaac 3009be51f97SToby Isaac Input Parameters: 3019be51f97SToby Isaac + dm - the forest 3029be51f97SToby Isaac - base - the base DM of the forest 3039be51f97SToby Isaac 3049be51f97SToby Isaac Level: intermediate 3059be51f97SToby Isaac 3069be51f97SToby Isaac .seealso(): DMForestGetBaseDM() 3079be51f97SToby Isaac @*/ 308dd8e54a2SToby Isaac PetscErrorCode DMForestSetBaseDM(DM dm, DM base) 309dd8e54a2SToby Isaac { 310dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 311dd8e54a2SToby Isaac PetscInt dim, dimEmbed; 312dd8e54a2SToby Isaac PetscErrorCode ierr; 313dd8e54a2SToby Isaac 314dd8e54a2SToby Isaac PetscFunctionBegin; 315dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 316ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the base after setup"); 317dd8e54a2SToby Isaac ierr = PetscObjectReference((PetscObject)base);CHKERRQ(ierr); 318dd8e54a2SToby Isaac ierr = DMDestroy(&forest->base);CHKERRQ(ierr); 319dd8e54a2SToby Isaac forest->base = base; 320a0452a8eSToby Isaac if (base) { 321a0452a8eSToby Isaac PetscValidHeaderSpecific(base, DM_CLASSID, 2); 322dd8e54a2SToby Isaac ierr = DMGetDimension(base,&dim);CHKERRQ(ierr); 323dd8e54a2SToby Isaac ierr = DMSetDimension(dm,dim);CHKERRQ(ierr); 324dd8e54a2SToby Isaac ierr = DMGetCoordinateDim(base,&dimEmbed);CHKERRQ(ierr); 325dd8e54a2SToby Isaac ierr = DMSetCoordinateDim(dm,dimEmbed);CHKERRQ(ierr); 326a0452a8eSToby Isaac } 327dd8e54a2SToby Isaac PetscFunctionReturn(0); 328dd8e54a2SToby Isaac } 329dd8e54a2SToby Isaac 330dd8e54a2SToby Isaac #undef __FUNCT__ 331dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetBaseDM" 3329be51f97SToby Isaac /*@ 3339be51f97SToby Isaac DMForestGetBaseDM - Get the base DM of a DMForest forest. The forest will be hierarchically refined from the base, 3349be51f97SToby Isaac and all refinements/coarsenings of the forest will share its base. In general, two forest must share a bse to be 3359be51f97SToby Isaac comparable, to do things like construct interpolators. 3369be51f97SToby Isaac 3379be51f97SToby Isaac Not collective 3389be51f97SToby Isaac 3399be51f97SToby Isaac Input Parameter: 3409be51f97SToby Isaac . dm - the forest 3419be51f97SToby Isaac 3429be51f97SToby Isaac Output Parameter: 3439be51f97SToby Isaac . base - the base DM of the forest 3449be51f97SToby Isaac 3459be51f97SToby Isaac Level: intermediate 3469be51f97SToby Isaac 3479be51f97SToby Isaac .seealso(); DMForestSetBaseDM() 3489be51f97SToby Isaac @*/ 349dd8e54a2SToby Isaac PetscErrorCode DMForestGetBaseDM(DM dm, DM *base) 350dd8e54a2SToby Isaac { 351dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 352dd8e54a2SToby Isaac 353dd8e54a2SToby Isaac PetscFunctionBegin; 354dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 355dd8e54a2SToby Isaac PetscValidPointer(base, 2); 356dd8e54a2SToby Isaac *base = forest->base; 357dd8e54a2SToby Isaac PetscFunctionReturn(0); 358dd8e54a2SToby Isaac } 359dd8e54a2SToby Isaac 360dd8e54a2SToby Isaac #undef __FUNCT__ 361cf38a08cSToby Isaac #define __FUNCT__ "DMForestSetBaseCoordinateMapping" 36299478f86SToby Isaac PetscErrorCode DMForestSetBaseCoordinateMapping(DM dm, PetscErrorCode (*func)(DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx) 363cf38a08cSToby Isaac { 364cf38a08cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 365cf38a08cSToby Isaac 366cf38a08cSToby Isaac PetscFunctionBegin; 367cf38a08cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 368cf38a08cSToby Isaac forest->mapcoordinates = func; 369cf38a08cSToby Isaac forest->mapcoordinatesctx = ctx; 370cf38a08cSToby Isaac PetscFunctionReturn(0); 371cf38a08cSToby Isaac } 372cf38a08cSToby Isaac 373cf38a08cSToby Isaac #undef __FUNCT__ 374cf38a08cSToby Isaac #define __FUNCT__ "DMForestGetBaseCoordinateMapping" 37599478f86SToby Isaac PetscErrorCode DMForestGetBaseCoordinateMapping(DM dm, PetscErrorCode (**func) (DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx) 376cf38a08cSToby Isaac { 377cf38a08cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 378cf38a08cSToby Isaac 379cf38a08cSToby Isaac PetscFunctionBegin; 380cf38a08cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 381cf38a08cSToby Isaac if (func) *func = forest->mapcoordinates; 382cf38a08cSToby Isaac if (ctx) *((void**) ctx) = forest->mapcoordinatesctx; 383cf38a08cSToby Isaac PetscFunctionReturn(0); 384cf38a08cSToby Isaac } 385cf38a08cSToby Isaac 386cf38a08cSToby Isaac #undef __FUNCT__ 387ba936b91SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityForest" 3889be51f97SToby Isaac /*@ 3899be51f97SToby Isaac DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be 3909be51f97SToby Isaac adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) im DMSetUp(). Usually not needed 3919be51f97SToby Isaac by users directly: DMForestTemplate() constructs a new forest to be adapted from an old forest and calls this 3929be51f97SToby Isaac routine. 3939be51f97SToby Isaac 394dffe73a3SToby Isaac Note that this can be called after setup with adapt = NULL, which will clear all internal data related to the 395dffe73a3SToby Isaac adaptivity forest from dm. This way, repeatedly adapting does not leave stale DM objects in memory. 396dffe73a3SToby Isaac 3979be51f97SToby Isaac Logically collective on dm 3989be51f97SToby Isaac 3999be51f97SToby Isaac Input Parameter: 4009be51f97SToby Isaac + dm - the new forest, which will be constructed from adapt 4019be51f97SToby Isaac - adapt - the old forest 4029be51f97SToby Isaac 4039be51f97SToby Isaac Level: intermediate 4049be51f97SToby Isaac 4059be51f97SToby Isaac .seealso: DMForestGetAdaptivityForest(), DMForestSetAdaptivityPurpose() 4069be51f97SToby Isaac @*/ 407ba936b91SToby Isaac PetscErrorCode DMForestSetAdaptivityForest(DM dm,DM adapt) 408dd8e54a2SToby Isaac { 409dffe73a3SToby Isaac DM_Forest *forest, *adaptForest, *oldAdaptForest; 410dffe73a3SToby Isaac DM oldAdapt; 411dd8e54a2SToby Isaac PetscErrorCode ierr; 412dd8e54a2SToby Isaac 413dd8e54a2SToby Isaac PetscFunctionBegin; 414dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 415ba936b91SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 416ba936b91SToby Isaac forest = (DM_Forest*) dm->data; 417dffe73a3SToby Isaac ierr = DMForestGetAdaptivityForest(dm,&oldAdapt);CHKERRQ(ierr); 418dffe73a3SToby Isaac if (adapt != NULL && dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup"); 419193eb951SToby Isaac adaptForest = (DM_Forest*) (adapt ? adapt->data : NULL); 420193eb951SToby Isaac oldAdaptForest = (DM_Forest*) (oldAdapt ? oldAdapt->data : NULL); 421dffe73a3SToby Isaac if (adaptForest != oldAdaptForest) { 422dffe73a3SToby Isaac ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr); 423dffe73a3SToby Isaac ierr = PetscSFDestroy(&forest->coarseToPreFine);CHKERRQ(ierr); 424dffe73a3SToby Isaac if (forest->clearadaptivityforest) {ierr = (forest->clearadaptivityforest)(dm);CHKERRQ(ierr);} 425dffe73a3SToby Isaac } 42626d9498aSToby Isaac switch (forest->adaptPurpose) { 427*cd3c525cSToby Isaac case DM_ADAPT_DETERMINE: 428ba936b91SToby Isaac ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr); 429ba936b91SToby Isaac ierr = DMDestroy(&(forest->adapt));CHKERRQ(ierr); 430ba936b91SToby Isaac forest->adapt = adapt; 43126d9498aSToby Isaac break; 432a1b0c543SToby Isaac case DM_ADAPT_REFINE: 43326d9498aSToby Isaac ierr = DMSetCoarseDM(dm,adapt);CHKERRQ(ierr); 43426d9498aSToby Isaac break; 435a1b0c543SToby Isaac case DM_ADAPT_COARSEN: 43626d9498aSToby Isaac ierr = DMSetFineDM(dm,adapt);CHKERRQ(ierr); 43726d9498aSToby Isaac break; 43826d9498aSToby Isaac default: 43926d9498aSToby Isaac SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose"); 44026d9498aSToby Isaac } 441dd8e54a2SToby Isaac PetscFunctionReturn(0); 442dd8e54a2SToby Isaac } 443dd8e54a2SToby Isaac 444dd8e54a2SToby Isaac #undef __FUNCT__ 445ba936b91SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityForest" 4469be51f97SToby Isaac /*@ 4479be51f97SToby Isaac DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted. 4489be51f97SToby Isaac 4499be51f97SToby Isaac Not collective 4509be51f97SToby Isaac 4519be51f97SToby Isaac Input Parameter: 4529be51f97SToby Isaac . dm - the forest 4539be51f97SToby Isaac 4549be51f97SToby Isaac Output Parameter: 4559be51f97SToby Isaac . adapt - the forest from which dm is/was adapted 4569be51f97SToby Isaac 4579be51f97SToby Isaac Level: intermediate 4589be51f97SToby Isaac 4599be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose() 4609be51f97SToby Isaac @*/ 461ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt) 462dd8e54a2SToby Isaac { 463ba936b91SToby Isaac DM_Forest *forest; 46426d9498aSToby Isaac PetscErrorCode ierr; 465dd8e54a2SToby Isaac 466dd8e54a2SToby Isaac PetscFunctionBegin; 467dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 468ba936b91SToby Isaac forest = (DM_Forest*) dm->data; 46926d9498aSToby Isaac switch (forest->adaptPurpose) { 470*cd3c525cSToby Isaac case DM_ADAPT_DETERMINE: 471ba936b91SToby Isaac *adapt = forest->adapt; 47226d9498aSToby Isaac break; 473a1b0c543SToby Isaac case DM_ADAPT_REFINE: 47426d9498aSToby Isaac ierr = DMGetCoarseDM(dm,adapt);CHKERRQ(ierr); 47526d9498aSToby Isaac break; 476a1b0c543SToby Isaac case DM_ADAPT_COARSEN: 47726d9498aSToby Isaac ierr = DMGetFineDM(dm,adapt);CHKERRQ(ierr); 47826d9498aSToby Isaac break; 47926d9498aSToby Isaac default: 48026d9498aSToby Isaac SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose"); 48126d9498aSToby Isaac } 48226d9498aSToby Isaac PetscFunctionReturn(0); 48326d9498aSToby Isaac } 48426d9498aSToby Isaac 48526d9498aSToby Isaac #undef __FUNCT__ 48626d9498aSToby Isaac #define __FUNCT__ "DMForestSetAdaptivityPurpose" 4879be51f97SToby Isaac /*@ 4889be51f97SToby Isaac DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its 489a1b0c543SToby Isaac source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening 490*cd3c525cSToby Isaac (DM_ADAPT_COARSEN), or undefined (DM_ADAPT_DETERMINE). This only matters for the purposes of reference counting: 4919be51f97SToby Isaac during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse 4929be51f97SToby Isaac relationship (see DMSetFineDM()/DMSetCoarseDM()). If the purpose is not refinement or coarsening, and the user does 4939be51f97SToby Isaac not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can 4949be51f97SToby Isaac cause a memory leak. This method is used by subtypes of DMForest when automatically constructing mesh hierarchies. 4959be51f97SToby Isaac 4969be51f97SToby Isaac Logically collective on dm 4979be51f97SToby Isaac 4989be51f97SToby Isaac Input Parameters: 4999be51f97SToby Isaac + dm - the forest 500*cd3c525cSToby Isaac - purpose - the adaptivity purpose (DM_ADAPT_DETERMINE/DM_ADAPT_REFINE/DM_ADAPT_COARSEN) 5019be51f97SToby Isaac 5029be51f97SToby Isaac Level: advanced 5039be51f97SToby Isaac 5049be51f97SToby Isaac .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest() 5059be51f97SToby Isaac @*/ 506a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMAdaptFlag purpose) 50726d9498aSToby Isaac { 50826d9498aSToby Isaac DM_Forest *forest; 50926d9498aSToby Isaac PetscErrorCode ierr; 51026d9498aSToby Isaac 51126d9498aSToby Isaac PetscFunctionBegin; 51226d9498aSToby Isaac forest = (DM_Forest*) dm->data; 5139be51f97SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup"); 51426d9498aSToby Isaac if (purpose != forest->adaptPurpose) { 51526d9498aSToby Isaac DM adapt; 51626d9498aSToby Isaac 51726d9498aSToby Isaac ierr = DMForestGetAdaptivityForest(dm,&adapt);CHKERRQ(ierr); 51826d9498aSToby Isaac ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr); 51926d9498aSToby Isaac ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr); 520f885a11aSToby Isaac 52126d9498aSToby Isaac forest->adaptPurpose = purpose; 522f885a11aSToby Isaac 52326d9498aSToby Isaac ierr = DMForestSetAdaptivityForest(dm,adapt);CHKERRQ(ierr); 52426d9498aSToby Isaac ierr = DMDestroy(&adapt);CHKERRQ(ierr); 52526d9498aSToby Isaac } 526dd8e54a2SToby Isaac PetscFunctionReturn(0); 527dd8e54a2SToby Isaac } 528dd8e54a2SToby Isaac 529dd8e54a2SToby Isaac #undef __FUNCT__ 53056c0450aSToby Isaac #define __FUNCT__ "DMForestGetAdaptivityPurpose" 53156c0450aSToby Isaac /*@ 53256c0450aSToby Isaac DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with 533a1b0c543SToby Isaac DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening (DM_ADAPT_COARSEN), or 534*cd3c525cSToby Isaac undefined (DM_ADAPT_DETERMINE). This only matters for the purposes of reference counting: during DMDestroy(), cyclic 53556c0450aSToby Isaac references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see 53656c0450aSToby Isaac DMSetFineDM()/DMSetCoarseDM()). If the purpose is not refinement or coarsening, and the user does not maintain a 53756c0450aSToby Isaac reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory 53856c0450aSToby Isaac leak. This method is used by subtypes of DMForest when automatically constructing mesh hierarchies. 53956c0450aSToby Isaac 54056c0450aSToby Isaac Not collective 54156c0450aSToby Isaac 54256c0450aSToby Isaac Input Parameter: 54356c0450aSToby Isaac . dm - the forest 54456c0450aSToby Isaac 54556c0450aSToby Isaac Output Parameter: 546*cd3c525cSToby Isaac . purpose - the adaptivity purpose (DM_ADAPT_DETERMINE/DM_ADAPT_REFINE/DM_ADAPT_COARSEN) 54756c0450aSToby Isaac 54856c0450aSToby Isaac Level: advanced 54956c0450aSToby Isaac 55056c0450aSToby Isaac .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest() 55156c0450aSToby Isaac @*/ 552a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMAdaptFlag *purpose) 55356c0450aSToby Isaac { 55456c0450aSToby Isaac DM_Forest *forest; 55556c0450aSToby Isaac 55656c0450aSToby Isaac PetscFunctionBegin; 55756c0450aSToby Isaac forest = (DM_Forest*) dm->data; 55856c0450aSToby Isaac *purpose = forest->adaptPurpose; 55956c0450aSToby Isaac PetscFunctionReturn(0); 56056c0450aSToby Isaac } 56156c0450aSToby Isaac 56256c0450aSToby Isaac #undef __FUNCT__ 563dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetAdjacencyDimension" 5649be51f97SToby Isaac /*@ 5659be51f97SToby Isaac DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine 5669be51f97SToby Isaac cell adjacency (for the purposes of partitioning and overlap). 5679be51f97SToby Isaac 5689be51f97SToby Isaac Logically collective on dm 5699be51f97SToby Isaac 5709be51f97SToby Isaac Input Parameters: 5719be51f97SToby Isaac + dm - the forest 5729be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency) 5739be51f97SToby Isaac 5749be51f97SToby Isaac Level: intermediate 5759be51f97SToby Isaac 5769be51f97SToby Isaac .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap() 5779be51f97SToby Isaac @*/ 578dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim) 579dd8e54a2SToby Isaac { 580dd8e54a2SToby Isaac PetscInt dim; 581dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 582dd8e54a2SToby Isaac PetscErrorCode ierr; 583dd8e54a2SToby Isaac 584dd8e54a2SToby Isaac PetscFunctionBegin; 585dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 586ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup"); 587dd8e54a2SToby Isaac if (adjDim < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %d", adjDim); 588dd8e54a2SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 589dd8e54a2SToby Isaac if (adjDim > dim) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %d: %d", dim, adjDim); 590dd8e54a2SToby Isaac forest->adjDim = adjDim; 591dd8e54a2SToby Isaac PetscFunctionReturn(0); 592dd8e54a2SToby Isaac } 593dd8e54a2SToby Isaac 594dd8e54a2SToby Isaac #undef __FUNCT__ 595dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetAdjacencyCodimension" 5969be51f97SToby Isaac /*@ 5979be51f97SToby Isaac DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that, 5989be51f97SToby Isaac e.g., adjacency based on facets can be specified by codimension 1 in all cases) 5999be51f97SToby Isaac 6009be51f97SToby Isaac Logically collective on dm 6019be51f97SToby Isaac 6029be51f97SToby Isaac Input Parameters: 6039be51f97SToby Isaac + dm - the forest 6049be51f97SToby Isaac - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices 6059be51f97SToby Isaac 6069be51f97SToby Isaac Level: intermediate 6079be51f97SToby Isaac 6089be51f97SToby Isaac .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension() 6099be51f97SToby Isaac @*/ 610dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim) 611dd8e54a2SToby Isaac { 612dd8e54a2SToby Isaac PetscInt dim; 613dd8e54a2SToby Isaac PetscErrorCode ierr; 614dd8e54a2SToby Isaac 615dd8e54a2SToby Isaac PetscFunctionBegin; 616dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 617dd8e54a2SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 618dd8e54a2SToby Isaac ierr = DMForestSetAdjacencyDimension(dm,dim-adjCodim);CHKERRQ(ierr); 619dd8e54a2SToby Isaac PetscFunctionReturn(0); 620dd8e54a2SToby Isaac } 621dd8e54a2SToby Isaac 622dd8e54a2SToby Isaac #undef __FUNCT__ 623dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetAdjacencyDimension" 6249be51f97SToby Isaac /*@ 6259be51f97SToby Isaac DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the 6269be51f97SToby Isaac purposes of partitioning and overlap). 6279be51f97SToby Isaac 6289be51f97SToby Isaac Not collective 6299be51f97SToby Isaac 6309be51f97SToby Isaac Input Parameter: 6319be51f97SToby Isaac . dm - the forest 6329be51f97SToby Isaac 6339be51f97SToby Isaac Output Parameter: 6349be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency) 6359be51f97SToby Isaac 6369be51f97SToby Isaac Level: intermediate 6379be51f97SToby Isaac 6389be51f97SToby Isaac .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap() 6399be51f97SToby Isaac @*/ 640dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim) 641dd8e54a2SToby Isaac { 642dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 643dd8e54a2SToby Isaac 644dd8e54a2SToby Isaac PetscFunctionBegin; 645dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 646dd8e54a2SToby Isaac PetscValidIntPointer(adjDim,2); 647dd8e54a2SToby Isaac *adjDim = forest->adjDim; 648dd8e54a2SToby Isaac PetscFunctionReturn(0); 649dd8e54a2SToby Isaac } 650dd8e54a2SToby Isaac 651dd8e54a2SToby Isaac #undef __FUNCT__ 652dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetAdjacencyCodimension" 6539be51f97SToby Isaac /*@ 6549be51f97SToby Isaac DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that, 6559be51f97SToby Isaac e.g., adjacency based on facets can be specified by codimension 1 in all cases) 6569be51f97SToby Isaac 6579be51f97SToby Isaac Not collective 6589be51f97SToby Isaac 6599be51f97SToby Isaac Input Parameter: 6609be51f97SToby Isaac . dm - the forest 6619be51f97SToby Isaac 6629be51f97SToby Isaac Output Parameter: 6639be51f97SToby Isaac . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices 6649be51f97SToby Isaac 6659be51f97SToby Isaac Level: intermediate 6669be51f97SToby Isaac 6679be51f97SToby Isaac .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension() 6689be51f97SToby Isaac @*/ 669dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim) 670dd8e54a2SToby Isaac { 671dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 672dd8e54a2SToby Isaac PetscInt dim; 673dd8e54a2SToby Isaac PetscErrorCode ierr; 674dd8e54a2SToby Isaac 675dd8e54a2SToby Isaac PetscFunctionBegin; 676dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 677dd8e54a2SToby Isaac PetscValidIntPointer(adjCodim,2); 678dd8e54a2SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 679dd8e54a2SToby Isaac *adjCodim = dim - forest->adjDim; 680dd8e54a2SToby Isaac PetscFunctionReturn(0); 681dd8e54a2SToby Isaac } 682dd8e54a2SToby Isaac 683dd8e54a2SToby Isaac #undef __FUNCT__ 684ef51cf95SToby Isaac #define __FUNCT__ "DMForestSetPartitionOverlap" 6859be51f97SToby Isaac /*@ 6869be51f97SToby Isaac DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel 6879be51f97SToby Isaac partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding 6889be51f97SToby Isaac adjacent cells 6899be51f97SToby Isaac 6909be51f97SToby Isaac Logically collective on dm 6919be51f97SToby Isaac 6929be51f97SToby Isaac Input Parameters: 6939be51f97SToby Isaac + dm - the forest 6949be51f97SToby Isaac - overlap - default 0 6959be51f97SToby Isaac 6969be51f97SToby Isaac Level: intermediate 6979be51f97SToby Isaac 6989be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension() 6999be51f97SToby Isaac @*/ 700dd8e54a2SToby Isaac PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap) 701dd8e54a2SToby Isaac { 702dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 703dd8e54a2SToby Isaac 704dd8e54a2SToby Isaac PetscFunctionBegin; 705dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 706ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup"); 707dd8e54a2SToby Isaac if (overlap < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %d", overlap); 708dd8e54a2SToby Isaac forest->overlap = overlap; 709dd8e54a2SToby Isaac PetscFunctionReturn(0); 710dd8e54a2SToby Isaac } 711dd8e54a2SToby Isaac 712dd8e54a2SToby Isaac #undef __FUNCT__ 713dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetPartitionOverlap" 7149be51f97SToby Isaac /*@ 7159be51f97SToby Isaac DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values 7169be51f97SToby Isaac > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells 7179be51f97SToby Isaac 7189be51f97SToby Isaac Not collective 7199be51f97SToby Isaac 7209be51f97SToby Isaac Input Parameter: 7219be51f97SToby Isaac . dm - the forest 7229be51f97SToby Isaac 7239be51f97SToby Isaac Output Parameter: 7249be51f97SToby Isaac . overlap - default 0 7259be51f97SToby Isaac 7269be51f97SToby Isaac Level: intermediate 7279be51f97SToby Isaac 7289be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension() 7299be51f97SToby Isaac @*/ 730dd8e54a2SToby Isaac PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap) 731dd8e54a2SToby Isaac { 732dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 733dd8e54a2SToby Isaac 734dd8e54a2SToby Isaac PetscFunctionBegin; 735dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 736dd8e54a2SToby Isaac PetscValidIntPointer(overlap,2); 737dd8e54a2SToby Isaac *overlap = forest->overlap; 738dd8e54a2SToby Isaac PetscFunctionReturn(0); 739dd8e54a2SToby Isaac } 740dd8e54a2SToby Isaac 741dd8e54a2SToby Isaac #undef __FUNCT__ 742dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetMinimumRefinement" 7439be51f97SToby Isaac /*@ 7449be51f97SToby Isaac DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base 7459be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. If the forest is being created by coarsening a previous forest 7469be51f97SToby Isaac (see DMForestGetAdaptivityForest()) this limits the amount of coarsening. 7479be51f97SToby Isaac 7489be51f97SToby Isaac Logically collective on dm 7499be51f97SToby Isaac 7509be51f97SToby Isaac Input Parameters: 7519be51f97SToby Isaac + dm - the forest 7529be51f97SToby Isaac - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7539be51f97SToby Isaac 7549be51f97SToby Isaac Level: intermediate 7559be51f97SToby Isaac 7569be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 7579be51f97SToby Isaac @*/ 758dd8e54a2SToby Isaac PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement) 759dd8e54a2SToby Isaac { 760dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 761dd8e54a2SToby Isaac 762dd8e54a2SToby Isaac PetscFunctionBegin; 763dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 764ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup"); 765dd8e54a2SToby Isaac forest->minRefinement = minRefinement; 766dd8e54a2SToby Isaac PetscFunctionReturn(0); 767dd8e54a2SToby Isaac } 768dd8e54a2SToby Isaac 769dd8e54a2SToby Isaac #undef __FUNCT__ 770dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetMinimumRefinement" 7719be51f97SToby Isaac /*@ 7729be51f97SToby Isaac DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see 7739be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. If the forest is being created by coarsening a previous forest (see 7749be51f97SToby Isaac DMForestGetAdaptivityForest()), this limits the amount of coarsening. 7759be51f97SToby Isaac 7769be51f97SToby Isaac Not collective 7779be51f97SToby Isaac 7789be51f97SToby Isaac Input Parameter: 7799be51f97SToby Isaac . dm - the forest 7809be51f97SToby Isaac 7819be51f97SToby Isaac Output Parameter: 7829be51f97SToby Isaac . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7839be51f97SToby Isaac 7849be51f97SToby Isaac Level: intermediate 7859be51f97SToby Isaac 7869be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestGetMaximumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 7879be51f97SToby Isaac @*/ 788dd8e54a2SToby Isaac PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement) 789dd8e54a2SToby Isaac { 790dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 791dd8e54a2SToby Isaac 792dd8e54a2SToby Isaac PetscFunctionBegin; 793dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 794dd8e54a2SToby Isaac PetscValidIntPointer(minRefinement,2); 795dd8e54a2SToby Isaac *minRefinement = forest->minRefinement; 796dd8e54a2SToby Isaac PetscFunctionReturn(0); 797dd8e54a2SToby Isaac } 798dd8e54a2SToby Isaac 799dd8e54a2SToby Isaac #undef __FUNCT__ 80056ba9f64SToby Isaac #define __FUNCT__ "DMForestSetInitialRefinement" 8019be51f97SToby Isaac /*@ 8029be51f97SToby Isaac DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base 8039be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. 8049be51f97SToby Isaac 8059be51f97SToby Isaac Logically collective on dm 8069be51f97SToby Isaac 8079be51f97SToby Isaac Input Parameters: 8089be51f97SToby Isaac + dm - the forest 8099be51f97SToby Isaac - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 8109be51f97SToby Isaac 8119be51f97SToby Isaac Level: intermediate 8129be51f97SToby Isaac 8139be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM() 8149be51f97SToby Isaac @*/ 81556ba9f64SToby Isaac PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement) 81656ba9f64SToby Isaac { 81756ba9f64SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 81856ba9f64SToby Isaac 81956ba9f64SToby Isaac PetscFunctionBegin; 82056ba9f64SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 82156ba9f64SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup"); 82256ba9f64SToby Isaac forest->initRefinement = initRefinement; 82356ba9f64SToby Isaac PetscFunctionReturn(0); 82456ba9f64SToby Isaac } 82556ba9f64SToby Isaac 82656ba9f64SToby Isaac #undef __FUNCT__ 82756ba9f64SToby Isaac #define __FUNCT__ "DMForestGetInitialRefinement" 8289be51f97SToby Isaac /*@ 8299be51f97SToby Isaac DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see 8309be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. 8319be51f97SToby Isaac 8329be51f97SToby Isaac Not collective 8339be51f97SToby Isaac 8349be51f97SToby Isaac Input Parameter: 8359be51f97SToby Isaac . dm - the forest 8369be51f97SToby Isaac 8379be51f97SToby Isaac Output Paramater: 8389be51f97SToby Isaac . initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 8399be51f97SToby Isaac 8409be51f97SToby Isaac Level: intermediate 8419be51f97SToby Isaac 8429be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM() 8439be51f97SToby Isaac @*/ 84456ba9f64SToby Isaac PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement) 84556ba9f64SToby Isaac { 84656ba9f64SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 84756ba9f64SToby Isaac 84856ba9f64SToby Isaac PetscFunctionBegin; 84956ba9f64SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 85056ba9f64SToby Isaac PetscValidIntPointer(initRefinement,2); 85156ba9f64SToby Isaac *initRefinement = forest->initRefinement; 85256ba9f64SToby Isaac PetscFunctionReturn(0); 85356ba9f64SToby Isaac } 85456ba9f64SToby Isaac 85556ba9f64SToby Isaac #undef __FUNCT__ 856c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetMaximumRefinement" 8579be51f97SToby Isaac /*@ 8589be51f97SToby Isaac DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base 8599be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. If the forest is being created by refining a previous forest 8609be51f97SToby Isaac (see DMForestGetAdaptivityForest()), this limits the amount of refinement. 8619be51f97SToby Isaac 8629be51f97SToby Isaac Logically collective on dm 8639be51f97SToby Isaac 8649be51f97SToby Isaac Input Parameters: 8659be51f97SToby Isaac + dm - the forest 8669be51f97SToby Isaac - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 8679be51f97SToby Isaac 8689be51f97SToby Isaac Level: intermediate 8699be51f97SToby Isaac 8709be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityDM() 8719be51f97SToby Isaac @*/ 872c7eeac06SToby Isaac PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement) 873dd8e54a2SToby Isaac { 874dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 875dd8e54a2SToby Isaac 876dd8e54a2SToby Isaac PetscFunctionBegin; 877dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 878ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup"); 879c7eeac06SToby Isaac forest->maxRefinement = maxRefinement; 880dd8e54a2SToby Isaac PetscFunctionReturn(0); 881dd8e54a2SToby Isaac } 882dd8e54a2SToby Isaac 883dd8e54a2SToby Isaac #undef __FUNCT__ 884c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetMaximumRefinement" 8859be51f97SToby Isaac /*@ 8869be51f97SToby Isaac DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see 8879be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. If the forest is being created by refining a previous forest (see 8889be51f97SToby Isaac DMForestGetAdaptivityForest()), this limits the amount of refinement. 8899be51f97SToby Isaac 8909be51f97SToby Isaac Not collective 8919be51f97SToby Isaac 8929be51f97SToby Isaac Input Parameter: 8939be51f97SToby Isaac . dm - the forest 8949be51f97SToby Isaac 8959be51f97SToby Isaac Output Parameter: 8969be51f97SToby Isaac . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 8979be51f97SToby Isaac 8989be51f97SToby Isaac Level: intermediate 8999be51f97SToby Isaac 9009be51f97SToby Isaac .seealso: DMForestSetMaximumRefinement(), DMForestGetMinimumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 9019be51f97SToby Isaac @*/ 902c7eeac06SToby Isaac PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement) 903dd8e54a2SToby Isaac { 904dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 905dd8e54a2SToby Isaac 906dd8e54a2SToby Isaac PetscFunctionBegin; 907dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 908c7eeac06SToby Isaac PetscValidIntPointer(maxRefinement,2); 909c7eeac06SToby Isaac *maxRefinement = forest->maxRefinement; 910dd8e54a2SToby Isaac PetscFunctionReturn(0); 911dd8e54a2SToby Isaac } 912c7eeac06SToby Isaac 913c7eeac06SToby Isaac #undef __FUNCT__ 914c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityStrategy" 9159be51f97SToby Isaac /*@C 9169be51f97SToby Isaac DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes. 9179be51f97SToby Isaac Subtypes of DMForest may define their own strategies. Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree 9189be51f97SToby Isaac for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to 9199be51f97SToby Isaac specify refinement/coarsening. 9209be51f97SToby Isaac 9219be51f97SToby Isaac Logically collective on dm 9229be51f97SToby Isaac 9239be51f97SToby Isaac Input Parameters: 9249be51f97SToby Isaac + dm - the forest 9259be51f97SToby Isaac - adaptStrategy - default DMFORESTADAPTALL 9269be51f97SToby Isaac 9279be51f97SToby Isaac Level: advanced 9289be51f97SToby Isaac 9299be51f97SToby Isaac .seealso: DMForestGetAdaptivityStrategy() 9309be51f97SToby Isaac @*/ 931c7eeac06SToby Isaac PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy) 932c7eeac06SToby Isaac { 933c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 934c7eeac06SToby Isaac PetscErrorCode ierr; 935c7eeac06SToby Isaac 936c7eeac06SToby Isaac PetscFunctionBegin; 937c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 938c7eeac06SToby Isaac ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr); 939a73e2921SToby Isaac ierr = PetscStrallocpy((const char*) adaptStrategy,(char**)&forest->adaptStrategy);CHKERRQ(ierr); 940c7eeac06SToby Isaac PetscFunctionReturn(0); 941c7eeac06SToby Isaac } 942c7eeac06SToby Isaac 943c7eeac06SToby Isaac #undef __FUNCT__ 944c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityStrategy" 9459be51f97SToby Isaac /*@C 9469be51f97SToby Isaac DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes. Subtypes 9479be51f97SToby Isaac of DMForest may define their own strategies. Two default strategies are DMFORESTADAPTALL, which indicates that all 9489be51f97SToby Isaac processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only 9499be51f97SToby Isaac one process needs to specify refinement/coarsening. 9509be51f97SToby Isaac 9519be51f97SToby Isaac Not collective 9529be51f97SToby Isaac 9539be51f97SToby Isaac Input Parameter: 9549be51f97SToby Isaac . dm - the forest 9559be51f97SToby Isaac 9569be51f97SToby Isaac Output Parameter: 9579be51f97SToby Isaac . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL) 9589be51f97SToby Isaac 9599be51f97SToby Isaac Level: advanced 9609be51f97SToby Isaac 9619be51f97SToby Isaac .seealso: DMForestSetAdaptivityStrategy() 9629be51f97SToby Isaac @*/ 963c7eeac06SToby Isaac PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy) 964c7eeac06SToby Isaac { 965c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 966c7eeac06SToby Isaac 967c7eeac06SToby Isaac PetscFunctionBegin; 968c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 969c7eeac06SToby Isaac PetscValidPointer(adaptStrategy,2); 970c7eeac06SToby Isaac *adaptStrategy = forest->adaptStrategy; 971c7eeac06SToby Isaac PetscFunctionReturn(0); 972c7eeac06SToby Isaac } 973c7eeac06SToby Isaac 974c7eeac06SToby Isaac #undef __FUNCT__ 9752a133e43SToby Isaac #define __FUNCT__ "DMForestGetAdaptivitySuccess" 9762a133e43SToby Isaac /*@ 9772a133e43SToby Isaac DMForestGetAdaptivitySuccess - Return whether the requested adaptation (refinement, coarsening, repartitioning, 9782a133e43SToby Isaac etc.) was successful. PETSC_FALSE indicates that the post-adaptation forest is the same as the pre-adpatation 9792a133e43SToby Isaac forest. A requested adaptation may have been unsuccessful if, for example, the requested refinement would have 9802a133e43SToby Isaac exceeded the maximum refinement level. 9812a133e43SToby Isaac 9822a133e43SToby Isaac Collective on dm 9832a133e43SToby Isaac 9842a133e43SToby Isaac Input Parameter: 9852a133e43SToby Isaac 9862a133e43SToby Isaac . dm - the post-adaptation forest 9872a133e43SToby Isaac 9882a133e43SToby Isaac Output Parameter: 9892a133e43SToby Isaac 9902a133e43SToby Isaac . success - PETSC_TRUE if the post-adaptation forest is different from the pre-adaptation forest. 9912a133e43SToby Isaac 9922a133e43SToby Isaac Level: intermediate 9932a133e43SToby Isaac 9942a133e43SToby Isaac .see 9952a133e43SToby Isaac @*/ 9962a133e43SToby Isaac PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success) 9972a133e43SToby Isaac { 9982a133e43SToby Isaac DM_Forest *forest; 9992a133e43SToby Isaac PetscErrorCode ierr; 10002a133e43SToby Isaac 10012a133e43SToby Isaac PetscFunctionBegin; 10022a133e43SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 10032a133e43SToby Isaac if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() has not been called yet."); 10042a133e43SToby Isaac forest = (DM_Forest *) dm->data; 10052a133e43SToby Isaac ierr = (forest->getadaptivitysuccess)(dm,success);CHKERRQ(ierr); 10062a133e43SToby Isaac PetscFunctionReturn(0); 10072a133e43SToby Isaac } 10082a133e43SToby Isaac 10092a133e43SToby Isaac #undef __FUNCT__ 1010bf9b5d84SToby Isaac #define __FUNCT__ "DMForestSetComputeAdaptivitySF" 1011bf9b5d84SToby Isaac /*@ 1012bf9b5d84SToby Isaac DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed 1013bf9b5d84SToby 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(). 1014bf9b5d84SToby Isaac 1015bf9b5d84SToby Isaac Logically collective on dm 1016bf9b5d84SToby Isaac 1017bf9b5d84SToby Isaac Input Parameters: 1018bf9b5d84SToby Isaac + dm - the post-adaptation forest 1019bf9b5d84SToby Isaac - computeSF - default PETSC_TRUE 1020bf9b5d84SToby Isaac 1021bf9b5d84SToby Isaac Level: advanced 1022bf9b5d84SToby Isaac 1023bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF() 1024bf9b5d84SToby Isaac @*/ 1025bf9b5d84SToby Isaac PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF) 1026bf9b5d84SToby Isaac { 1027bf9b5d84SToby Isaac DM_Forest *forest; 1028bf9b5d84SToby Isaac 1029bf9b5d84SToby Isaac PetscFunctionBegin; 1030bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1031bf9b5d84SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called"); 1032bf9b5d84SToby Isaac forest = (DM_Forest*) dm->data; 1033bf9b5d84SToby Isaac forest->computeAdaptSF = computeSF; 1034bf9b5d84SToby Isaac PetscFunctionReturn(0); 1035bf9b5d84SToby Isaac } 1036bf9b5d84SToby Isaac 1037bf9b5d84SToby Isaac #undef __FUNCT__ 103880b27e07SToby Isaac #define __FUNCT__ "DMForestTransferVec" 10390eb7e1eaSToby Isaac PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time) 104080b27e07SToby Isaac { 104180b27e07SToby Isaac DM_Forest *forest; 104280b27e07SToby Isaac PetscErrorCode ierr; 104380b27e07SToby Isaac 104480b27e07SToby Isaac PetscFunctionBegin; 104580b27e07SToby Isaac PetscValidHeaderSpecific(dmIn ,DM_CLASSID ,1); 104680b27e07SToby Isaac PetscValidHeaderSpecific(vecIn ,VEC_CLASSID ,2); 104780b27e07SToby Isaac PetscValidHeaderSpecific(dmOut ,DM_CLASSID ,3); 104880b27e07SToby Isaac PetscValidHeaderSpecific(vecOut ,VEC_CLASSID ,4); 104980b27e07SToby Isaac forest = (DM_Forest *) dmIn->data; 105080b27e07SToby Isaac if (!forest->transfervec) SETERRQ(PetscObjectComm((PetscObject)dmIn),PETSC_ERR_SUP,"DMForestTransferVec() not implemented"); 10510eb7e1eaSToby Isaac ierr = (forest->transfervec)(dmIn,vecIn,dmOut,vecOut,useBCs,time);CHKERRQ(ierr); 105280b27e07SToby Isaac PetscFunctionReturn(0); 105380b27e07SToby Isaac } 105480b27e07SToby Isaac 105580b27e07SToby Isaac #undef __FUNCT__ 1056bf9b5d84SToby Isaac #define __FUNCT__ "DMForestGetComputeAdaptivitySF" 1057bf9b5d84SToby Isaac /*@ 1058bf9b5d84SToby Isaac DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the 1059bf9b5d84SToby Isaac pre-adaptation forest to the post-adaptiation forest. After DMSetUp() is called, these transfer PetscSFs can be 1060bf9b5d84SToby Isaac accessed with DMForestGetAdaptivitySF(). 1061bf9b5d84SToby Isaac 1062bf9b5d84SToby Isaac Not collective 1063bf9b5d84SToby Isaac 1064bf9b5d84SToby Isaac Input Parameter: 1065bf9b5d84SToby Isaac . dm - the post-adaptation forest 1066bf9b5d84SToby Isaac 1067bf9b5d84SToby Isaac Output Parameter: 1068bf9b5d84SToby Isaac . computeSF - default PETSC_TRUE 1069bf9b5d84SToby Isaac 1070bf9b5d84SToby Isaac Level: advanced 1071bf9b5d84SToby Isaac 1072bf9b5d84SToby Isaac .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF() 1073bf9b5d84SToby Isaac @*/ 1074bf9b5d84SToby Isaac PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF) 1075bf9b5d84SToby Isaac { 1076bf9b5d84SToby Isaac DM_Forest *forest; 1077bf9b5d84SToby Isaac 1078bf9b5d84SToby Isaac PetscFunctionBegin; 1079bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1080bf9b5d84SToby Isaac forest = (DM_Forest*) dm->data; 1081bf9b5d84SToby Isaac *computeSF = forest->computeAdaptSF; 1082bf9b5d84SToby Isaac PetscFunctionReturn(0); 1083bf9b5d84SToby Isaac } 1084bf9b5d84SToby Isaac 1085bf9b5d84SToby Isaac #undef __FUNCT__ 1086bf9b5d84SToby Isaac #define __FUNCT__ "DMForestGetAdaptivitySF" 1087bf9b5d84SToby Isaac /*@ 1088bf9b5d84SToby Isaac DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest. 1089bf9b5d84SToby Isaac Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be 1090bf9b5d84SToby Isaac some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa. Therefore there are two 1091bf9b5d84SToby Isaac PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates 1092bf9b5d84SToby Isaac pre-adaptation fine cells to post-adaptation coarse cells. 1093bf9b5d84SToby Isaac 1094bf9b5d84SToby Isaac Not collective 1095bf9b5d84SToby Isaac 1096bf9b5d84SToby Isaac Input Parameter: 1097bf9b5d84SToby Isaac dm - the post-adaptation forest 1098bf9b5d84SToby Isaac 1099bf9b5d84SToby Isaac Output Parameter: 11000f17b9e3SToby Isaac preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post- 11010f17b9e3SToby Isaac coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre- 1102bf9b5d84SToby Isaac 1103bf9b5d84SToby Isaac Level: advanced 1104bf9b5d84SToby Isaac 1105bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF() 1106bf9b5d84SToby Isaac @*/ 11070f17b9e3SToby Isaac PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine) 1108bf9b5d84SToby Isaac { 1109bf9b5d84SToby Isaac DM_Forest *forest; 1110bf9b5d84SToby Isaac PetscErrorCode ierr; 1111bf9b5d84SToby Isaac 1112bf9b5d84SToby Isaac PetscFunctionBegin; 1113bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1114bf9b5d84SToby Isaac ierr = DMSetUp(dm);CHKERRQ(ierr); 1115bf9b5d84SToby Isaac forest = (DM_Forest*) dm->data; 1116f885a11aSToby Isaac if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine; 1117f885a11aSToby Isaac if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine; 1118bf9b5d84SToby Isaac PetscFunctionReturn(0); 1119bf9b5d84SToby Isaac } 1120bf9b5d84SToby Isaac 1121bf9b5d84SToby Isaac #undef __FUNCT__ 1122c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetGradeFactor" 11239be51f97SToby Isaac /*@ 11249be51f97SToby Isaac DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to 11259be51f97SToby Isaac indicate that the diameter of neighboring cells should differ by at most a factor of 2. Subtypes of DMForest may 11269be51f97SToby Isaac only support one particular choice of grading factor. 11279be51f97SToby Isaac 11289be51f97SToby Isaac Logically collective on dm 11299be51f97SToby Isaac 11309be51f97SToby Isaac Input Parameters: 11319be51f97SToby Isaac + dm - the forest 11329be51f97SToby Isaac - grade - the grading factor 11339be51f97SToby Isaac 11349be51f97SToby Isaac Level: advanced 11359be51f97SToby Isaac 11369be51f97SToby Isaac .seealso: DMForestGetGradeFactor() 11379be51f97SToby Isaac @*/ 1138c7eeac06SToby Isaac PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade) 1139c7eeac06SToby Isaac { 1140c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1141c7eeac06SToby Isaac 1142c7eeac06SToby Isaac PetscFunctionBegin; 1143c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1144ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup"); 1145c7eeac06SToby Isaac forest->gradeFactor = grade; 1146c7eeac06SToby Isaac PetscFunctionReturn(0); 1147c7eeac06SToby Isaac } 1148c7eeac06SToby Isaac 1149c7eeac06SToby Isaac #undef __FUNCT__ 1150c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetGradeFactor" 11519be51f97SToby Isaac /*@ 11529be51f97SToby Isaac DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of 11539be51f97SToby Isaac neighboring cells should differ by at most a factor of 2. Subtypes of DMForest may only support one particular 11549be51f97SToby Isaac choice of grading factor. 11559be51f97SToby Isaac 11569be51f97SToby Isaac Not collective 11579be51f97SToby Isaac 11589be51f97SToby Isaac Input Parameter: 11599be51f97SToby Isaac . dm - the forest 11609be51f97SToby Isaac 11619be51f97SToby Isaac Output Parameter: 11629be51f97SToby Isaac . grade - the grading factor 11639be51f97SToby Isaac 11649be51f97SToby Isaac Level: advanced 11659be51f97SToby Isaac 11669be51f97SToby Isaac .seealso: DMForestSetGradeFactor() 11679be51f97SToby Isaac @*/ 1168c7eeac06SToby Isaac PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade) 1169c7eeac06SToby Isaac { 1170c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1171c7eeac06SToby Isaac 1172c7eeac06SToby Isaac PetscFunctionBegin; 1173c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1174c7eeac06SToby Isaac PetscValidIntPointer(grade,2); 1175c7eeac06SToby Isaac *grade = forest->gradeFactor; 1176c7eeac06SToby Isaac PetscFunctionReturn(0); 1177c7eeac06SToby Isaac } 1178c7eeac06SToby Isaac 1179c7eeac06SToby Isaac #undef __FUNCT__ 1180ef51cf95SToby Isaac #define __FUNCT__ "DMForestSetCellWeightFactor" 11819be51f97SToby Isaac /*@ 11829be51f97SToby Isaac DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes 11839be51f97SToby Isaac the cell weight (see DMForestSetCellWeights()) when calculating partitions. The final weight of a cell will be 11849be51f97SToby Isaac (cellWeight) * (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on 11859be51f97SToby Isaac its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the 11869be51f97SToby Isaac computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement. 11879be51f97SToby Isaac 11889be51f97SToby Isaac Logically collective on dm 11899be51f97SToby Isaac 11909be51f97SToby Isaac Input Parameters: 11919be51f97SToby Isaac + dm - the forest 11929be51f97SToby Isaac - weightsFactors - default 1. 11939be51f97SToby Isaac 11949be51f97SToby Isaac Level: advanced 11959be51f97SToby Isaac 11969be51f97SToby Isaac .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights() 11979be51f97SToby Isaac @*/ 1198ef51cf95SToby Isaac PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor) 1199c7eeac06SToby Isaac { 1200c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1201c7eeac06SToby Isaac 1202c7eeac06SToby Isaac PetscFunctionBegin; 1203c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1204ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup"); 1205c7eeac06SToby Isaac forest->weightsFactor = weightsFactor; 1206c7eeac06SToby Isaac PetscFunctionReturn(0); 1207c7eeac06SToby Isaac } 1208c7eeac06SToby Isaac 1209c7eeac06SToby Isaac #undef __FUNCT__ 1210ef51cf95SToby Isaac #define __FUNCT__ "DMForestGetCellWeightFactor" 12119be51f97SToby Isaac /*@ 12129be51f97SToby Isaac DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see 12139be51f97SToby Isaac DMForestSetCellWeights()) when calculating partitions. The final weight of a cell will be (cellWeight) * 12149be51f97SToby Isaac (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on its level; a 12159be51f97SToby Isaac factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation 12169be51f97SToby Isaac associated with a cell is multiplied by a factor of 2 for each additional level of refinement. 12179be51f97SToby Isaac 12189be51f97SToby Isaac Not collective 12199be51f97SToby Isaac 12209be51f97SToby Isaac Input Parameter: 12219be51f97SToby Isaac . dm - the forest 12229be51f97SToby Isaac 12239be51f97SToby Isaac Output Parameter: 12249be51f97SToby Isaac . weightsFactors - default 1. 12259be51f97SToby Isaac 12269be51f97SToby Isaac Level: advanced 12279be51f97SToby Isaac 12289be51f97SToby Isaac .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights() 12299be51f97SToby Isaac @*/ 1230ef51cf95SToby Isaac PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor) 1231c7eeac06SToby Isaac { 1232c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1233c7eeac06SToby Isaac 1234c7eeac06SToby Isaac PetscFunctionBegin; 1235c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1236c7eeac06SToby Isaac PetscValidRealPointer(weightsFactor,2); 1237c7eeac06SToby Isaac *weightsFactor = forest->weightsFactor; 1238c7eeac06SToby Isaac PetscFunctionReturn(0); 1239c7eeac06SToby Isaac } 1240c7eeac06SToby Isaac 1241c7eeac06SToby Isaac #undef __FUNCT__ 1242c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellChart" 12439be51f97SToby Isaac /*@ 12449be51f97SToby Isaac DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process 12459be51f97SToby Isaac 12469be51f97SToby Isaac Not collective 12479be51f97SToby Isaac 12489be51f97SToby Isaac Input Parameter: 12499be51f97SToby Isaac . dm - the forest 12509be51f97SToby Isaac 12519be51f97SToby Isaac Output Parameters: 12529be51f97SToby Isaac + cStart - the first cell on this process 12539be51f97SToby Isaac - cEnd - one after the final cell on this process 12549be51f97SToby Isaac 12551a244344SSatish Balay Level: intermediate 12569be51f97SToby Isaac 12579be51f97SToby Isaac .seealso: DMForestGetCellSF() 12589be51f97SToby Isaac @*/ 1259c7eeac06SToby Isaac PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd) 1260c7eeac06SToby Isaac { 1261c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1262c7eeac06SToby Isaac PetscErrorCode ierr; 1263c7eeac06SToby Isaac 1264c7eeac06SToby Isaac PetscFunctionBegin; 1265c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1266c7eeac06SToby Isaac PetscValidIntPointer(cStart,2); 1267c7eeac06SToby Isaac PetscValidIntPointer(cEnd,2); 1268c7eeac06SToby Isaac if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) { 1269c7eeac06SToby Isaac ierr = forest->createcellchart(dm,&forest->cStart,&forest->cEnd);CHKERRQ(ierr); 1270c7eeac06SToby Isaac } 1271c7eeac06SToby Isaac *cStart = forest->cStart; 1272c7eeac06SToby Isaac *cEnd = forest->cEnd; 1273c7eeac06SToby Isaac PetscFunctionReturn(0); 1274c7eeac06SToby Isaac } 1275c7eeac06SToby Isaac 1276c7eeac06SToby Isaac #undef __FUNCT__ 1277c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellSF" 12789be51f97SToby Isaac /*@ 12799be51f97SToby Isaac DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes 12809be51f97SToby Isaac 12819be51f97SToby Isaac Not collective 12829be51f97SToby Isaac 12839be51f97SToby Isaac Input Parameter: 12849be51f97SToby Isaac . dm - the forest 12859be51f97SToby Isaac 12869be51f97SToby Isaac Output Parameter: 12879be51f97SToby Isaac . cellSF - the PetscSF 12889be51f97SToby Isaac 12891a244344SSatish Balay Level: intermediate 12909be51f97SToby Isaac 12919be51f97SToby Isaac .seealso: DMForestGetCellChart() 12929be51f97SToby Isaac @*/ 1293c7eeac06SToby Isaac PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF) 1294c7eeac06SToby Isaac { 1295c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1296c7eeac06SToby Isaac PetscErrorCode ierr; 1297c7eeac06SToby Isaac 1298c7eeac06SToby Isaac PetscFunctionBegin; 1299c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1300c7eeac06SToby Isaac PetscValidPointer(cellSF,2); 1301c7eeac06SToby Isaac if ((!forest->cellSF) && forest->createcellsf) { 1302c7eeac06SToby Isaac ierr = forest->createcellsf(dm,&forest->cellSF);CHKERRQ(ierr); 1303c7eeac06SToby Isaac } 1304c7eeac06SToby Isaac *cellSF = forest->cellSF; 1305c7eeac06SToby Isaac PetscFunctionReturn(0); 1306c7eeac06SToby Isaac } 1307c7eeac06SToby Isaac 1308c7eeac06SToby Isaac #undef __FUNCT__ 1309ebdf65a2SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityLabel" 13109be51f97SToby Isaac /*@C 13119be51f97SToby Isaac DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see 13129be51f97SToby Isaac DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination). The 1313*cd3c525cSToby Isaac interpretation of the label values is up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, 1314*cd3c525cSToby Isaac DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have been reserved as choices that should be accepted by all subtypes. 13159be51f97SToby Isaac 13169be51f97SToby Isaac Logically collective on dm 13179be51f97SToby Isaac 13189be51f97SToby Isaac Input Parameters: 13199be51f97SToby Isaac - dm - the forest 1320a1b0c543SToby Isaac + adaptLabel - the label in the pre-adaptation forest 13219be51f97SToby Isaac 13229be51f97SToby Isaac Level: intermediate 13239be51f97SToby Isaac 13249be51f97SToby Isaac .seealso DMForestGetAdaptivityLabel() 13259be51f97SToby Isaac @*/ 1326a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityLabel(DM dm, DMLabel adaptLabel) 1327c7eeac06SToby Isaac { 1328c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1329c7eeac06SToby Isaac PetscErrorCode ierr; 1330c7eeac06SToby Isaac 1331c7eeac06SToby Isaac PetscFunctionBegin; 1332c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1333a1b0c543SToby Isaac adaptLabel->refct++; 1334a1b0c543SToby Isaac if (forest->adaptLabel) {ierr = DMLabelDestroy(&forest->adaptLabel);CHKERRQ(ierr);} 1335a1b0c543SToby Isaac forest->adaptLabel = adaptLabel; 1336c7eeac06SToby Isaac PetscFunctionReturn(0); 1337c7eeac06SToby Isaac } 1338c7eeac06SToby Isaac 1339c7eeac06SToby Isaac #undef __FUNCT__ 1340ebdf65a2SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityLabel" 13419be51f97SToby Isaac /*@C 13429be51f97SToby Isaac DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that 13439be51f97SToby Isaac holds the adaptation flags (refinement, coarsening, or some combination). The interpretation of the label values is 1344*cd3c525cSToby Isaac up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have 1345*cd3c525cSToby Isaac been reserved as choices that should be accepted by all subtypes. 13469be51f97SToby Isaac 13479be51f97SToby Isaac Not collective 13489be51f97SToby Isaac 13499be51f97SToby Isaac Input Parameter: 13509be51f97SToby Isaac . dm - the forest 13519be51f97SToby Isaac 13529be51f97SToby Isaac Output Parameter: 13539be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest 13549be51f97SToby Isaac 13559be51f97SToby Isaac Level: intermediate 13569be51f97SToby Isaac 13579be51f97SToby Isaac .seealso DMForestSetAdaptivityLabel() 13589be51f97SToby Isaac @*/ 1359a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityLabel(DM dm, DMLabel *adaptLabel) 1360c7eeac06SToby Isaac { 1361c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1362c7eeac06SToby Isaac 1363c7eeac06SToby Isaac PetscFunctionBegin; 1364c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1365ba936b91SToby Isaac *adaptLabel = forest->adaptLabel; 1366c7eeac06SToby Isaac PetscFunctionReturn(0); 1367c7eeac06SToby Isaac } 1368c7eeac06SToby Isaac 1369c7eeac06SToby Isaac #undef __FUNCT__ 1370c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetCellWeights" 13719be51f97SToby Isaac /*@ 13729be51f97SToby Isaac DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current 13739be51f97SToby Isaac process: weights are used to determine parallel partitioning. Partitions will be created so that each process's 13749be51f97SToby Isaac ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight 13759be51f97SToby Isaac of 1. 13769be51f97SToby Isaac 13779be51f97SToby Isaac Logically collective on dm 13789be51f97SToby Isaac 13799be51f97SToby Isaac Input Parameters: 13809be51f97SToby Isaac + dm - the forest 13819be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1. 13829be51f97SToby Isaac - copyMode - how weights should reference weights 13839be51f97SToby Isaac 13849be51f97SToby Isaac Level: advanced 13859be51f97SToby Isaac 13869be51f97SToby Isaac .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity() 13879be51f97SToby Isaac @*/ 1388c7eeac06SToby Isaac PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode) 1389c7eeac06SToby Isaac { 1390c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1391c7eeac06SToby Isaac PetscInt cStart, cEnd; 1392c7eeac06SToby Isaac PetscErrorCode ierr; 1393c7eeac06SToby Isaac 1394c7eeac06SToby Isaac PetscFunctionBegin; 1395c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1396c7eeac06SToby Isaac ierr = DMForestGetCellChart(dm,&cStart,&cEnd);CHKERRQ(ierr); 1397c7eeac06SToby Isaac if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd); 1398c7eeac06SToby Isaac if (copyMode == PETSC_COPY_VALUES) { 1399c7eeac06SToby Isaac if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) { 1400c7eeac06SToby Isaac ierr = PetscMalloc1(cEnd-cStart,&forest->cellWeights);CHKERRQ(ierr); 1401c7eeac06SToby Isaac } 1402c7eeac06SToby Isaac ierr = PetscMemcpy(forest->cellWeights,weights,(cEnd-cStart)*sizeof(*weights));CHKERRQ(ierr); 1403c7eeac06SToby Isaac forest->cellWeightsCopyMode = PETSC_OWN_POINTER; 1404c7eeac06SToby Isaac PetscFunctionReturn(0); 1405c7eeac06SToby Isaac } 1406c7eeac06SToby Isaac if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) { 1407c7eeac06SToby Isaac ierr = PetscFree(forest->cellWeights);CHKERRQ(ierr); 1408c7eeac06SToby Isaac } 1409c7eeac06SToby Isaac forest->cellWeights = weights; 1410c7eeac06SToby Isaac forest->cellWeightsCopyMode = copyMode; 1411c7eeac06SToby Isaac PetscFunctionReturn(0); 1412c7eeac06SToby Isaac } 1413c7eeac06SToby Isaac 1414c7eeac06SToby Isaac #undef __FUNCT__ 1415c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellWeights" 14169be51f97SToby Isaac /*@ 14179be51f97SToby Isaac DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current 14189be51f97SToby Isaac process: weights are used to determine parallel partitioning. Partitions will be created so that each process's 14199be51f97SToby Isaac ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight 14209be51f97SToby Isaac of 1. 14219be51f97SToby Isaac 14229be51f97SToby Isaac Not collective 14239be51f97SToby Isaac 14249be51f97SToby Isaac Input Parameter: 14259be51f97SToby Isaac . dm - the forest 14269be51f97SToby Isaac 14279be51f97SToby Isaac Output Parameter: 14289be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1. 14299be51f97SToby Isaac 14309be51f97SToby Isaac Level: advanced 14319be51f97SToby Isaac 14329be51f97SToby Isaac .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity() 14339be51f97SToby Isaac @*/ 1434c7eeac06SToby Isaac PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights) 1435c7eeac06SToby Isaac { 1436c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1437c7eeac06SToby Isaac 1438c7eeac06SToby Isaac PetscFunctionBegin; 1439c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1440c7eeac06SToby Isaac PetscValidPointer(weights,2); 1441c7eeac06SToby Isaac *weights = forest->cellWeights; 1442c7eeac06SToby Isaac PetscFunctionReturn(0); 1443c7eeac06SToby Isaac } 1444c7eeac06SToby Isaac 1445c7eeac06SToby Isaac #undef __FUNCT__ 1446c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetWeightCapacity" 14479be51f97SToby Isaac /*@ 14489be51f97SToby Isaac DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning 14499be51f97SToby Isaac a pre-adaptation forest (see DMForestGetAdaptivityForest()). After partitioning, the ratio of the weight of each 14509be51f97SToby Isaac process's cells to the process's capacity will be roughly equal for all processes. A capacity of 0 indicates that 14519be51f97SToby Isaac the current process should not have any cells after repartitioning. 14529be51f97SToby Isaac 14539be51f97SToby Isaac Logically Collective on dm 14549be51f97SToby Isaac 14559be51f97SToby Isaac Input parameters: 14569be51f97SToby Isaac + dm - the forest 14579be51f97SToby Isaac - capacity - this process's capacity 14589be51f97SToby Isaac 14599be51f97SToby Isaac Level: advanced 14609be51f97SToby Isaac 14619be51f97SToby Isaac .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor() 14629be51f97SToby Isaac @*/ 1463c7eeac06SToby Isaac PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity) 1464c7eeac06SToby Isaac { 1465c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1466c7eeac06SToby Isaac 1467c7eeac06SToby Isaac PetscFunctionBegin; 1468c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1469ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup"); 1470c7eeac06SToby Isaac if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity); 1471c7eeac06SToby Isaac forest->weightCapacity = capacity; 1472c7eeac06SToby Isaac PetscFunctionReturn(0); 1473c7eeac06SToby Isaac } 1474c7eeac06SToby Isaac 1475c7eeac06SToby Isaac #undef __FUNCT__ 1476c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetWeightCapacity" 14779be51f97SToby Isaac /*@ 14789be51f97SToby Isaac DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see 14799be51f97SToby Isaac DMForestGetAdaptivityForest()). After partitioning, the ratio of the weight of each process's cells to the 14809be51f97SToby Isaac process's capacity will be roughly equal for all processes. A capacity of 0 indicates that the current process 14819be51f97SToby Isaac should not have any cells after repartitioning. 14829be51f97SToby Isaac 14839be51f97SToby Isaac Not collective 14849be51f97SToby Isaac 14859be51f97SToby Isaac Input parameter: 14869be51f97SToby Isaac . dm - the forest 14879be51f97SToby Isaac 14889be51f97SToby Isaac Output parameter: 14899be51f97SToby Isaac . capacity - this process's capacity 14909be51f97SToby Isaac 14919be51f97SToby Isaac Level: advanced 14929be51f97SToby Isaac 14939be51f97SToby Isaac .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor() 14949be51f97SToby Isaac @*/ 1495c7eeac06SToby Isaac PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity) 1496c7eeac06SToby Isaac { 1497c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1498c7eeac06SToby Isaac 1499c7eeac06SToby Isaac PetscFunctionBegin; 1500c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1501c7eeac06SToby Isaac PetscValidRealPointer(capacity,2); 1502c7eeac06SToby Isaac *capacity = forest->weightCapacity; 1503c7eeac06SToby Isaac PetscFunctionReturn(0); 1504c7eeac06SToby Isaac } 1505c7eeac06SToby Isaac 1506dd8e54a2SToby Isaac #undef __FUNCT__ 1507db4d5e8cSToby Isaac #define __FUNCT__ "DMSetFromOptions_Forest" 15080709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm) 1509db4d5e8cSToby Isaac { 1510db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 151156ba9f64SToby Isaac PetscBool flg, flg1, flg2, flg3, flg4; 1512dd8e54a2SToby Isaac DMForestTopology oldTopo; 1513c7eeac06SToby Isaac char stringBuffer[256]; 1514dd8e54a2SToby Isaac PetscViewer viewer; 1515dd8e54a2SToby Isaac PetscViewerFormat format; 151656ba9f64SToby Isaac PetscInt adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade; 1517c7eeac06SToby Isaac PetscReal weightsFactor; 1518c7eeac06SToby Isaac DMForestAdaptivityStrategy adaptStrategy; 1519db4d5e8cSToby Isaac PetscErrorCode ierr; 1520db4d5e8cSToby Isaac 1521db4d5e8cSToby Isaac PetscFunctionBegin; 1522db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 152358762b62SToby Isaac forest->setfromoptionscalled = PETSC_TRUE; 1524dd8e54a2SToby Isaac ierr = DMForestGetTopology(dm, &oldTopo);CHKERRQ(ierr); 1525a3eda1eaSToby Isaac ierr = PetscOptionsHead(PetscOptionsObject,"DMForest Options");CHKERRQ(ierr); 152656ba9f64SToby Isaac ierr = PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);CHKERRQ(ierr); 152756ba9f64SToby Isaac ierr = PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);CHKERRQ(ierr); 152856ba9f64SToby Isaac ierr = PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);CHKERRQ(ierr); 152956ba9f64SToby Isaac ierr = PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);CHKERRQ(ierr); 1530f885a11aSToby Isaac if ((PetscInt) flg1 + (PetscInt) flg2 + (PetscInt) flg3 + (PetscInt) flg4 > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Specify only one of -dm_forest_{topology,base_dm,coarse_forest,fine_forest}"); 153156ba9f64SToby Isaac if (flg1) { 153256ba9f64SToby Isaac ierr = DMForestSetTopology(dm,(DMForestTopology)stringBuffer);CHKERRQ(ierr); 153356ba9f64SToby Isaac ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr); 153420e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr); 153556ba9f64SToby Isaac } 153656ba9f64SToby Isaac if (flg2) { 1537dd8e54a2SToby Isaac DM base; 1538dd8e54a2SToby Isaac 1539dd8e54a2SToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),&base);CHKERRQ(ierr); 1540dd8e54a2SToby Isaac ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 1541dd8e54a2SToby Isaac ierr = DMLoad(base,viewer);CHKERRQ(ierr); 1542dd8e54a2SToby Isaac ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1543dd8e54a2SToby Isaac ierr = DMForestSetBaseDM(dm,base);CHKERRQ(ierr); 1544dd8e54a2SToby Isaac ierr = DMDestroy(&base);CHKERRQ(ierr); 154556ba9f64SToby Isaac ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr); 154620e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr); 1547dd8e54a2SToby Isaac } 154856ba9f64SToby Isaac if (flg3) { 1549dd8e54a2SToby Isaac DM coarse; 1550dd8e54a2SToby Isaac 1551dd8e54a2SToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),&coarse);CHKERRQ(ierr); 1552dd8e54a2SToby Isaac ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 1553dd8e54a2SToby Isaac ierr = DMLoad(coarse,viewer);CHKERRQ(ierr); 1554dd8e54a2SToby Isaac ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 155520e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,coarse);CHKERRQ(ierr); 1556dd8e54a2SToby Isaac ierr = DMDestroy(&coarse);CHKERRQ(ierr); 155756ba9f64SToby Isaac ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr); 155856ba9f64SToby Isaac ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr); 1559dd8e54a2SToby Isaac } 156056ba9f64SToby Isaac if (flg4) { 1561dd8e54a2SToby Isaac DM fine; 1562dd8e54a2SToby Isaac 1563dd8e54a2SToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),&fine);CHKERRQ(ierr); 1564dd8e54a2SToby Isaac ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 1565dd8e54a2SToby Isaac ierr = DMLoad(fine,viewer);CHKERRQ(ierr); 1566dd8e54a2SToby Isaac ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 156720e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,fine);CHKERRQ(ierr); 1568dd8e54a2SToby Isaac ierr = DMDestroy(&fine);CHKERRQ(ierr); 156956ba9f64SToby Isaac ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr); 157056ba9f64SToby Isaac ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr); 1571dd8e54a2SToby Isaac } 1572dd8e54a2SToby Isaac ierr = DMForestGetAdjacencyDimension(dm,&adjDim);CHKERRQ(ierr); 1573dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg);CHKERRQ(ierr); 1574dd8e54a2SToby Isaac if (flg) { 1575dd8e54a2SToby Isaac ierr = DMForestSetAdjacencyDimension(dm,adjDim);CHKERRQ(ierr); 1576f885a11aSToby Isaac } else { 1577dd8e54a2SToby Isaac ierr = DMForestGetAdjacencyCodimension(dm,&adjCodim);CHKERRQ(ierr); 1578dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg);CHKERRQ(ierr); 1579dd8e54a2SToby Isaac if (flg) { 1580dd8e54a2SToby Isaac ierr = DMForestSetAdjacencyCodimension(dm,adjCodim);CHKERRQ(ierr); 1581dd8e54a2SToby Isaac } 1582dd8e54a2SToby Isaac } 1583dd8e54a2SToby Isaac ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr); 1584dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg);CHKERRQ(ierr); 1585dd8e54a2SToby Isaac if (flg) { 1586dd8e54a2SToby Isaac ierr = DMForestSetPartitionOverlap(dm,overlap);CHKERRQ(ierr); 1587dd8e54a2SToby Isaac } 1588a6121fbdSMatthew G. Knepley #if 0 1589a6121fbdSMatthew 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); 1590a6121fbdSMatthew G. Knepley if (flg) { 1591a6121fbdSMatthew G. Knepley ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr); 1592a6121fbdSMatthew G. Knepley ierr = DMForestSetInitialRefinement(dm,minRefinement);CHKERRQ(ierr); 1593a6121fbdSMatthew G. Knepley } 1594a6121fbdSMatthew 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); 1595a6121fbdSMatthew G. Knepley if (flg) { 1596a6121fbdSMatthew G. Knepley ierr = DMForestSetMinimumRefinement(dm,0);CHKERRQ(ierr); 1597a6121fbdSMatthew G. Knepley ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr); 1598a6121fbdSMatthew G. Knepley } 1599a6121fbdSMatthew G. Knepley #endif 1600dd8e54a2SToby Isaac ierr = DMForestGetMinimumRefinement(dm,&minRefinement);CHKERRQ(ierr); 1601dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg);CHKERRQ(ierr); 1602dd8e54a2SToby Isaac if (flg) { 1603dd8e54a2SToby Isaac ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr); 1604db4d5e8cSToby Isaac } 160556ba9f64SToby Isaac ierr = DMForestGetInitialRefinement(dm,&initRefinement);CHKERRQ(ierr); 160656ba9f64SToby Isaac ierr = PetscOptionsInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg);CHKERRQ(ierr); 160756ba9f64SToby Isaac if (flg) { 160856ba9f64SToby Isaac ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr); 160956ba9f64SToby Isaac } 1610c7eeac06SToby Isaac ierr = DMForestGetMaximumRefinement(dm,&maxRefinement);CHKERRQ(ierr); 1611c7eeac06SToby Isaac ierr = PetscOptionsInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg);CHKERRQ(ierr); 1612c7eeac06SToby Isaac if (flg) { 1613c7eeac06SToby Isaac ierr = DMForestSetMaximumRefinement(dm,maxRefinement);CHKERRQ(ierr); 1614c7eeac06SToby Isaac } 1615c7eeac06SToby Isaac ierr = DMForestGetAdaptivityStrategy(dm,&adaptStrategy);CHKERRQ(ierr); 1616c7eeac06SToby Isaac ierr = PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,256,&flg);CHKERRQ(ierr); 1617c7eeac06SToby Isaac if (flg) { 1618c7eeac06SToby Isaac ierr = DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);CHKERRQ(ierr); 1619c7eeac06SToby Isaac } 1620c7eeac06SToby Isaac ierr = DMForestGetGradeFactor(dm,&grade);CHKERRQ(ierr); 1621c7eeac06SToby Isaac ierr = PetscOptionsInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg);CHKERRQ(ierr); 1622c7eeac06SToby Isaac if (flg) { 1623c7eeac06SToby Isaac ierr = DMForestSetGradeFactor(dm,grade);CHKERRQ(ierr); 1624c7eeac06SToby Isaac } 1625c7eeac06SToby Isaac ierr = DMForestGetCellWeightFactor(dm,&weightsFactor);CHKERRQ(ierr); 1626c7eeac06SToby Isaac ierr = PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);CHKERRQ(ierr); 1627c7eeac06SToby Isaac if (flg) { 1628c7eeac06SToby Isaac ierr = DMForestSetCellWeightFactor(dm,weightsFactor);CHKERRQ(ierr); 1629c7eeac06SToby Isaac } 1630db4d5e8cSToby Isaac ierr = PetscOptionsTail();CHKERRQ(ierr); 1631db4d5e8cSToby Isaac PetscFunctionReturn(0); 1632db4d5e8cSToby Isaac } 1633db4d5e8cSToby Isaac 1634db4d5e8cSToby Isaac #undef __FUNCT__ 1635d8984e3bSMatthew G. Knepley #define __FUNCT__ "DMCreateSubDM_Forest" 1636d8984e3bSMatthew G. Knepley PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm) 1637d8984e3bSMatthew G. Knepley { 1638d8984e3bSMatthew G. Knepley PetscErrorCode ierr; 1639d8984e3bSMatthew G. Knepley 1640d8984e3bSMatthew G. Knepley PetscFunctionBegin; 1641d8984e3bSMatthew G. Knepley if (subdm) {ierr = DMClone(dm, subdm);CHKERRQ(ierr);} 1642d8984e3bSMatthew G. Knepley ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr); 1643d8984e3bSMatthew G. Knepley PetscFunctionReturn(0); 1644d8984e3bSMatthew G. Knepley } 1645d8984e3bSMatthew G. Knepley 1646d8984e3bSMatthew G. Knepley #undef __FUNCT__ 16475421bac9SToby Isaac #define __FUNCT__ "DMRefine_Forest" 16485421bac9SToby Isaac PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined) 16495421bac9SToby Isaac { 16505421bac9SToby Isaac DMLabel refine; 16515421bac9SToby Isaac DM fineDM; 16525421bac9SToby Isaac PetscErrorCode ierr; 16535421bac9SToby Isaac 16545421bac9SToby Isaac PetscFunctionBegin; 16555421bac9SToby Isaac ierr = DMGetFineDM(dm,&fineDM);CHKERRQ(ierr); 16565421bac9SToby Isaac if (fineDM) { 16575421bac9SToby Isaac ierr = PetscObjectReference((PetscObject)fineDM);CHKERRQ(ierr); 16585421bac9SToby Isaac *dmRefined = fineDM; 16595421bac9SToby Isaac PetscFunctionReturn(0); 16605421bac9SToby Isaac } 16615421bac9SToby Isaac ierr = DMForestTemplate(dm,comm,dmRefined);CHKERRQ(ierr); 16625421bac9SToby Isaac ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr); 16635421bac9SToby Isaac if (!refine) { 1664a1b0c543SToby Isaac ierr = DMLabelCreate("refine",&refine);CHKERRQ(ierr); 1665a1b0c543SToby Isaac ierr = DMLabelSetDefaultValue(refine,DM_ADAPT_REFINE);CHKERRQ(ierr); 16665421bac9SToby Isaac } 1667a1b0c543SToby Isaac else { 1668a1b0c543SToby Isaac refine->refct++; 1669a1b0c543SToby Isaac } 1670a1b0c543SToby Isaac ierr = DMForestSetAdaptivityLabel(*dmRefined,refine);CHKERRQ(ierr); 1671a1b0c543SToby Isaac ierr = DMLabelDestroy(&refine);CHKERRQ(ierr); 16725421bac9SToby Isaac PetscFunctionReturn(0); 16735421bac9SToby Isaac } 16745421bac9SToby Isaac 16755421bac9SToby Isaac #undef __FUNCT__ 16765421bac9SToby Isaac #define __FUNCT__ "DMCoarsen_Forest" 16775421bac9SToby Isaac PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened) 16785421bac9SToby Isaac { 16795421bac9SToby Isaac DMLabel coarsen; 16805421bac9SToby Isaac DM coarseDM; 16815421bac9SToby Isaac PetscErrorCode ierr; 16825421bac9SToby Isaac 16835421bac9SToby Isaac PetscFunctionBegin; 16844098eed7SToby Isaac { 16854098eed7SToby Isaac PetscMPIInt mpiComparison; 16864098eed7SToby Isaac MPI_Comm dmcomm = PetscObjectComm((PetscObject)dm); 16874098eed7SToby Isaac 16884098eed7SToby Isaac ierr = MPI_Comm_compare(comm, dmcomm, &mpiComparison);CHKERRQ(ierr); 1689f885a11aSToby Isaac if (mpiComparison != MPI_IDENT && mpiComparison != MPI_CONGRUENT) SETERRQ(dmcomm,PETSC_ERR_SUP,"No support for different communicators yet"); 16904098eed7SToby Isaac } 16915421bac9SToby Isaac ierr = DMGetCoarseDM(dm,&coarseDM);CHKERRQ(ierr); 16925421bac9SToby Isaac if (coarseDM) { 16935421bac9SToby Isaac ierr = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr); 16945421bac9SToby Isaac *dmCoarsened = coarseDM; 16955421bac9SToby Isaac PetscFunctionReturn(0); 16965421bac9SToby Isaac } 16975421bac9SToby Isaac ierr = DMForestTemplate(dm,comm,dmCoarsened);CHKERRQ(ierr); 1698a1b0c543SToby Isaac ierr = DMForestSetAdaptivityPurpose(coarseDM,DM_ADAPT_COARSEN);CHKERRQ(ierr); 16995421bac9SToby Isaac ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr); 17005421bac9SToby Isaac if (!coarsen) { 1701a1b0c543SToby Isaac ierr = DMLabelCreate("coarsen",&coarsen);CHKERRQ(ierr); 1702a1b0c543SToby Isaac ierr = DMLabelSetDefaultValue(coarsen,DM_ADAPT_COARSEN);CHKERRQ(ierr); 1703a1b0c543SToby Isaac } else { 1704a1b0c543SToby Isaac coarsen->refct++; 17055421bac9SToby Isaac } 1706a1b0c543SToby Isaac ierr = DMForestSetAdaptivityLabel(*dmCoarsened,coarsen);CHKERRQ(ierr); 1707a1b0c543SToby Isaac ierr = DMLabelDestroy(&coarsen);CHKERRQ(ierr); 17085421bac9SToby Isaac PetscFunctionReturn(0); 17095421bac9SToby Isaac } 17105421bac9SToby Isaac 17115421bac9SToby Isaac #undef __FUNCT__ 1712a1b0c543SToby Isaac #define __FUNCT__ "DMAdaptLabel_Forest" 1713a1b0c543SToby Isaac static PetscErrorCode DMAdaptLabel_Forest(DM dm, DMLabel label, DM *adaptedDM) 171409350103SToby Isaac { 171509350103SToby Isaac PetscBool success; 171609350103SToby Isaac PetscErrorCode ierr; 171709350103SToby Isaac 171809350103SToby Isaac PetscFunctionBegin; 171909350103SToby Isaac ierr = DMForestTemplate(dm,PetscObjectComm((PetscObject)dm),adaptedDM);CHKERRQ(ierr); 1720a1b0c543SToby Isaac ierr = DMForestSetAdaptivityLabel(*adaptedDM,label);CHKERRQ(ierr); 172109350103SToby Isaac ierr = DMSetUp(*adaptedDM);CHKERRQ(ierr); 172209350103SToby Isaac ierr = DMForestGetAdaptivitySuccess(*adaptedDM,&success);CHKERRQ(ierr); 172309350103SToby Isaac if (!success) { 172409350103SToby Isaac ierr = DMDestroy(adaptedDM);CHKERRQ(ierr); 172509350103SToby Isaac *adaptedDM = NULL; 172609350103SToby Isaac } 172709350103SToby Isaac PetscFunctionReturn(0); 172809350103SToby Isaac } 172909350103SToby Isaac 173009350103SToby Isaac #undef __FUNCT__ 1731d222f98bSToby Isaac #define __FUNCT__ "DMInitialize_Forest" 1732d222f98bSToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm) 1733d222f98bSToby Isaac { 1734d222f98bSToby Isaac PetscErrorCode ierr; 1735d222f98bSToby Isaac 1736d222f98bSToby Isaac PetscFunctionBegin; 1737d222f98bSToby Isaac ierr = PetscMemzero(dm->ops,sizeof(*(dm->ops)));CHKERRQ(ierr); 1738d222f98bSToby Isaac 1739d222f98bSToby Isaac dm->ops->clone = DMClone_Forest; 1740d222f98bSToby Isaac dm->ops->setfromoptions = DMSetFromOptions_Forest; 1741d222f98bSToby Isaac dm->ops->destroy = DMDestroy_Forest; 1742d8984e3bSMatthew G. Knepley dm->ops->createsubdm = DMCreateSubDM_Forest; 17435421bac9SToby Isaac dm->ops->refine = DMRefine_Forest; 17445421bac9SToby Isaac dm->ops->coarsen = DMCoarsen_Forest; 1745a1b0c543SToby Isaac ierr = PetscObjectComposeFunction((PetscObject)dm,"DMAdaptLabel_C",DMAdaptLabel_Forest);CHKERRQ(ierr); 1746d222f98bSToby Isaac PetscFunctionReturn(0); 1747d222f98bSToby Isaac } 1748d222f98bSToby Isaac 17499be51f97SToby Isaac /*MC 17509be51f97SToby 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. 17519be51f97SToby Isaac 1752a1b0c543SToby 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_ADAPT_REFINE) or coarsened (DM_ADAPT_COARSEN) and how (subtypes are free to allow additional values for things like anisotropic refinement). The label should be given to the *new* mesh with DMForestSetAdaptivityLabel(). 17539be51f97SToby Isaac 17549be51f97SToby Isaac Level: advanced 17559be51f97SToby Isaac 17569be51f97SToby Isaac .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel() 17579be51f97SToby Isaac M*/ 17589be51f97SToby Isaac 1759d222f98bSToby Isaac #undef __FUNCT__ 1760db4d5e8cSToby Isaac #define __FUNCT__ "DMCreate_Forest" 1761db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm) 1762db4d5e8cSToby Isaac { 1763db4d5e8cSToby Isaac DM_Forest *forest; 1764db4d5e8cSToby Isaac PetscErrorCode ierr; 1765db4d5e8cSToby Isaac 1766db4d5e8cSToby Isaac PetscFunctionBegin; 1767db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1768db4d5e8cSToby Isaac ierr = PetscNewLog(dm,&forest);CHKERRQ(ierr); 1769db4d5e8cSToby Isaac dm->dim = 0; 1770db4d5e8cSToby Isaac dm->data = forest; 1771db4d5e8cSToby Isaac forest->refct = 1; 1772db4d5e8cSToby Isaac forest->data = NULL; 177358762b62SToby Isaac forest->setfromoptionscalled = PETSC_FALSE; 1774db4d5e8cSToby Isaac forest->topology = NULL; 1775*cd3c525cSToby Isaac forest->adapt = NULL; 1776db4d5e8cSToby Isaac forest->base = NULL; 1777*cd3c525cSToby Isaac forest->adaptPurpose = PETSC_DETERMINE; 1778db4d5e8cSToby Isaac forest->adjDim = PETSC_DEFAULT; 1779db4d5e8cSToby Isaac forest->overlap = PETSC_DEFAULT; 1780db4d5e8cSToby Isaac forest->minRefinement = PETSC_DEFAULT; 1781db4d5e8cSToby Isaac forest->maxRefinement = PETSC_DEFAULT; 178256ba9f64SToby Isaac forest->initRefinement = PETSC_DEFAULT; 1783c7eeac06SToby Isaac forest->cStart = PETSC_DETERMINE; 1784c7eeac06SToby Isaac forest->cEnd = PETSC_DETERMINE; 1785*cd3c525cSToby Isaac forest->cellSF = NULL; 1786ebdf65a2SToby Isaac forest->adaptLabel = NULL; 1787db4d5e8cSToby Isaac forest->gradeFactor = 2; 1788db4d5e8cSToby Isaac forest->cellWeights = NULL; 1789db4d5e8cSToby Isaac forest->cellWeightsCopyMode = PETSC_USE_POINTER; 1790db4d5e8cSToby Isaac forest->weightsFactor = 1.; 1791db4d5e8cSToby Isaac forest->weightCapacity = 1.; 1792a73e2921SToby Isaac ierr = DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);CHKERRQ(ierr); 1793d222f98bSToby Isaac ierr = DMInitialize_Forest(dm);CHKERRQ(ierr); 1794db4d5e8cSToby Isaac PetscFunctionReturn(0); 1795db4d5e8cSToby Isaac } 1796db4d5e8cSToby Isaac 1797