19be51f97SToby Isaac #include <petsc/private/dmforestimpl.h> /*I "petscdmforest.h" I*/ 29be51f97SToby Isaac #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 3ef19d27cSToby Isaac #include <petscsf.h> 4db4d5e8cSToby Isaac 527d4645fSToby Isaac PetscBool DMForestPackageInitialized = PETSC_FALSE; 627d4645fSToby Isaac 727d4645fSToby Isaac typedef struct _DMForestTypeLink*DMForestTypeLink; 827d4645fSToby Isaac 927d4645fSToby Isaac struct _DMForestTypeLink 1027d4645fSToby Isaac { 1127d4645fSToby Isaac char *name; 1227d4645fSToby Isaac DMForestTypeLink next; 1327d4645fSToby Isaac }; 1427d4645fSToby Isaac 1527d4645fSToby Isaac DMForestTypeLink DMForestTypeList; 1627d4645fSToby Isaac 1727d4645fSToby Isaac #undef __FUNCT__ 1827d4645fSToby Isaac #define __FUNCT__ "DMForestPackageFinalize" 1927d4645fSToby Isaac static PetscErrorCode DMForestPackageFinalize(void) 2027d4645fSToby Isaac { 2127d4645fSToby Isaac DMForestTypeLink oldLink, link = DMForestTypeList; 2227d4645fSToby Isaac PetscErrorCode ierr; 2327d4645fSToby Isaac 2427d4645fSToby Isaac PetscFunctionBegin; 2527d4645fSToby Isaac while (link) { 2627d4645fSToby Isaac oldLink = link; 27f416af30SBarry Smith ierr = PetscFree(oldLink->name);CHKERRQ(ierr); 2827d4645fSToby Isaac link = oldLink->next; 2927d4645fSToby Isaac ierr = PetscFree(oldLink);CHKERRQ(ierr); 3027d4645fSToby Isaac } 3127d4645fSToby Isaac PetscFunctionReturn(0); 3227d4645fSToby Isaac } 3327d4645fSToby Isaac 3427d4645fSToby Isaac #undef __FUNCT__ 3527d4645fSToby Isaac #define __FUNCT__ "DMForestPackageInitialize" 3627d4645fSToby Isaac static PetscErrorCode DMForestPackageInitialize(void) 3727d4645fSToby Isaac { 3827d4645fSToby Isaac PetscErrorCode ierr; 3927d4645fSToby Isaac 4027d4645fSToby Isaac PetscFunctionBegin; 4127d4645fSToby Isaac if (DMForestPackageInitialized) PetscFunctionReturn(0); 4227d4645fSToby Isaac DMForestPackageInitialized = PETSC_TRUE; 43f885a11aSToby Isaac 4427d4645fSToby Isaac ierr = DMForestRegisterType(DMFOREST);CHKERRQ(ierr); 4527d4645fSToby Isaac ierr = PetscRegisterFinalize(DMForestPackageFinalize);CHKERRQ(ierr); 4627d4645fSToby Isaac PetscFunctionReturn(0); 4727d4645fSToby Isaac } 4827d4645fSToby Isaac 4927d4645fSToby Isaac #undef __FUNCT__ 5027d4645fSToby Isaac #define __FUNCT__ "DMForestRegisterType" 519be51f97SToby Isaac /*@C 529be51f97SToby Isaac DMForestRegisterType - Registers a DMType as a subtype of DMFOREST (so that DMIsForest() will be correct) 539be51f97SToby Isaac 549be51f97SToby Isaac Not Collective 559be51f97SToby Isaac 569be51f97SToby Isaac Input parameter: 579be51f97SToby Isaac . name - the name of the type 589be51f97SToby Isaac 599be51f97SToby Isaac Level: advanced 609be51f97SToby Isaac 619be51f97SToby Isaac .seealso: DMFOREST, DMIsForest() 629be51f97SToby Isaac @*/ 6327d4645fSToby Isaac PetscErrorCode DMForestRegisterType(DMType name) 6427d4645fSToby Isaac { 6527d4645fSToby Isaac DMForestTypeLink link; 6627d4645fSToby Isaac PetscErrorCode ierr; 6727d4645fSToby Isaac 6827d4645fSToby Isaac PetscFunctionBegin; 6927d4645fSToby Isaac ierr = DMForestPackageInitialize();CHKERRQ(ierr); 7027d4645fSToby Isaac ierr = PetscNew(&link);CHKERRQ(ierr); 7127d4645fSToby Isaac ierr = PetscStrallocpy(name,&link->name);CHKERRQ(ierr); 7227d4645fSToby Isaac link->next = DMForestTypeList; 7327d4645fSToby Isaac DMForestTypeList = link; 7427d4645fSToby Isaac PetscFunctionReturn(0); 7527d4645fSToby Isaac } 7627d4645fSToby Isaac 7727d4645fSToby Isaac #undef __FUNCT__ 7827d4645fSToby Isaac #define __FUNCT__ "DMIsForest" 799be51f97SToby Isaac /*@ 809be51f97SToby Isaac DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes 819be51f97SToby Isaac 829be51f97SToby Isaac Not Collective 839be51f97SToby Isaac 849be51f97SToby Isaac Input parameter: 859be51f97SToby Isaac . dm - the DM object 869be51f97SToby Isaac 879be51f97SToby Isaac Output parameter: 889be51f97SToby Isaac . isForest - whether dm is a subtype of DMFOREST 899be51f97SToby Isaac 909be51f97SToby Isaac Level: intermediate 919be51f97SToby Isaac 929be51f97SToby Isaac .seealso: DMFOREST, DMForestRegisterType() 939be51f97SToby Isaac @*/ 9427d4645fSToby Isaac PetscErrorCode DMIsForest(DM dm, PetscBool *isForest) 9527d4645fSToby Isaac { 9627d4645fSToby Isaac DMForestTypeLink link = DMForestTypeList; 9727d4645fSToby Isaac PetscErrorCode ierr; 9827d4645fSToby Isaac 9927d4645fSToby Isaac PetscFunctionBegin; 10027d4645fSToby Isaac while (link) { 10127d4645fSToby Isaac PetscBool sameType; 10227d4645fSToby Isaac ierr = PetscObjectTypeCompare((PetscObject)dm,link->name,&sameType);CHKERRQ(ierr); 10327d4645fSToby Isaac if (sameType) { 10427d4645fSToby Isaac *isForest = PETSC_TRUE; 10527d4645fSToby Isaac PetscFunctionReturn(0); 10627d4645fSToby Isaac } 10727d4645fSToby Isaac link = link->next; 10827d4645fSToby Isaac } 10927d4645fSToby Isaac *isForest = PETSC_FALSE; 11027d4645fSToby Isaac PetscFunctionReturn(0); 11127d4645fSToby Isaac } 11227d4645fSToby Isaac 113db4d5e8cSToby Isaac #undef __FUNCT__ 114a0452a8eSToby Isaac #define __FUNCT__ "DMForestTemplate" 1159be51f97SToby Isaac /*@ 1169be51f97SToby Isaac DMForestTemplate - Create a new DM that will be adapted from a source DM. The new DM reproduces the configuration 1179be51f97SToby 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 1189be51f97SToby Isaac (by, e.g., refinement or repartitioning). The source DM is also set as the adaptivity source DM of the new DM (see 1199be51f97SToby Isaac DMForestSetAdaptivityForest()). 1209be51f97SToby Isaac 1219be51f97SToby Isaac Collective on dm 1229be51f97SToby Isaac 1239be51f97SToby Isaac Input Parameters: 1249be51f97SToby Isaac + dm - the source DM object 1259be51f97SToby Isaac - comm - the communicator for the new DM (this communicator is currently ignored, but is present so that DMForestTemplate() can be used within DMCoarsen()) 1269be51f97SToby Isaac 1279be51f97SToby Isaac Output Parameter: 1289be51f97SToby Isaac . tdm - the new DM object 1299be51f97SToby Isaac 1309be51f97SToby Isaac Level: intermediate 1319be51f97SToby Isaac 1329be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest() 1339be51f97SToby Isaac @*/ 1349be51f97SToby Isaac PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm) 135a0452a8eSToby Isaac { 136a0452a8eSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 13720e8089bSToby Isaac DMType type; 138a0452a8eSToby Isaac DM base; 139a0452a8eSToby Isaac DMForestTopology topology; 140a0452a8eSToby Isaac PetscInt dim, overlap, ref, factor; 141a0452a8eSToby Isaac DMForestAdaptivityStrategy strat; 142795844e7SToby Isaac PetscDS ds; 143795844e7SToby Isaac void *ctx; 14449fc9a2fSToby Isaac PetscErrorCode (*map)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void*); 1453e58adeeSToby Isaac void *mapCtx; 146a0452a8eSToby Isaac PetscErrorCode ierr; 147a0452a8eSToby Isaac 148a0452a8eSToby Isaac PetscFunctionBegin; 149a0452a8eSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 15020e8089bSToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),tdm);CHKERRQ(ierr); 15120e8089bSToby Isaac ierr = DMGetType(dm,&type);CHKERRQ(ierr); 15220e8089bSToby Isaac ierr = DMSetType(*tdm,type);CHKERRQ(ierr); 153a0452a8eSToby Isaac ierr = DMForestGetBaseDM(dm,&base);CHKERRQ(ierr); 15420e8089bSToby Isaac ierr = DMForestSetBaseDM(*tdm,base);CHKERRQ(ierr); 155a0452a8eSToby Isaac ierr = DMForestGetTopology(dm,&topology);CHKERRQ(ierr); 15620e8089bSToby Isaac ierr = DMForestSetTopology(*tdm,topology);CHKERRQ(ierr); 157a0452a8eSToby Isaac ierr = DMForestGetAdjacencyDimension(dm,&dim);CHKERRQ(ierr); 15820e8089bSToby Isaac ierr = DMForestSetAdjacencyDimension(*tdm,dim);CHKERRQ(ierr); 159a0452a8eSToby Isaac ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr); 16020e8089bSToby Isaac ierr = DMForestSetPartitionOverlap(*tdm,overlap);CHKERRQ(ierr); 161a0452a8eSToby Isaac ierr = DMForestGetMinimumRefinement(dm,&ref);CHKERRQ(ierr); 16220e8089bSToby Isaac ierr = DMForestSetMinimumRefinement(*tdm,ref);CHKERRQ(ierr); 163a0452a8eSToby Isaac ierr = DMForestGetMaximumRefinement(dm,&ref);CHKERRQ(ierr); 16420e8089bSToby Isaac ierr = DMForestSetMaximumRefinement(*tdm,ref);CHKERRQ(ierr); 165a0452a8eSToby Isaac ierr = DMForestGetAdaptivityStrategy(dm,&strat);CHKERRQ(ierr); 16620e8089bSToby Isaac ierr = DMForestSetAdaptivityStrategy(*tdm,strat);CHKERRQ(ierr); 167a0452a8eSToby Isaac ierr = DMForestGetGradeFactor(dm,&factor);CHKERRQ(ierr); 16820e8089bSToby Isaac ierr = DMForestSetGradeFactor(*tdm,factor);CHKERRQ(ierr); 1693e58adeeSToby Isaac ierr = DMForestGetBaseCoordinateMapping(dm,&map,&mapCtx);CHKERRQ(ierr); 1705e8c540aSToby Isaac ierr = DMForestSetBaseCoordinateMapping(*tdm,map,mapCtx);CHKERRQ(ierr); 171a0452a8eSToby Isaac if (forest->ftemplate) { 17220e8089bSToby Isaac ierr = (forest->ftemplate)(dm, *tdm);CHKERRQ(ierr); 173a0452a8eSToby Isaac } 17420e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(*tdm,dm);CHKERRQ(ierr); 175795844e7SToby Isaac ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 176795844e7SToby Isaac ierr = DMSetDS(*tdm,ds);CHKERRQ(ierr); 177795844e7SToby Isaac ierr = DMGetApplicationContext(dm,&ctx);CHKERRQ(ierr); 178795844e7SToby Isaac ierr = DMSetApplicationContext(*tdm,&ctx);CHKERRQ(ierr); 179795844e7SToby Isaac if (dm->maxCell) { 180795844e7SToby Isaac const PetscReal *maxCell, *L; 181795844e7SToby Isaac const DMBoundaryType *bd; 182795844e7SToby Isaac 183795844e7SToby Isaac ierr = DMGetPeriodicity(dm,&maxCell,&L,&bd);CHKERRQ(ierr); 184795844e7SToby Isaac ierr = DMSetPeriodicity(*tdm,maxCell,L,bd);CHKERRQ(ierr); 185795844e7SToby Isaac } 186a0452a8eSToby Isaac PetscFunctionReturn(0); 187a0452a8eSToby Isaac } 188a0452a8eSToby Isaac 18901d9d024SToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm); 19001d9d024SToby Isaac 191a0452a8eSToby Isaac #undef __FUNCT__ 192db4d5e8cSToby Isaac #define __FUNCT__ "DMClone_Forest" 193db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm) 194db4d5e8cSToby Isaac { 195db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 196db4d5e8cSToby Isaac const char *type; 197db4d5e8cSToby Isaac PetscErrorCode ierr; 198db4d5e8cSToby Isaac 199db4d5e8cSToby Isaac PetscFunctionBegin; 200db4d5e8cSToby Isaac forest->refct++; 201db4d5e8cSToby Isaac (*newdm)->data = forest; 202db4d5e8cSToby Isaac ierr = PetscObjectGetType((PetscObject) dm, &type);CHKERRQ(ierr); 203db4d5e8cSToby Isaac ierr = PetscObjectChangeTypeName((PetscObject) *newdm, type);CHKERRQ(ierr); 20401d9d024SToby Isaac ierr = DMInitialize_Forest(*newdm);CHKERRQ(ierr); 205db4d5e8cSToby Isaac PetscFunctionReturn(0); 206db4d5e8cSToby Isaac } 207db4d5e8cSToby Isaac 208db4d5e8cSToby Isaac #undef __FUNCT__ 209db4d5e8cSToby Isaac #define __FUNCT__ "DMDestroy_Forest" 210d222f98bSToby Isaac static PetscErrorCode DMDestroy_Forest(DM dm) 211db4d5e8cSToby Isaac { 212db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 213db4d5e8cSToby Isaac PetscErrorCode ierr; 214db4d5e8cSToby Isaac 215db4d5e8cSToby Isaac PetscFunctionBegin; 216db4d5e8cSToby Isaac if (--forest->refct > 0) PetscFunctionReturn(0); 217d222f98bSToby Isaac if (forest->destroy) {ierr = forest->destroy(dm);CHKERRQ(ierr);} 218db4d5e8cSToby Isaac ierr = PetscSFDestroy(&forest->cellSF);CHKERRQ(ierr); 2190f17b9e3SToby Isaac ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr); 2200f17b9e3SToby Isaac ierr = PetscSFDestroy(&forest->coarseToPreFine);CHKERRQ(ierr); 221ebdf65a2SToby Isaac ierr = PetscFree(forest->adaptLabel);CHKERRQ(ierr); 2229a81d013SToby Isaac ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr); 22356ba9f64SToby Isaac ierr = DMDestroy(&forest->base);CHKERRQ(ierr); 224ba936b91SToby Isaac ierr = DMDestroy(&forest->adapt);CHKERRQ(ierr); 22530f902e7SToby Isaac ierr = PetscFree(forest->topology);CHKERRQ(ierr); 22630f902e7SToby Isaac ierr = PetscFree(forest);CHKERRQ(ierr); 227db4d5e8cSToby Isaac PetscFunctionReturn(0); 228db4d5e8cSToby Isaac } 229db4d5e8cSToby Isaac 230db4d5e8cSToby Isaac #undef __FUNCT__ 231dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetTopology" 2329be51f97SToby Isaac /*@C 2339be51f97SToby Isaac DMForestSetTopology - Set the topology of a DMForest during the pre-setup phase. The topology is a string (e.g. 2349be51f97SToby Isaac "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest durint 2359be51f97SToby Isaac DMSetUp(). 2369be51f97SToby Isaac 2379be51f97SToby Isaac Logically collective on dm 2389be51f97SToby Isaac 2399be51f97SToby Isaac Input parameters: 2409be51f97SToby Isaac + dm - the forest 2419be51f97SToby Isaac - topology - the topology of the forest 2429be51f97SToby Isaac 2439be51f97SToby Isaac Level: intermediate 2449be51f97SToby Isaac 2459be51f97SToby Isaac .seealso(): DMForestGetTopology(), DMForestSetBaseDM() 2469be51f97SToby Isaac @*/ 247dd8e54a2SToby Isaac PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology) 248db4d5e8cSToby Isaac { 249db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 250db4d5e8cSToby Isaac PetscErrorCode ierr; 251db4d5e8cSToby Isaac 252db4d5e8cSToby Isaac PetscFunctionBegin; 253db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 254ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the topology after setup"); 255dd8e54a2SToby Isaac ierr = PetscFree(forest->topology);CHKERRQ(ierr); 256dd8e54a2SToby Isaac ierr = PetscStrallocpy((const char*)topology,(char**) &forest->topology);CHKERRQ(ierr); 257db4d5e8cSToby Isaac PetscFunctionReturn(0); 258db4d5e8cSToby Isaac } 259db4d5e8cSToby Isaac 260db4d5e8cSToby Isaac #undef __FUNCT__ 261dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetTopology" 2629be51f97SToby Isaac /*@C 2639be51f97SToby Isaac DMForestGetTopology - Get a string describing the topology of a DMForest. 2649be51f97SToby Isaac 2659be51f97SToby Isaac Not collective 2669be51f97SToby Isaac 2679be51f97SToby Isaac Input parameter: 2689be51f97SToby Isaac . dm - the forest 2699be51f97SToby Isaac 2709be51f97SToby Isaac Output parameter: 2719be51f97SToby Isaac . topology - the topology of the forest (e.g., 'cube', 'shell') 2729be51f97SToby Isaac 2739be51f97SToby Isaac Level: intermediate 2749be51f97SToby Isaac 2759be51f97SToby Isaac .seealso: DMForestSetTopology() 2769be51f97SToby Isaac @*/ 277dd8e54a2SToby Isaac PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology) 278dd8e54a2SToby Isaac { 279dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 280dd8e54a2SToby Isaac 281dd8e54a2SToby Isaac PetscFunctionBegin; 282dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 283dd8e54a2SToby Isaac PetscValidPointer(topology,2); 284dd8e54a2SToby Isaac *topology = forest->topology; 285dd8e54a2SToby Isaac PetscFunctionReturn(0); 286dd8e54a2SToby Isaac } 287dd8e54a2SToby Isaac 288dd8e54a2SToby Isaac #undef __FUNCT__ 289dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetBaseDM" 2909be51f97SToby Isaac /*@ 2919be51f97SToby Isaac DMForestSetBaseDM - During the pre-setup phase, set the DM that defines the base mesh of a DMForest forest. The 2929be51f97SToby Isaac forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its 2939be51f97SToby Isaac base. In general, two forest must share a bse to be comparable, to do things like construct interpolators. 2949be51f97SToby Isaac 2959be51f97SToby Isaac Logically collective on dm 2969be51f97SToby Isaac 2979be51f97SToby Isaac Input Parameters: 2989be51f97SToby Isaac + dm - the forest 2999be51f97SToby Isaac - base - the base DM of the forest 3009be51f97SToby Isaac 3019be51f97SToby Isaac Level: intermediate 3029be51f97SToby Isaac 3039be51f97SToby Isaac .seealso(): DMForestGetBaseDM() 3049be51f97SToby Isaac @*/ 305dd8e54a2SToby Isaac PetscErrorCode DMForestSetBaseDM(DM dm, DM base) 306dd8e54a2SToby Isaac { 307dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 308dd8e54a2SToby Isaac PetscInt dim, dimEmbed; 309dd8e54a2SToby Isaac PetscErrorCode ierr; 310dd8e54a2SToby Isaac 311dd8e54a2SToby Isaac PetscFunctionBegin; 312dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 313ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the base after setup"); 314dd8e54a2SToby Isaac ierr = PetscObjectReference((PetscObject)base);CHKERRQ(ierr); 315dd8e54a2SToby Isaac ierr = DMDestroy(&forest->base);CHKERRQ(ierr); 316dd8e54a2SToby Isaac forest->base = base; 317a0452a8eSToby Isaac if (base) { 318a0452a8eSToby Isaac PetscValidHeaderSpecific(base, DM_CLASSID, 2); 319dd8e54a2SToby Isaac ierr = DMGetDimension(base,&dim);CHKERRQ(ierr); 320dd8e54a2SToby Isaac ierr = DMSetDimension(dm,dim);CHKERRQ(ierr); 321dd8e54a2SToby Isaac ierr = DMGetCoordinateDim(base,&dimEmbed);CHKERRQ(ierr); 322dd8e54a2SToby Isaac ierr = DMSetCoordinateDim(dm,dimEmbed);CHKERRQ(ierr); 323a0452a8eSToby Isaac } 324dd8e54a2SToby Isaac PetscFunctionReturn(0); 325dd8e54a2SToby Isaac } 326dd8e54a2SToby Isaac 327dd8e54a2SToby Isaac #undef __FUNCT__ 328dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetBaseDM" 3299be51f97SToby Isaac /*@ 3309be51f97SToby Isaac DMForestGetBaseDM - Get the base DM of a DMForest forest. The forest will be hierarchically refined from the base, 3319be51f97SToby Isaac and all refinements/coarsenings of the forest will share its base. In general, two forest must share a bse to be 3329be51f97SToby Isaac comparable, to do things like construct interpolators. 3339be51f97SToby Isaac 3349be51f97SToby Isaac Not collective 3359be51f97SToby Isaac 3369be51f97SToby Isaac Input Parameter: 3379be51f97SToby Isaac . dm - the forest 3389be51f97SToby Isaac 3399be51f97SToby Isaac Output Parameter: 3409be51f97SToby Isaac . base - the base DM of the forest 3419be51f97SToby Isaac 3429be51f97SToby Isaac Level: intermediate 3439be51f97SToby Isaac 3449be51f97SToby Isaac .seealso(); DMForestSetBaseDM() 3459be51f97SToby Isaac @*/ 346dd8e54a2SToby Isaac PetscErrorCode DMForestGetBaseDM(DM dm, DM *base) 347dd8e54a2SToby Isaac { 348dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 349dd8e54a2SToby Isaac 350dd8e54a2SToby Isaac PetscFunctionBegin; 351dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 352dd8e54a2SToby Isaac PetscValidPointer(base, 2); 353dd8e54a2SToby Isaac *base = forest->base; 354dd8e54a2SToby Isaac PetscFunctionReturn(0); 355dd8e54a2SToby Isaac } 356dd8e54a2SToby Isaac 357dd8e54a2SToby Isaac #undef __FUNCT__ 358cf38a08cSToby Isaac #define __FUNCT__ "DMForestSetBaseCoordinateMapping" 35999478f86SToby Isaac PetscErrorCode DMForestSetBaseCoordinateMapping(DM dm, PetscErrorCode (*func)(DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx) 360cf38a08cSToby Isaac { 361cf38a08cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 362cf38a08cSToby Isaac 363cf38a08cSToby Isaac PetscFunctionBegin; 364cf38a08cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 365cf38a08cSToby Isaac forest->mapcoordinates = func; 366cf38a08cSToby Isaac forest->mapcoordinatesctx = ctx; 367cf38a08cSToby Isaac PetscFunctionReturn(0); 368cf38a08cSToby Isaac } 369cf38a08cSToby Isaac 370cf38a08cSToby Isaac #undef __FUNCT__ 371cf38a08cSToby Isaac #define __FUNCT__ "DMForestGetBaseCoordinateMapping" 37299478f86SToby Isaac PetscErrorCode DMForestGetBaseCoordinateMapping(DM dm, PetscErrorCode (**func) (DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx) 373cf38a08cSToby Isaac { 374cf38a08cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 375cf38a08cSToby Isaac 376cf38a08cSToby Isaac PetscFunctionBegin; 377cf38a08cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 378cf38a08cSToby Isaac if (func) *func = forest->mapcoordinates; 379cf38a08cSToby Isaac if (ctx) *((void**) ctx) = forest->mapcoordinatesctx; 380cf38a08cSToby Isaac PetscFunctionReturn(0); 381cf38a08cSToby Isaac } 382cf38a08cSToby Isaac 383cf38a08cSToby Isaac #undef __FUNCT__ 384ba936b91SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityForest" 3859be51f97SToby Isaac /*@ 3869be51f97SToby Isaac DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be 3879be51f97SToby Isaac adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) im DMSetUp(). Usually not needed 3889be51f97SToby Isaac by users directly: DMForestTemplate() constructs a new forest to be adapted from an old forest and calls this 3899be51f97SToby Isaac routine. 3909be51f97SToby Isaac 391dffe73a3SToby Isaac Note that this can be called after setup with adapt = NULL, which will clear all internal data related to the 392dffe73a3SToby Isaac adaptivity forest from dm. This way, repeatedly adapting does not leave stale DM objects in memory. 393dffe73a3SToby Isaac 3949be51f97SToby Isaac Logically collective on dm 3959be51f97SToby Isaac 3969be51f97SToby Isaac Input Parameter: 3979be51f97SToby Isaac + dm - the new forest, which will be constructed from adapt 3989be51f97SToby Isaac - adapt - the old forest 3999be51f97SToby Isaac 4009be51f97SToby Isaac Level: intermediate 4019be51f97SToby Isaac 4029be51f97SToby Isaac .seealso: DMForestGetAdaptivityForest(), DMForestSetAdaptivityPurpose() 4039be51f97SToby Isaac @*/ 404ba936b91SToby Isaac PetscErrorCode DMForestSetAdaptivityForest(DM dm,DM adapt) 405dd8e54a2SToby Isaac { 406dffe73a3SToby Isaac DM_Forest *forest, *adaptForest, *oldAdaptForest; 407dffe73a3SToby Isaac DM oldAdapt; 408dd8e54a2SToby Isaac PetscErrorCode ierr; 409dd8e54a2SToby Isaac 410dd8e54a2SToby Isaac PetscFunctionBegin; 411dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 412ba936b91SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 413ba936b91SToby Isaac forest = (DM_Forest*) dm->data; 414dffe73a3SToby Isaac ierr = DMForestGetAdaptivityForest(dm,&oldAdapt);CHKERRQ(ierr); 415dffe73a3SToby Isaac if (adapt != NULL && dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup"); 416dffe73a3SToby Isaac adaptForest = (DM_Forest*) adapt ? adapt->data : NULL; 417dffe73a3SToby Isaac oldAdaptForest = (DM_Forest*) oldAdapt ? oldAdapt->data : NULL; 418dffe73a3SToby Isaac if (adaptForest != oldAdaptForest) { 419dffe73a3SToby Isaac ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr); 420dffe73a3SToby Isaac ierr = PetscSFDestroy(&forest->coarseToPreFine);CHKERRQ(ierr); 421dffe73a3SToby Isaac if (forest->clearadaptivityforest) {ierr = (forest->clearadaptivityforest)(dm);CHKERRQ(ierr);} 422dffe73a3SToby Isaac } 42326d9498aSToby Isaac switch (forest->adaptPurpose) { 42426d9498aSToby Isaac case DM_FOREST_KEEP: 425ba936b91SToby Isaac ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr); 426ba936b91SToby Isaac ierr = DMDestroy(&(forest->adapt));CHKERRQ(ierr); 427ba936b91SToby Isaac forest->adapt = adapt; 42826d9498aSToby Isaac break; 42926d9498aSToby Isaac case DM_FOREST_REFINE: 43026d9498aSToby Isaac ierr = DMSetCoarseDM(dm,adapt);CHKERRQ(ierr); 43126d9498aSToby Isaac break; 43226d9498aSToby Isaac case DM_FOREST_COARSEN: 43326d9498aSToby Isaac ierr = DMSetFineDM(dm,adapt);CHKERRQ(ierr); 43426d9498aSToby Isaac break; 43526d9498aSToby Isaac default: 43626d9498aSToby Isaac SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose"); 43726d9498aSToby Isaac } 438dd8e54a2SToby Isaac PetscFunctionReturn(0); 439dd8e54a2SToby Isaac } 440dd8e54a2SToby Isaac 441dd8e54a2SToby Isaac #undef __FUNCT__ 442ba936b91SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityForest" 4439be51f97SToby Isaac /*@ 4449be51f97SToby Isaac DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted. 4459be51f97SToby Isaac 4469be51f97SToby Isaac Not collective 4479be51f97SToby Isaac 4489be51f97SToby Isaac Input Parameter: 4499be51f97SToby Isaac . dm - the forest 4509be51f97SToby Isaac 4519be51f97SToby Isaac Output Parameter: 4529be51f97SToby Isaac . adapt - the forest from which dm is/was adapted 4539be51f97SToby Isaac 4549be51f97SToby Isaac Level: intermediate 4559be51f97SToby Isaac 4569be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose() 4579be51f97SToby Isaac @*/ 458ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt) 459dd8e54a2SToby Isaac { 460ba936b91SToby Isaac DM_Forest *forest; 46126d9498aSToby Isaac PetscErrorCode ierr; 462dd8e54a2SToby Isaac 463dd8e54a2SToby Isaac PetscFunctionBegin; 464dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 465ba936b91SToby Isaac forest = (DM_Forest*) dm->data; 46626d9498aSToby Isaac switch (forest->adaptPurpose) { 46726d9498aSToby Isaac case DM_FOREST_KEEP: 468ba936b91SToby Isaac *adapt = forest->adapt; 46926d9498aSToby Isaac break; 47026d9498aSToby Isaac case DM_FOREST_REFINE: 47126d9498aSToby Isaac ierr = DMGetCoarseDM(dm,adapt);CHKERRQ(ierr); 47226d9498aSToby Isaac break; 47326d9498aSToby Isaac case DM_FOREST_COARSEN: 47426d9498aSToby Isaac ierr = DMGetFineDM(dm,adapt);CHKERRQ(ierr); 47526d9498aSToby Isaac break; 47626d9498aSToby Isaac default: 47726d9498aSToby Isaac SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose"); 47826d9498aSToby Isaac } 47926d9498aSToby Isaac PetscFunctionReturn(0); 48026d9498aSToby Isaac } 48126d9498aSToby Isaac 48226d9498aSToby Isaac #undef __FUNCT__ 48326d9498aSToby Isaac #define __FUNCT__ "DMForestSetAdaptivityPurpose" 4849be51f97SToby Isaac /*@ 4859be51f97SToby Isaac DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its 4869be51f97SToby Isaac source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_FOREST_REFINE), coarsening 4879be51f97SToby Isaac (DM_FOREST_COARSEN), or undefined (DM_FOREST_NONE). This only matters for the purposes of reference counting: 4889be51f97SToby Isaac during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse 4899be51f97SToby Isaac relationship (see DMSetFineDM()/DMSetCoarseDM()). If the purpose is not refinement or coarsening, and the user does 4909be51f97SToby Isaac not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can 4919be51f97SToby Isaac cause a memory leak. This method is used by subtypes of DMForest when automatically constructing mesh hierarchies. 4929be51f97SToby Isaac 4939be51f97SToby Isaac Logically collective on dm 4949be51f97SToby Isaac 4959be51f97SToby Isaac Input Parameters: 4969be51f97SToby Isaac + dm - the forest 4979be51f97SToby Isaac - purpose - the adaptivity purpose (DM_FOREST_NONE/DM_FOREST_REFINE/DM_FOREST_COARSEN) 4989be51f97SToby Isaac 4999be51f97SToby Isaac Level: advanced 5009be51f97SToby Isaac 5019be51f97SToby Isaac .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest() 5029be51f97SToby Isaac @*/ 50326d9498aSToby Isaac PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMForestAdaptivityPurpose purpose) 50426d9498aSToby Isaac { 50526d9498aSToby Isaac DM_Forest *forest; 50626d9498aSToby Isaac PetscErrorCode ierr; 50726d9498aSToby Isaac 50826d9498aSToby Isaac PetscFunctionBegin; 50926d9498aSToby Isaac forest = (DM_Forest*) dm->data; 5109be51f97SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup"); 51126d9498aSToby Isaac if (purpose != forest->adaptPurpose) { 51226d9498aSToby Isaac DM adapt; 51326d9498aSToby Isaac 51426d9498aSToby Isaac ierr = DMForestGetAdaptivityForest(dm,&adapt);CHKERRQ(ierr); 51526d9498aSToby Isaac ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr); 51626d9498aSToby Isaac ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr); 517f885a11aSToby Isaac 51826d9498aSToby Isaac forest->adaptPurpose = purpose; 519f885a11aSToby Isaac 52026d9498aSToby Isaac ierr = DMForestSetAdaptivityForest(dm,adapt);CHKERRQ(ierr); 52126d9498aSToby Isaac ierr = DMDestroy(&adapt);CHKERRQ(ierr); 52226d9498aSToby Isaac } 523dd8e54a2SToby Isaac PetscFunctionReturn(0); 524dd8e54a2SToby Isaac } 525dd8e54a2SToby Isaac 526dd8e54a2SToby Isaac #undef __FUNCT__ 52756c0450aSToby Isaac #define __FUNCT__ "DMForestGetAdaptivityPurpose" 52856c0450aSToby Isaac /*@ 52956c0450aSToby Isaac DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with 53056c0450aSToby Isaac DMForestSetAdaptivityForest()) for the purpose of refinement (DM_FOREST_REFINE), coarsening (DM_FOREST_COARSEN), or 53156c0450aSToby Isaac undefined (DM_FOREST_NONE). This only matters for the purposes of reference counting: during DMDestroy(), cyclic 53256c0450aSToby Isaac references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see 53356c0450aSToby Isaac DMSetFineDM()/DMSetCoarseDM()). If the purpose is not refinement or coarsening, and the user does not maintain a 53456c0450aSToby Isaac reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory 53556c0450aSToby Isaac leak. This method is used by subtypes of DMForest when automatically constructing mesh hierarchies. 53656c0450aSToby Isaac 53756c0450aSToby Isaac Not collective 53856c0450aSToby Isaac 53956c0450aSToby Isaac Input Parameter: 54056c0450aSToby Isaac . dm - the forest 54156c0450aSToby Isaac 54256c0450aSToby Isaac Output Parameter: 54356c0450aSToby Isaac . purpose - the adaptivity purpose (DM_FOREST_NONE/DM_FOREST_REFINE/DM_FOREST_COARSEN) 54456c0450aSToby Isaac 54556c0450aSToby Isaac Level: advanced 54656c0450aSToby Isaac 54756c0450aSToby Isaac .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest() 54856c0450aSToby Isaac @*/ 54956c0450aSToby Isaac PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMForestAdaptivityPurpose *purpose) 55056c0450aSToby Isaac { 55156c0450aSToby Isaac DM_Forest *forest; 55256c0450aSToby Isaac 55356c0450aSToby Isaac PetscFunctionBegin; 55456c0450aSToby Isaac forest = (DM_Forest*) dm->data; 55556c0450aSToby Isaac *purpose = forest->adaptPurpose; 55656c0450aSToby Isaac PetscFunctionReturn(0); 55756c0450aSToby Isaac } 55856c0450aSToby Isaac 55956c0450aSToby Isaac #undef __FUNCT__ 560dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetAdjacencyDimension" 5619be51f97SToby Isaac /*@ 5629be51f97SToby Isaac DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine 5639be51f97SToby Isaac cell adjacency (for the purposes of partitioning and overlap). 5649be51f97SToby Isaac 5659be51f97SToby Isaac Logically collective on dm 5669be51f97SToby Isaac 5679be51f97SToby Isaac Input Parameters: 5689be51f97SToby Isaac + dm - the forest 5699be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency) 5709be51f97SToby Isaac 5719be51f97SToby Isaac Level: intermediate 5729be51f97SToby Isaac 5739be51f97SToby Isaac .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap() 5749be51f97SToby Isaac @*/ 575dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim) 576dd8e54a2SToby Isaac { 577dd8e54a2SToby Isaac PetscInt dim; 578dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 579dd8e54a2SToby Isaac PetscErrorCode ierr; 580dd8e54a2SToby Isaac 581dd8e54a2SToby Isaac PetscFunctionBegin; 582dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 583ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup"); 584dd8e54a2SToby Isaac if (adjDim < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %d", adjDim); 585dd8e54a2SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 586dd8e54a2SToby Isaac if (adjDim > dim) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %d: %d", dim, adjDim); 587dd8e54a2SToby Isaac forest->adjDim = adjDim; 588dd8e54a2SToby Isaac PetscFunctionReturn(0); 589dd8e54a2SToby Isaac } 590dd8e54a2SToby Isaac 591dd8e54a2SToby Isaac #undef __FUNCT__ 592dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetAdjacencyCodimension" 5939be51f97SToby Isaac /*@ 5949be51f97SToby Isaac DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that, 5959be51f97SToby Isaac e.g., adjacency based on facets can be specified by codimension 1 in all cases) 5969be51f97SToby Isaac 5979be51f97SToby Isaac Logically collective on dm 5989be51f97SToby Isaac 5999be51f97SToby Isaac Input Parameters: 6009be51f97SToby Isaac + dm - the forest 6019be51f97SToby Isaac - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices 6029be51f97SToby Isaac 6039be51f97SToby Isaac Level: intermediate 6049be51f97SToby Isaac 6059be51f97SToby Isaac .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension() 6069be51f97SToby Isaac @*/ 607dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim) 608dd8e54a2SToby Isaac { 609dd8e54a2SToby Isaac PetscInt dim; 610dd8e54a2SToby Isaac PetscErrorCode ierr; 611dd8e54a2SToby Isaac 612dd8e54a2SToby Isaac PetscFunctionBegin; 613dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 614dd8e54a2SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 615dd8e54a2SToby Isaac ierr = DMForestSetAdjacencyDimension(dm,dim-adjCodim);CHKERRQ(ierr); 616dd8e54a2SToby Isaac PetscFunctionReturn(0); 617dd8e54a2SToby Isaac } 618dd8e54a2SToby Isaac 619dd8e54a2SToby Isaac #undef __FUNCT__ 620dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetAdjacencyDimension" 6219be51f97SToby Isaac /*@ 6229be51f97SToby Isaac DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the 6239be51f97SToby Isaac purposes of partitioning and overlap). 6249be51f97SToby Isaac 6259be51f97SToby Isaac Not collective 6269be51f97SToby Isaac 6279be51f97SToby Isaac Input Parameter: 6289be51f97SToby Isaac . dm - the forest 6299be51f97SToby Isaac 6309be51f97SToby Isaac Output Parameter: 6319be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency) 6329be51f97SToby Isaac 6339be51f97SToby Isaac Level: intermediate 6349be51f97SToby Isaac 6359be51f97SToby Isaac .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap() 6369be51f97SToby Isaac @*/ 637dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim) 638dd8e54a2SToby Isaac { 639dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 640dd8e54a2SToby Isaac 641dd8e54a2SToby Isaac PetscFunctionBegin; 642dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 643dd8e54a2SToby Isaac PetscValidIntPointer(adjDim,2); 644dd8e54a2SToby Isaac *adjDim = forest->adjDim; 645dd8e54a2SToby Isaac PetscFunctionReturn(0); 646dd8e54a2SToby Isaac } 647dd8e54a2SToby Isaac 648dd8e54a2SToby Isaac #undef __FUNCT__ 649dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetAdjacencyCodimension" 6509be51f97SToby Isaac /*@ 6519be51f97SToby Isaac DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that, 6529be51f97SToby Isaac e.g., adjacency based on facets can be specified by codimension 1 in all cases) 6539be51f97SToby Isaac 6549be51f97SToby Isaac Not collective 6559be51f97SToby Isaac 6569be51f97SToby Isaac Input Parameter: 6579be51f97SToby Isaac . dm - the forest 6589be51f97SToby Isaac 6599be51f97SToby Isaac Output Parameter: 6609be51f97SToby Isaac . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices 6619be51f97SToby Isaac 6629be51f97SToby Isaac Level: intermediate 6639be51f97SToby Isaac 6649be51f97SToby Isaac .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension() 6659be51f97SToby Isaac @*/ 666dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim) 667dd8e54a2SToby Isaac { 668dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 669dd8e54a2SToby Isaac PetscInt dim; 670dd8e54a2SToby Isaac PetscErrorCode ierr; 671dd8e54a2SToby Isaac 672dd8e54a2SToby Isaac PetscFunctionBegin; 673dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 674dd8e54a2SToby Isaac PetscValidIntPointer(adjCodim,2); 675dd8e54a2SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 676dd8e54a2SToby Isaac *adjCodim = dim - forest->adjDim; 677dd8e54a2SToby Isaac PetscFunctionReturn(0); 678dd8e54a2SToby Isaac } 679dd8e54a2SToby Isaac 680dd8e54a2SToby Isaac #undef __FUNCT__ 681ef51cf95SToby Isaac #define __FUNCT__ "DMForestSetPartitionOverlap" 6829be51f97SToby Isaac /*@ 6839be51f97SToby Isaac DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel 6849be51f97SToby Isaac partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding 6859be51f97SToby Isaac adjacent cells 6869be51f97SToby Isaac 6879be51f97SToby Isaac Logically collective on dm 6889be51f97SToby Isaac 6899be51f97SToby Isaac Input Parameters: 6909be51f97SToby Isaac + dm - the forest 6919be51f97SToby Isaac - overlap - default 0 6929be51f97SToby Isaac 6939be51f97SToby Isaac Level: intermediate 6949be51f97SToby Isaac 6959be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension() 6969be51f97SToby Isaac @*/ 697dd8e54a2SToby Isaac PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap) 698dd8e54a2SToby Isaac { 699dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 700dd8e54a2SToby Isaac 701dd8e54a2SToby Isaac PetscFunctionBegin; 702dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 703ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup"); 704dd8e54a2SToby Isaac if (overlap < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %d", overlap); 705dd8e54a2SToby Isaac forest->overlap = overlap; 706dd8e54a2SToby Isaac PetscFunctionReturn(0); 707dd8e54a2SToby Isaac } 708dd8e54a2SToby Isaac 709dd8e54a2SToby Isaac #undef __FUNCT__ 710dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetPartitionOverlap" 7119be51f97SToby Isaac /*@ 7129be51f97SToby Isaac DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values 7139be51f97SToby Isaac > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells 7149be51f97SToby Isaac 7159be51f97SToby Isaac Not collective 7169be51f97SToby Isaac 7179be51f97SToby Isaac Input Parameter: 7189be51f97SToby Isaac . dm - the forest 7199be51f97SToby Isaac 7209be51f97SToby Isaac Output Parameter: 7219be51f97SToby Isaac . overlap - default 0 7229be51f97SToby Isaac 7239be51f97SToby Isaac Level: intermediate 7249be51f97SToby Isaac 7259be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension() 7269be51f97SToby Isaac @*/ 727dd8e54a2SToby Isaac PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap) 728dd8e54a2SToby Isaac { 729dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 730dd8e54a2SToby Isaac 731dd8e54a2SToby Isaac PetscFunctionBegin; 732dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 733dd8e54a2SToby Isaac PetscValidIntPointer(overlap,2); 734dd8e54a2SToby Isaac *overlap = forest->overlap; 735dd8e54a2SToby Isaac PetscFunctionReturn(0); 736dd8e54a2SToby Isaac } 737dd8e54a2SToby Isaac 738dd8e54a2SToby Isaac #undef __FUNCT__ 739dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetMinimumRefinement" 7409be51f97SToby Isaac /*@ 7419be51f97SToby Isaac DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base 7429be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. If the forest is being created by coarsening a previous forest 7439be51f97SToby Isaac (see DMForestGetAdaptivityForest()) this limits the amount of coarsening. 7449be51f97SToby Isaac 7459be51f97SToby Isaac Logically collective on dm 7469be51f97SToby Isaac 7479be51f97SToby Isaac Input Parameters: 7489be51f97SToby Isaac + dm - the forest 7499be51f97SToby Isaac - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7509be51f97SToby Isaac 7519be51f97SToby Isaac Level: intermediate 7529be51f97SToby Isaac 7539be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 7549be51f97SToby Isaac @*/ 755dd8e54a2SToby Isaac PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement) 756dd8e54a2SToby Isaac { 757dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 758dd8e54a2SToby Isaac 759dd8e54a2SToby Isaac PetscFunctionBegin; 760dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 761ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup"); 762dd8e54a2SToby Isaac forest->minRefinement = minRefinement; 763dd8e54a2SToby Isaac PetscFunctionReturn(0); 764dd8e54a2SToby Isaac } 765dd8e54a2SToby Isaac 766dd8e54a2SToby Isaac #undef __FUNCT__ 767dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetMinimumRefinement" 7689be51f97SToby Isaac /*@ 7699be51f97SToby Isaac DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see 7709be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. If the forest is being created by coarsening a previous forest (see 7719be51f97SToby Isaac DMForestGetAdaptivityForest()), this limits the amount of coarsening. 7729be51f97SToby Isaac 7739be51f97SToby Isaac Not collective 7749be51f97SToby Isaac 7759be51f97SToby Isaac Input Parameter: 7769be51f97SToby Isaac . dm - the forest 7779be51f97SToby Isaac 7789be51f97SToby Isaac Output Parameter: 7799be51f97SToby Isaac . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7809be51f97SToby Isaac 7819be51f97SToby Isaac Level: intermediate 7829be51f97SToby Isaac 7839be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestGetMaximumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 7849be51f97SToby Isaac @*/ 785dd8e54a2SToby Isaac PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement) 786dd8e54a2SToby Isaac { 787dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 788dd8e54a2SToby Isaac 789dd8e54a2SToby Isaac PetscFunctionBegin; 790dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 791dd8e54a2SToby Isaac PetscValidIntPointer(minRefinement,2); 792dd8e54a2SToby Isaac *minRefinement = forest->minRefinement; 793dd8e54a2SToby Isaac PetscFunctionReturn(0); 794dd8e54a2SToby Isaac } 795dd8e54a2SToby Isaac 796dd8e54a2SToby Isaac #undef __FUNCT__ 79756ba9f64SToby Isaac #define __FUNCT__ "DMForestSetInitialRefinement" 7989be51f97SToby Isaac /*@ 7999be51f97SToby Isaac DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base 8009be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. 8019be51f97SToby Isaac 8029be51f97SToby Isaac Logically collective on dm 8039be51f97SToby Isaac 8049be51f97SToby Isaac Input Parameters: 8059be51f97SToby Isaac + dm - the forest 8069be51f97SToby Isaac - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 8079be51f97SToby Isaac 8089be51f97SToby Isaac Level: intermediate 8099be51f97SToby Isaac 8109be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM() 8119be51f97SToby Isaac @*/ 81256ba9f64SToby Isaac PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement) 81356ba9f64SToby Isaac { 81456ba9f64SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 81556ba9f64SToby Isaac 81656ba9f64SToby Isaac PetscFunctionBegin; 81756ba9f64SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 81856ba9f64SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup"); 81956ba9f64SToby Isaac forest->initRefinement = initRefinement; 82056ba9f64SToby Isaac PetscFunctionReturn(0); 82156ba9f64SToby Isaac } 82256ba9f64SToby Isaac 82356ba9f64SToby Isaac #undef __FUNCT__ 82456ba9f64SToby Isaac #define __FUNCT__ "DMForestGetInitialRefinement" 8259be51f97SToby Isaac /*@ 8269be51f97SToby Isaac DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see 8279be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. 8289be51f97SToby Isaac 8299be51f97SToby Isaac Not collective 8309be51f97SToby Isaac 8319be51f97SToby Isaac Input Parameter: 8329be51f97SToby Isaac . dm - the forest 8339be51f97SToby Isaac 8349be51f97SToby Isaac Output Paramater: 8359be51f97SToby Isaac . initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 8369be51f97SToby Isaac 8379be51f97SToby Isaac Level: intermediate 8389be51f97SToby Isaac 8399be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM() 8409be51f97SToby Isaac @*/ 84156ba9f64SToby Isaac PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement) 84256ba9f64SToby Isaac { 84356ba9f64SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 84456ba9f64SToby Isaac 84556ba9f64SToby Isaac PetscFunctionBegin; 84656ba9f64SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 84756ba9f64SToby Isaac PetscValidIntPointer(initRefinement,2); 84856ba9f64SToby Isaac *initRefinement = forest->initRefinement; 84956ba9f64SToby Isaac PetscFunctionReturn(0); 85056ba9f64SToby Isaac } 85156ba9f64SToby Isaac 85256ba9f64SToby Isaac #undef __FUNCT__ 853c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetMaximumRefinement" 8549be51f97SToby Isaac /*@ 8559be51f97SToby Isaac DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base 8569be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. If the forest is being created by refining a previous forest 8579be51f97SToby Isaac (see DMForestGetAdaptivityForest()), this limits the amount of refinement. 8589be51f97SToby Isaac 8599be51f97SToby Isaac Logically collective on dm 8609be51f97SToby Isaac 8619be51f97SToby Isaac Input Parameters: 8629be51f97SToby Isaac + dm - the forest 8639be51f97SToby Isaac - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 8649be51f97SToby Isaac 8659be51f97SToby Isaac Level: intermediate 8669be51f97SToby Isaac 8679be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityDM() 8689be51f97SToby Isaac @*/ 869c7eeac06SToby Isaac PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement) 870dd8e54a2SToby Isaac { 871dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 872dd8e54a2SToby Isaac 873dd8e54a2SToby Isaac PetscFunctionBegin; 874dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 875ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup"); 876c7eeac06SToby Isaac forest->maxRefinement = maxRefinement; 877dd8e54a2SToby Isaac PetscFunctionReturn(0); 878dd8e54a2SToby Isaac } 879dd8e54a2SToby Isaac 880dd8e54a2SToby Isaac #undef __FUNCT__ 881c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetMaximumRefinement" 8829be51f97SToby Isaac /*@ 8839be51f97SToby Isaac DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see 8849be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. If the forest is being created by refining a previous forest (see 8859be51f97SToby Isaac DMForestGetAdaptivityForest()), this limits the amount of refinement. 8869be51f97SToby Isaac 8879be51f97SToby Isaac Not collective 8889be51f97SToby Isaac 8899be51f97SToby Isaac Input Parameter: 8909be51f97SToby Isaac . dm - the forest 8919be51f97SToby Isaac 8929be51f97SToby Isaac Output Parameter: 8939be51f97SToby Isaac . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 8949be51f97SToby Isaac 8959be51f97SToby Isaac Level: intermediate 8969be51f97SToby Isaac 8979be51f97SToby Isaac .seealso: DMForestSetMaximumRefinement(), DMForestGetMinimumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 8989be51f97SToby Isaac @*/ 899c7eeac06SToby Isaac PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement) 900dd8e54a2SToby Isaac { 901dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 902dd8e54a2SToby Isaac 903dd8e54a2SToby Isaac PetscFunctionBegin; 904dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 905c7eeac06SToby Isaac PetscValidIntPointer(maxRefinement,2); 906c7eeac06SToby Isaac *maxRefinement = forest->maxRefinement; 907dd8e54a2SToby Isaac PetscFunctionReturn(0); 908dd8e54a2SToby Isaac } 909c7eeac06SToby Isaac 910c7eeac06SToby Isaac #undef __FUNCT__ 911c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityStrategy" 9129be51f97SToby Isaac /*@C 9139be51f97SToby Isaac DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes. 9149be51f97SToby Isaac Subtypes of DMForest may define their own strategies. Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree 9159be51f97SToby Isaac for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to 9169be51f97SToby Isaac specify refinement/coarsening. 9179be51f97SToby Isaac 9189be51f97SToby Isaac Logically collective on dm 9199be51f97SToby Isaac 9209be51f97SToby Isaac Input Parameters: 9219be51f97SToby Isaac + dm - the forest 9229be51f97SToby Isaac - adaptStrategy - default DMFORESTADAPTALL 9239be51f97SToby Isaac 9249be51f97SToby Isaac Level: advanced 9259be51f97SToby Isaac 9269be51f97SToby Isaac .seealso: DMForestGetAdaptivityStrategy() 9279be51f97SToby Isaac @*/ 928c7eeac06SToby Isaac PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy) 929c7eeac06SToby Isaac { 930c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 931c7eeac06SToby Isaac PetscErrorCode ierr; 932c7eeac06SToby Isaac 933c7eeac06SToby Isaac PetscFunctionBegin; 934c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 935c7eeac06SToby Isaac ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr); 936a73e2921SToby Isaac ierr = PetscStrallocpy((const char*) adaptStrategy,(char**)&forest->adaptStrategy);CHKERRQ(ierr); 937c7eeac06SToby Isaac PetscFunctionReturn(0); 938c7eeac06SToby Isaac } 939c7eeac06SToby Isaac 940c7eeac06SToby Isaac #undef __FUNCT__ 941c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityStrategy" 9429be51f97SToby Isaac /*@C 9439be51f97SToby Isaac DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes. Subtypes 9449be51f97SToby Isaac of DMForest may define their own strategies. Two default strategies are DMFORESTADAPTALL, which indicates that all 9459be51f97SToby Isaac processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only 9469be51f97SToby Isaac one process needs to specify refinement/coarsening. 9479be51f97SToby Isaac 9489be51f97SToby Isaac Not collective 9499be51f97SToby Isaac 9509be51f97SToby Isaac Input Parameter: 9519be51f97SToby Isaac . dm - the forest 9529be51f97SToby Isaac 9539be51f97SToby Isaac Output Parameter: 9549be51f97SToby Isaac . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL) 9559be51f97SToby Isaac 9569be51f97SToby Isaac Level: advanced 9579be51f97SToby Isaac 9589be51f97SToby Isaac .seealso: DMForestSetAdaptivityStrategy() 9599be51f97SToby Isaac @*/ 960c7eeac06SToby Isaac PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy) 961c7eeac06SToby Isaac { 962c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 963c7eeac06SToby Isaac 964c7eeac06SToby Isaac PetscFunctionBegin; 965c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 966c7eeac06SToby Isaac PetscValidPointer(adaptStrategy,2); 967c7eeac06SToby Isaac *adaptStrategy = forest->adaptStrategy; 968c7eeac06SToby Isaac PetscFunctionReturn(0); 969c7eeac06SToby Isaac } 970c7eeac06SToby Isaac 971c7eeac06SToby Isaac #undef __FUNCT__ 972bf9b5d84SToby Isaac #define __FUNCT__ "DMForestSetComputeAdaptivitySF" 973bf9b5d84SToby Isaac /*@ 974bf9b5d84SToby Isaac DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed 975bf9b5d84SToby 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(). 976bf9b5d84SToby Isaac 977bf9b5d84SToby Isaac Logically collective on dm 978bf9b5d84SToby Isaac 979bf9b5d84SToby Isaac Input Parameters: 980bf9b5d84SToby Isaac + dm - the post-adaptation forest 981bf9b5d84SToby Isaac - computeSF - default PETSC_TRUE 982bf9b5d84SToby Isaac 983bf9b5d84SToby Isaac Level: advanced 984bf9b5d84SToby Isaac 985bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF() 986bf9b5d84SToby Isaac @*/ 987bf9b5d84SToby Isaac PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF) 988bf9b5d84SToby Isaac { 989bf9b5d84SToby Isaac DM_Forest *forest; 990bf9b5d84SToby Isaac 991bf9b5d84SToby Isaac PetscFunctionBegin; 992bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 993bf9b5d84SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called"); 994bf9b5d84SToby Isaac forest = (DM_Forest*) dm->data; 995bf9b5d84SToby Isaac forest->computeAdaptSF = computeSF; 996bf9b5d84SToby Isaac PetscFunctionReturn(0); 997bf9b5d84SToby Isaac } 998bf9b5d84SToby Isaac 999bf9b5d84SToby Isaac #undef __FUNCT__ 100080b27e07SToby Isaac #define __FUNCT__ "DMForestTransferVec" 1001*0eb7e1eaSToby Isaac PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time) 100280b27e07SToby Isaac { 100380b27e07SToby Isaac DM_Forest *forest; 100480b27e07SToby Isaac PetscErrorCode ierr; 100580b27e07SToby Isaac 100680b27e07SToby Isaac PetscFunctionBegin; 100780b27e07SToby Isaac PetscValidHeaderSpecific(dmIn ,DM_CLASSID ,1); 100880b27e07SToby Isaac PetscValidHeaderSpecific(vecIn ,VEC_CLASSID ,2); 100980b27e07SToby Isaac PetscValidHeaderSpecific(dmOut ,DM_CLASSID ,3); 101080b27e07SToby Isaac PetscValidHeaderSpecific(vecOut ,VEC_CLASSID ,4); 101180b27e07SToby Isaac forest = (DM_Forest *) dmIn->data; 101280b27e07SToby Isaac if (!forest->transfervec) SETERRQ(PetscObjectComm((PetscObject)dmIn),PETSC_ERR_SUP,"DMForestTransferVec() not implemented"); 1013*0eb7e1eaSToby Isaac ierr = (forest->transfervec)(dmIn,vecIn,dmOut,vecOut,useBCs,time);CHKERRQ(ierr); 101480b27e07SToby Isaac PetscFunctionReturn(0); 101580b27e07SToby Isaac } 101680b27e07SToby Isaac 101780b27e07SToby Isaac #undef __FUNCT__ 1018bf9b5d84SToby Isaac #define __FUNCT__ "DMForestGetComputeAdaptivitySF" 1019bf9b5d84SToby Isaac /*@ 1020bf9b5d84SToby Isaac DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the 1021bf9b5d84SToby Isaac pre-adaptation forest to the post-adaptiation forest. After DMSetUp() is called, these transfer PetscSFs can be 1022bf9b5d84SToby Isaac accessed with DMForestGetAdaptivitySF(). 1023bf9b5d84SToby Isaac 1024bf9b5d84SToby Isaac Not collective 1025bf9b5d84SToby Isaac 1026bf9b5d84SToby Isaac Input Parameter: 1027bf9b5d84SToby Isaac . dm - the post-adaptation forest 1028bf9b5d84SToby Isaac 1029bf9b5d84SToby Isaac Output Parameter: 1030bf9b5d84SToby Isaac . computeSF - default PETSC_TRUE 1031bf9b5d84SToby Isaac 1032bf9b5d84SToby Isaac Level: advanced 1033bf9b5d84SToby Isaac 1034bf9b5d84SToby Isaac .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF() 1035bf9b5d84SToby Isaac @*/ 1036bf9b5d84SToby Isaac PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF) 1037bf9b5d84SToby Isaac { 1038bf9b5d84SToby Isaac DM_Forest *forest; 1039bf9b5d84SToby Isaac 1040bf9b5d84SToby Isaac PetscFunctionBegin; 1041bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1042bf9b5d84SToby Isaac forest = (DM_Forest*) dm->data; 1043bf9b5d84SToby Isaac *computeSF = forest->computeAdaptSF; 1044bf9b5d84SToby Isaac PetscFunctionReturn(0); 1045bf9b5d84SToby Isaac } 1046bf9b5d84SToby Isaac 1047bf9b5d84SToby Isaac #undef __FUNCT__ 1048bf9b5d84SToby Isaac #define __FUNCT__ "DMForestGetAdaptivitySF" 1049bf9b5d84SToby Isaac /*@ 1050bf9b5d84SToby Isaac DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest. 1051bf9b5d84SToby Isaac Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be 1052bf9b5d84SToby Isaac some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa. Therefore there are two 1053bf9b5d84SToby Isaac PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates 1054bf9b5d84SToby Isaac pre-adaptation fine cells to post-adaptation coarse cells. 1055bf9b5d84SToby Isaac 1056bf9b5d84SToby Isaac Not collective 1057bf9b5d84SToby Isaac 1058bf9b5d84SToby Isaac Input Parameter: 1059bf9b5d84SToby Isaac dm - the post-adaptation forest 1060bf9b5d84SToby Isaac 1061bf9b5d84SToby Isaac Output Parameter: 10620f17b9e3SToby Isaac preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post- 10630f17b9e3SToby Isaac coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre- 1064bf9b5d84SToby Isaac 1065bf9b5d84SToby Isaac Level: advanced 1066bf9b5d84SToby Isaac 1067bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF() 1068bf9b5d84SToby Isaac @*/ 10690f17b9e3SToby Isaac PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine) 1070bf9b5d84SToby Isaac { 1071bf9b5d84SToby Isaac DM_Forest *forest; 1072bf9b5d84SToby Isaac PetscErrorCode ierr; 1073bf9b5d84SToby Isaac 1074bf9b5d84SToby Isaac PetscFunctionBegin; 1075bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1076bf9b5d84SToby Isaac ierr = DMSetUp(dm);CHKERRQ(ierr); 1077bf9b5d84SToby Isaac forest = (DM_Forest*) dm->data; 1078f885a11aSToby Isaac if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine; 1079f885a11aSToby Isaac if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine; 1080bf9b5d84SToby Isaac PetscFunctionReturn(0); 1081bf9b5d84SToby Isaac } 1082bf9b5d84SToby Isaac 1083bf9b5d84SToby Isaac #undef __FUNCT__ 1084c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetGradeFactor" 10859be51f97SToby Isaac /*@ 10869be51f97SToby Isaac DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to 10879be51f97SToby Isaac indicate that the diameter of neighboring cells should differ by at most a factor of 2. Subtypes of DMForest may 10889be51f97SToby Isaac only support one particular choice of grading factor. 10899be51f97SToby Isaac 10909be51f97SToby Isaac Logically collective on dm 10919be51f97SToby Isaac 10929be51f97SToby Isaac Input Parameters: 10939be51f97SToby Isaac + dm - the forest 10949be51f97SToby Isaac - grade - the grading factor 10959be51f97SToby Isaac 10969be51f97SToby Isaac Level: advanced 10979be51f97SToby Isaac 10989be51f97SToby Isaac .seealso: DMForestGetGradeFactor() 10999be51f97SToby Isaac @*/ 1100c7eeac06SToby Isaac PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade) 1101c7eeac06SToby Isaac { 1102c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1103c7eeac06SToby Isaac 1104c7eeac06SToby Isaac PetscFunctionBegin; 1105c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1106ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup"); 1107c7eeac06SToby Isaac forest->gradeFactor = grade; 1108c7eeac06SToby Isaac PetscFunctionReturn(0); 1109c7eeac06SToby Isaac } 1110c7eeac06SToby Isaac 1111c7eeac06SToby Isaac #undef __FUNCT__ 1112c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetGradeFactor" 11139be51f97SToby Isaac /*@ 11149be51f97SToby Isaac DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of 11159be51f97SToby Isaac neighboring cells should differ by at most a factor of 2. Subtypes of DMForest may only support one particular 11169be51f97SToby Isaac choice of grading factor. 11179be51f97SToby Isaac 11189be51f97SToby Isaac Not collective 11199be51f97SToby Isaac 11209be51f97SToby Isaac Input Parameter: 11219be51f97SToby Isaac . dm - the forest 11229be51f97SToby Isaac 11239be51f97SToby Isaac Output Parameter: 11249be51f97SToby Isaac . grade - the grading factor 11259be51f97SToby Isaac 11269be51f97SToby Isaac Level: advanced 11279be51f97SToby Isaac 11289be51f97SToby Isaac .seealso: DMForestSetGradeFactor() 11299be51f97SToby Isaac @*/ 1130c7eeac06SToby Isaac PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade) 1131c7eeac06SToby Isaac { 1132c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1133c7eeac06SToby Isaac 1134c7eeac06SToby Isaac PetscFunctionBegin; 1135c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1136c7eeac06SToby Isaac PetscValidIntPointer(grade,2); 1137c7eeac06SToby Isaac *grade = forest->gradeFactor; 1138c7eeac06SToby Isaac PetscFunctionReturn(0); 1139c7eeac06SToby Isaac } 1140c7eeac06SToby Isaac 1141c7eeac06SToby Isaac #undef __FUNCT__ 1142ef51cf95SToby Isaac #define __FUNCT__ "DMForestSetCellWeightFactor" 11439be51f97SToby Isaac /*@ 11449be51f97SToby Isaac DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes 11459be51f97SToby Isaac the cell weight (see DMForestSetCellWeights()) when calculating partitions. The final weight of a cell will be 11469be51f97SToby Isaac (cellWeight) * (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on 11479be51f97SToby Isaac its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the 11489be51f97SToby Isaac computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement. 11499be51f97SToby Isaac 11509be51f97SToby Isaac Logically collective on dm 11519be51f97SToby Isaac 11529be51f97SToby Isaac Input Parameters: 11539be51f97SToby Isaac + dm - the forest 11549be51f97SToby Isaac - weightsFactors - default 1. 11559be51f97SToby Isaac 11569be51f97SToby Isaac Level: advanced 11579be51f97SToby Isaac 11589be51f97SToby Isaac .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights() 11599be51f97SToby Isaac @*/ 1160ef51cf95SToby Isaac PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor) 1161c7eeac06SToby Isaac { 1162c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1163c7eeac06SToby Isaac 1164c7eeac06SToby Isaac PetscFunctionBegin; 1165c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1166ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup"); 1167c7eeac06SToby Isaac forest->weightsFactor = weightsFactor; 1168c7eeac06SToby Isaac PetscFunctionReturn(0); 1169c7eeac06SToby Isaac } 1170c7eeac06SToby Isaac 1171c7eeac06SToby Isaac #undef __FUNCT__ 1172ef51cf95SToby Isaac #define __FUNCT__ "DMForestGetCellWeightFactor" 11739be51f97SToby Isaac /*@ 11749be51f97SToby Isaac DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see 11759be51f97SToby Isaac DMForestSetCellWeights()) when calculating partitions. The final weight of a cell will be (cellWeight) * 11769be51f97SToby Isaac (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on its level; a 11779be51f97SToby Isaac factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation 11789be51f97SToby Isaac associated with a cell is multiplied by a factor of 2 for each additional level of refinement. 11799be51f97SToby Isaac 11809be51f97SToby Isaac Not collective 11819be51f97SToby Isaac 11829be51f97SToby Isaac Input Parameter: 11839be51f97SToby Isaac . dm - the forest 11849be51f97SToby Isaac 11859be51f97SToby Isaac Output Parameter: 11869be51f97SToby Isaac . weightsFactors - default 1. 11879be51f97SToby Isaac 11889be51f97SToby Isaac Level: advanced 11899be51f97SToby Isaac 11909be51f97SToby Isaac .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights() 11919be51f97SToby Isaac @*/ 1192ef51cf95SToby Isaac PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor) 1193c7eeac06SToby Isaac { 1194c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1195c7eeac06SToby Isaac 1196c7eeac06SToby Isaac PetscFunctionBegin; 1197c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1198c7eeac06SToby Isaac PetscValidRealPointer(weightsFactor,2); 1199c7eeac06SToby Isaac *weightsFactor = forest->weightsFactor; 1200c7eeac06SToby Isaac PetscFunctionReturn(0); 1201c7eeac06SToby Isaac } 1202c7eeac06SToby Isaac 1203c7eeac06SToby Isaac #undef __FUNCT__ 1204c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellChart" 12059be51f97SToby Isaac /*@ 12069be51f97SToby Isaac DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process 12079be51f97SToby Isaac 12089be51f97SToby Isaac Not collective 12099be51f97SToby Isaac 12109be51f97SToby Isaac Input Parameter: 12119be51f97SToby Isaac . dm - the forest 12129be51f97SToby Isaac 12139be51f97SToby Isaac Output Parameters: 12149be51f97SToby Isaac + cStart - the first cell on this process 12159be51f97SToby Isaac - cEnd - one after the final cell on this process 12169be51f97SToby Isaac 12171a244344SSatish Balay Level: intermediate 12189be51f97SToby Isaac 12199be51f97SToby Isaac .seealso: DMForestGetCellSF() 12209be51f97SToby Isaac @*/ 1221c7eeac06SToby Isaac PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd) 1222c7eeac06SToby Isaac { 1223c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1224c7eeac06SToby Isaac PetscErrorCode ierr; 1225c7eeac06SToby Isaac 1226c7eeac06SToby Isaac PetscFunctionBegin; 1227c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1228c7eeac06SToby Isaac PetscValidIntPointer(cStart,2); 1229c7eeac06SToby Isaac PetscValidIntPointer(cEnd,2); 1230c7eeac06SToby Isaac if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) { 1231c7eeac06SToby Isaac ierr = forest->createcellchart(dm,&forest->cStart,&forest->cEnd);CHKERRQ(ierr); 1232c7eeac06SToby Isaac } 1233c7eeac06SToby Isaac *cStart = forest->cStart; 1234c7eeac06SToby Isaac *cEnd = forest->cEnd; 1235c7eeac06SToby Isaac PetscFunctionReturn(0); 1236c7eeac06SToby Isaac } 1237c7eeac06SToby Isaac 1238c7eeac06SToby Isaac #undef __FUNCT__ 1239c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellSF" 12409be51f97SToby Isaac /*@ 12419be51f97SToby Isaac DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes 12429be51f97SToby Isaac 12439be51f97SToby Isaac Not collective 12449be51f97SToby Isaac 12459be51f97SToby Isaac Input Parameter: 12469be51f97SToby Isaac . dm - the forest 12479be51f97SToby Isaac 12489be51f97SToby Isaac Output Parameter: 12499be51f97SToby Isaac . cellSF - the PetscSF 12509be51f97SToby Isaac 12511a244344SSatish Balay Level: intermediate 12529be51f97SToby Isaac 12539be51f97SToby Isaac .seealso: DMForestGetCellChart() 12549be51f97SToby Isaac @*/ 1255c7eeac06SToby Isaac PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF) 1256c7eeac06SToby Isaac { 1257c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1258c7eeac06SToby Isaac PetscErrorCode ierr; 1259c7eeac06SToby Isaac 1260c7eeac06SToby Isaac PetscFunctionBegin; 1261c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1262c7eeac06SToby Isaac PetscValidPointer(cellSF,2); 1263c7eeac06SToby Isaac if ((!forest->cellSF) && forest->createcellsf) { 1264c7eeac06SToby Isaac ierr = forest->createcellsf(dm,&forest->cellSF);CHKERRQ(ierr); 1265c7eeac06SToby Isaac } 1266c7eeac06SToby Isaac *cellSF = forest->cellSF; 1267c7eeac06SToby Isaac PetscFunctionReturn(0); 1268c7eeac06SToby Isaac } 1269c7eeac06SToby Isaac 1270c7eeac06SToby Isaac #undef __FUNCT__ 1271ebdf65a2SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityLabel" 12729be51f97SToby Isaac /*@C 12739be51f97SToby Isaac DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see 12749be51f97SToby Isaac DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination). The 12759be51f97SToby Isaac interpretation of the label values is up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and 12769be51f97SToby Isaac DM_FOREST_COARSEN have been reserved as choices that should be accepted by all subtypes. 12779be51f97SToby Isaac 12789be51f97SToby Isaac Logically collective on dm 12799be51f97SToby Isaac 12809be51f97SToby Isaac Input Parameters: 12819be51f97SToby Isaac - dm - the forest 12829be51f97SToby Isaac + adaptLabel - the name of the label in the pre-adaptation forest 12839be51f97SToby Isaac 12849be51f97SToby Isaac Level: intermediate 12859be51f97SToby Isaac 12869be51f97SToby Isaac .seealso DMForestGetAdaptivityLabel() 12879be51f97SToby Isaac @*/ 1288ebdf65a2SToby Isaac PetscErrorCode DMForestSetAdaptivityLabel(DM dm, const char * adaptLabel) 1289c7eeac06SToby Isaac { 1290c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1291c7eeac06SToby Isaac PetscErrorCode ierr; 1292c7eeac06SToby Isaac 1293c7eeac06SToby Isaac PetscFunctionBegin; 1294c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1295ebdf65a2SToby Isaac ierr = PetscFree(forest->adaptLabel);CHKERRQ(ierr); 1296ebdf65a2SToby Isaac ierr = PetscStrallocpy(adaptLabel,&forest->adaptLabel);CHKERRQ(ierr); 1297c7eeac06SToby Isaac PetscFunctionReturn(0); 1298c7eeac06SToby Isaac } 1299c7eeac06SToby Isaac 1300c7eeac06SToby Isaac #undef __FUNCT__ 1301ebdf65a2SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityLabel" 13029be51f97SToby Isaac /*@C 13039be51f97SToby Isaac DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that 13049be51f97SToby Isaac holds the adaptation flags (refinement, coarsening, or some combination). The interpretation of the label values is 13059be51f97SToby Isaac up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and DM_FOREST_COARSEN have been reserved as 13069be51f97SToby Isaac choices that should be accepted by all subtypes. 13079be51f97SToby Isaac 13089be51f97SToby Isaac Not collective 13099be51f97SToby Isaac 13109be51f97SToby Isaac Input Parameter: 13119be51f97SToby Isaac . dm - the forest 13129be51f97SToby Isaac 13139be51f97SToby Isaac Output Parameter: 13149be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest 13159be51f97SToby Isaac 13169be51f97SToby Isaac Level: intermediate 13179be51f97SToby Isaac 13189be51f97SToby Isaac .seealso DMForestSetAdaptivityLabel() 13199be51f97SToby Isaac @*/ 1320ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityLabel(DM dm, const char ** adaptLabel) 1321c7eeac06SToby Isaac { 1322c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1323c7eeac06SToby Isaac 1324c7eeac06SToby Isaac PetscFunctionBegin; 1325c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1326ba936b91SToby Isaac *adaptLabel = forest->adaptLabel; 1327c7eeac06SToby Isaac PetscFunctionReturn(0); 1328c7eeac06SToby Isaac } 1329c7eeac06SToby Isaac 1330c7eeac06SToby Isaac #undef __FUNCT__ 1331c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetCellWeights" 13329be51f97SToby Isaac /*@ 13339be51f97SToby Isaac DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current 13349be51f97SToby Isaac process: weights are used to determine parallel partitioning. Partitions will be created so that each process's 13359be51f97SToby Isaac ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight 13369be51f97SToby Isaac of 1. 13379be51f97SToby Isaac 13389be51f97SToby Isaac Logically collective on dm 13399be51f97SToby Isaac 13409be51f97SToby Isaac Input Parameters: 13419be51f97SToby Isaac + dm - the forest 13429be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1. 13439be51f97SToby Isaac - copyMode - how weights should reference weights 13449be51f97SToby Isaac 13459be51f97SToby Isaac Level: advanced 13469be51f97SToby Isaac 13479be51f97SToby Isaac .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity() 13489be51f97SToby Isaac @*/ 1349c7eeac06SToby Isaac PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode) 1350c7eeac06SToby Isaac { 1351c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1352c7eeac06SToby Isaac PetscInt cStart, cEnd; 1353c7eeac06SToby Isaac PetscErrorCode ierr; 1354c7eeac06SToby Isaac 1355c7eeac06SToby Isaac PetscFunctionBegin; 1356c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1357c7eeac06SToby Isaac ierr = DMForestGetCellChart(dm,&cStart,&cEnd);CHKERRQ(ierr); 1358c7eeac06SToby Isaac if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd); 1359c7eeac06SToby Isaac if (copyMode == PETSC_COPY_VALUES) { 1360c7eeac06SToby Isaac if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) { 1361c7eeac06SToby Isaac ierr = PetscMalloc1(cEnd-cStart,&forest->cellWeights);CHKERRQ(ierr); 1362c7eeac06SToby Isaac } 1363c7eeac06SToby Isaac ierr = PetscMemcpy(forest->cellWeights,weights,(cEnd-cStart)*sizeof(*weights));CHKERRQ(ierr); 1364c7eeac06SToby Isaac forest->cellWeightsCopyMode = PETSC_OWN_POINTER; 1365c7eeac06SToby Isaac PetscFunctionReturn(0); 1366c7eeac06SToby Isaac } 1367c7eeac06SToby Isaac if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) { 1368c7eeac06SToby Isaac ierr = PetscFree(forest->cellWeights);CHKERRQ(ierr); 1369c7eeac06SToby Isaac } 1370c7eeac06SToby Isaac forest->cellWeights = weights; 1371c7eeac06SToby Isaac forest->cellWeightsCopyMode = copyMode; 1372c7eeac06SToby Isaac PetscFunctionReturn(0); 1373c7eeac06SToby Isaac } 1374c7eeac06SToby Isaac 1375c7eeac06SToby Isaac #undef __FUNCT__ 1376c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellWeights" 13779be51f97SToby Isaac /*@ 13789be51f97SToby Isaac DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current 13799be51f97SToby Isaac process: weights are used to determine parallel partitioning. Partitions will be created so that each process's 13809be51f97SToby Isaac ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight 13819be51f97SToby Isaac of 1. 13829be51f97SToby Isaac 13839be51f97SToby Isaac Not collective 13849be51f97SToby Isaac 13859be51f97SToby Isaac Input Parameter: 13869be51f97SToby Isaac . dm - the forest 13879be51f97SToby Isaac 13889be51f97SToby Isaac Output Parameter: 13899be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1. 13909be51f97SToby Isaac 13919be51f97SToby Isaac Level: advanced 13929be51f97SToby Isaac 13939be51f97SToby Isaac .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity() 13949be51f97SToby Isaac @*/ 1395c7eeac06SToby Isaac PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights) 1396c7eeac06SToby Isaac { 1397c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1398c7eeac06SToby Isaac 1399c7eeac06SToby Isaac PetscFunctionBegin; 1400c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1401c7eeac06SToby Isaac PetscValidPointer(weights,2); 1402c7eeac06SToby Isaac *weights = forest->cellWeights; 1403c7eeac06SToby Isaac PetscFunctionReturn(0); 1404c7eeac06SToby Isaac } 1405c7eeac06SToby Isaac 1406c7eeac06SToby Isaac #undef __FUNCT__ 1407c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetWeightCapacity" 14089be51f97SToby Isaac /*@ 14099be51f97SToby Isaac DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning 14109be51f97SToby Isaac a pre-adaptation forest (see DMForestGetAdaptivityForest()). After partitioning, the ratio of the weight of each 14119be51f97SToby Isaac process's cells to the process's capacity will be roughly equal for all processes. A capacity of 0 indicates that 14129be51f97SToby Isaac the current process should not have any cells after repartitioning. 14139be51f97SToby Isaac 14149be51f97SToby Isaac Logically Collective on dm 14159be51f97SToby Isaac 14169be51f97SToby Isaac Input parameters: 14179be51f97SToby Isaac + dm - the forest 14189be51f97SToby Isaac - capacity - this process's capacity 14199be51f97SToby Isaac 14209be51f97SToby Isaac Level: advanced 14219be51f97SToby Isaac 14229be51f97SToby Isaac .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor() 14239be51f97SToby Isaac @*/ 1424c7eeac06SToby Isaac PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity) 1425c7eeac06SToby Isaac { 1426c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1427c7eeac06SToby Isaac 1428c7eeac06SToby Isaac PetscFunctionBegin; 1429c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1430ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup"); 1431c7eeac06SToby Isaac if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity); 1432c7eeac06SToby Isaac forest->weightCapacity = capacity; 1433c7eeac06SToby Isaac PetscFunctionReturn(0); 1434c7eeac06SToby Isaac } 1435c7eeac06SToby Isaac 1436c7eeac06SToby Isaac #undef __FUNCT__ 1437c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetWeightCapacity" 14389be51f97SToby Isaac /*@ 14399be51f97SToby Isaac DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see 14409be51f97SToby Isaac DMForestGetAdaptivityForest()). After partitioning, the ratio of the weight of each process's cells to the 14419be51f97SToby Isaac process's capacity will be roughly equal for all processes. A capacity of 0 indicates that the current process 14429be51f97SToby Isaac should not have any cells after repartitioning. 14439be51f97SToby Isaac 14449be51f97SToby Isaac Not collective 14459be51f97SToby Isaac 14469be51f97SToby Isaac Input parameter: 14479be51f97SToby Isaac . dm - the forest 14489be51f97SToby Isaac 14499be51f97SToby Isaac Output parameter: 14509be51f97SToby Isaac . capacity - this process's capacity 14519be51f97SToby Isaac 14529be51f97SToby Isaac Level: advanced 14539be51f97SToby Isaac 14549be51f97SToby Isaac .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor() 14559be51f97SToby Isaac @*/ 1456c7eeac06SToby Isaac PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity) 1457c7eeac06SToby Isaac { 1458c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1459c7eeac06SToby Isaac 1460c7eeac06SToby Isaac PetscFunctionBegin; 1461c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1462c7eeac06SToby Isaac PetscValidRealPointer(capacity,2); 1463c7eeac06SToby Isaac *capacity = forest->weightCapacity; 1464c7eeac06SToby Isaac PetscFunctionReturn(0); 1465c7eeac06SToby Isaac } 1466c7eeac06SToby Isaac 1467dd8e54a2SToby Isaac #undef __FUNCT__ 1468db4d5e8cSToby Isaac #define __FUNCT__ "DMSetFromOptions_Forest" 14690709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm) 1470db4d5e8cSToby Isaac { 1471db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 147256ba9f64SToby Isaac PetscBool flg, flg1, flg2, flg3, flg4; 1473dd8e54a2SToby Isaac DMForestTopology oldTopo; 1474c7eeac06SToby Isaac char stringBuffer[256]; 1475dd8e54a2SToby Isaac PetscViewer viewer; 1476dd8e54a2SToby Isaac PetscViewerFormat format; 147756ba9f64SToby Isaac PetscInt adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade; 1478c7eeac06SToby Isaac PetscReal weightsFactor; 1479c7eeac06SToby Isaac DMForestAdaptivityStrategy adaptStrategy; 1480db4d5e8cSToby Isaac PetscErrorCode ierr; 1481db4d5e8cSToby Isaac 1482db4d5e8cSToby Isaac PetscFunctionBegin; 1483db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 148458762b62SToby Isaac forest->setfromoptionscalled = PETSC_TRUE; 1485dd8e54a2SToby Isaac ierr = DMForestGetTopology(dm, &oldTopo);CHKERRQ(ierr); 1486a3eda1eaSToby Isaac ierr = PetscOptionsHead(PetscOptionsObject,"DMForest Options");CHKERRQ(ierr); 148756ba9f64SToby Isaac ierr = PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);CHKERRQ(ierr); 148856ba9f64SToby Isaac ierr = PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);CHKERRQ(ierr); 148956ba9f64SToby Isaac ierr = PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);CHKERRQ(ierr); 149056ba9f64SToby Isaac ierr = PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);CHKERRQ(ierr); 1491f885a11aSToby 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}"); 149256ba9f64SToby Isaac if (flg1) { 149356ba9f64SToby Isaac ierr = DMForestSetTopology(dm,(DMForestTopology)stringBuffer);CHKERRQ(ierr); 149456ba9f64SToby Isaac ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr); 149520e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr); 149656ba9f64SToby Isaac } 149756ba9f64SToby Isaac if (flg2) { 1498dd8e54a2SToby Isaac DM base; 1499dd8e54a2SToby Isaac 1500dd8e54a2SToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),&base);CHKERRQ(ierr); 1501dd8e54a2SToby Isaac ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 1502dd8e54a2SToby Isaac ierr = DMLoad(base,viewer);CHKERRQ(ierr); 1503dd8e54a2SToby Isaac ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1504dd8e54a2SToby Isaac ierr = DMForestSetBaseDM(dm,base);CHKERRQ(ierr); 1505dd8e54a2SToby Isaac ierr = DMDestroy(&base);CHKERRQ(ierr); 150656ba9f64SToby Isaac ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr); 150720e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr); 1508dd8e54a2SToby Isaac } 150956ba9f64SToby Isaac if (flg3) { 1510dd8e54a2SToby Isaac DM coarse; 1511dd8e54a2SToby Isaac 1512dd8e54a2SToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),&coarse);CHKERRQ(ierr); 1513dd8e54a2SToby Isaac ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 1514dd8e54a2SToby Isaac ierr = DMLoad(coarse,viewer);CHKERRQ(ierr); 1515dd8e54a2SToby Isaac ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 151620e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,coarse);CHKERRQ(ierr); 1517dd8e54a2SToby Isaac ierr = DMDestroy(&coarse);CHKERRQ(ierr); 151856ba9f64SToby Isaac ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr); 151956ba9f64SToby Isaac ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr); 1520dd8e54a2SToby Isaac } 152156ba9f64SToby Isaac if (flg4) { 1522dd8e54a2SToby Isaac DM fine; 1523dd8e54a2SToby Isaac 1524dd8e54a2SToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),&fine);CHKERRQ(ierr); 1525dd8e54a2SToby Isaac ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 1526dd8e54a2SToby Isaac ierr = DMLoad(fine,viewer);CHKERRQ(ierr); 1527dd8e54a2SToby Isaac ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 152820e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,fine);CHKERRQ(ierr); 1529dd8e54a2SToby Isaac ierr = DMDestroy(&fine);CHKERRQ(ierr); 153056ba9f64SToby Isaac ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr); 153156ba9f64SToby Isaac ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr); 1532dd8e54a2SToby Isaac } 1533dd8e54a2SToby Isaac ierr = DMForestGetAdjacencyDimension(dm,&adjDim);CHKERRQ(ierr); 1534dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg);CHKERRQ(ierr); 1535dd8e54a2SToby Isaac if (flg) { 1536dd8e54a2SToby Isaac ierr = DMForestSetAdjacencyDimension(dm,adjDim);CHKERRQ(ierr); 1537f885a11aSToby Isaac } else { 1538dd8e54a2SToby Isaac ierr = DMForestGetAdjacencyCodimension(dm,&adjCodim);CHKERRQ(ierr); 1539dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg);CHKERRQ(ierr); 1540dd8e54a2SToby Isaac if (flg) { 1541dd8e54a2SToby Isaac ierr = DMForestSetAdjacencyCodimension(dm,adjCodim);CHKERRQ(ierr); 1542dd8e54a2SToby Isaac } 1543dd8e54a2SToby Isaac } 1544dd8e54a2SToby Isaac ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr); 1545dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg);CHKERRQ(ierr); 1546dd8e54a2SToby Isaac if (flg) { 1547dd8e54a2SToby Isaac ierr = DMForestSetPartitionOverlap(dm,overlap);CHKERRQ(ierr); 1548dd8e54a2SToby Isaac } 1549a6121fbdSMatthew G. Knepley #if 0 1550a6121fbdSMatthew 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); 1551a6121fbdSMatthew G. Knepley if (flg) { 1552a6121fbdSMatthew G. Knepley ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr); 1553a6121fbdSMatthew G. Knepley ierr = DMForestSetInitialRefinement(dm,minRefinement);CHKERRQ(ierr); 1554a6121fbdSMatthew G. Knepley } 1555a6121fbdSMatthew 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); 1556a6121fbdSMatthew G. Knepley if (flg) { 1557a6121fbdSMatthew G. Knepley ierr = DMForestSetMinimumRefinement(dm,0);CHKERRQ(ierr); 1558a6121fbdSMatthew G. Knepley ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr); 1559a6121fbdSMatthew G. Knepley } 1560a6121fbdSMatthew G. Knepley #endif 1561dd8e54a2SToby Isaac ierr = DMForestGetMinimumRefinement(dm,&minRefinement);CHKERRQ(ierr); 1562dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg);CHKERRQ(ierr); 1563dd8e54a2SToby Isaac if (flg) { 1564dd8e54a2SToby Isaac ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr); 1565db4d5e8cSToby Isaac } 156656ba9f64SToby Isaac ierr = DMForestGetInitialRefinement(dm,&initRefinement);CHKERRQ(ierr); 156756ba9f64SToby Isaac ierr = PetscOptionsInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg);CHKERRQ(ierr); 156856ba9f64SToby Isaac if (flg) { 156956ba9f64SToby Isaac ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr); 157056ba9f64SToby Isaac } 1571c7eeac06SToby Isaac ierr = DMForestGetMaximumRefinement(dm,&maxRefinement);CHKERRQ(ierr); 1572c7eeac06SToby Isaac ierr = PetscOptionsInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg);CHKERRQ(ierr); 1573c7eeac06SToby Isaac if (flg) { 1574c7eeac06SToby Isaac ierr = DMForestSetMaximumRefinement(dm,maxRefinement);CHKERRQ(ierr); 1575c7eeac06SToby Isaac } 1576c7eeac06SToby Isaac ierr = DMForestGetAdaptivityStrategy(dm,&adaptStrategy);CHKERRQ(ierr); 1577c7eeac06SToby Isaac ierr = PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,256,&flg);CHKERRQ(ierr); 1578c7eeac06SToby Isaac if (flg) { 1579c7eeac06SToby Isaac ierr = DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);CHKERRQ(ierr); 1580c7eeac06SToby Isaac } 1581c7eeac06SToby Isaac ierr = DMForestGetGradeFactor(dm,&grade);CHKERRQ(ierr); 1582c7eeac06SToby Isaac ierr = PetscOptionsInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg);CHKERRQ(ierr); 1583c7eeac06SToby Isaac if (flg) { 1584c7eeac06SToby Isaac ierr = DMForestSetGradeFactor(dm,grade);CHKERRQ(ierr); 1585c7eeac06SToby Isaac } 1586c7eeac06SToby Isaac ierr = DMForestGetCellWeightFactor(dm,&weightsFactor);CHKERRQ(ierr); 1587c7eeac06SToby Isaac ierr = PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);CHKERRQ(ierr); 1588c7eeac06SToby Isaac if (flg) { 1589c7eeac06SToby Isaac ierr = DMForestSetCellWeightFactor(dm,weightsFactor);CHKERRQ(ierr); 1590c7eeac06SToby Isaac } 1591db4d5e8cSToby Isaac ierr = PetscOptionsTail();CHKERRQ(ierr); 1592db4d5e8cSToby Isaac PetscFunctionReturn(0); 1593db4d5e8cSToby Isaac } 1594db4d5e8cSToby Isaac 1595db4d5e8cSToby Isaac #undef __FUNCT__ 1596d8984e3bSMatthew G. Knepley #define __FUNCT__ "DMCreateSubDM_Forest" 1597d8984e3bSMatthew G. Knepley PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm) 1598d8984e3bSMatthew G. Knepley { 1599d8984e3bSMatthew G. Knepley PetscErrorCode ierr; 1600d8984e3bSMatthew G. Knepley 1601d8984e3bSMatthew G. Knepley PetscFunctionBegin; 1602d8984e3bSMatthew G. Knepley if (subdm) {ierr = DMClone(dm, subdm);CHKERRQ(ierr);} 1603d8984e3bSMatthew G. Knepley ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr); 1604d8984e3bSMatthew G. Knepley PetscFunctionReturn(0); 1605d8984e3bSMatthew G. Knepley } 1606d8984e3bSMatthew G. Knepley 1607d8984e3bSMatthew G. Knepley #undef __FUNCT__ 16085421bac9SToby Isaac #define __FUNCT__ "DMRefine_Forest" 16095421bac9SToby Isaac PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined) 16105421bac9SToby Isaac { 16115421bac9SToby Isaac DMLabel refine; 16125421bac9SToby Isaac DM fineDM; 16135421bac9SToby Isaac PetscErrorCode ierr; 16145421bac9SToby Isaac 16155421bac9SToby Isaac PetscFunctionBegin; 16165421bac9SToby Isaac ierr = DMGetFineDM(dm,&fineDM);CHKERRQ(ierr); 16175421bac9SToby Isaac if (fineDM) { 16185421bac9SToby Isaac ierr = PetscObjectReference((PetscObject)fineDM);CHKERRQ(ierr); 16195421bac9SToby Isaac *dmRefined = fineDM; 16205421bac9SToby Isaac PetscFunctionReturn(0); 16215421bac9SToby Isaac } 16225421bac9SToby Isaac ierr = DMForestTemplate(dm,comm,dmRefined);CHKERRQ(ierr); 16235421bac9SToby Isaac ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr); 16245421bac9SToby Isaac if (!refine) { 16255421bac9SToby Isaac ierr = DMCreateLabel(dm,"refine");CHKERRQ(ierr); 16265421bac9SToby Isaac ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr); 16275421bac9SToby Isaac ierr = DMLabelSetDefaultValue(refine,DM_FOREST_REFINE);CHKERRQ(ierr); 16285421bac9SToby Isaac } 16295421bac9SToby Isaac ierr = DMForestSetAdaptivityLabel(*dmRefined,"refine");CHKERRQ(ierr); 16305421bac9SToby Isaac PetscFunctionReturn(0); 16315421bac9SToby Isaac } 16325421bac9SToby Isaac 16335421bac9SToby Isaac #undef __FUNCT__ 16345421bac9SToby Isaac #define __FUNCT__ "DMCoarsen_Forest" 16355421bac9SToby Isaac PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened) 16365421bac9SToby Isaac { 16375421bac9SToby Isaac DMLabel coarsen; 16385421bac9SToby Isaac DM coarseDM; 16395421bac9SToby Isaac PetscErrorCode ierr; 16405421bac9SToby Isaac 16415421bac9SToby Isaac PetscFunctionBegin; 16424098eed7SToby Isaac { 16434098eed7SToby Isaac PetscMPIInt mpiComparison; 16444098eed7SToby Isaac MPI_Comm dmcomm = PetscObjectComm((PetscObject)dm); 16454098eed7SToby Isaac 16464098eed7SToby Isaac ierr = MPI_Comm_compare(comm, dmcomm, &mpiComparison);CHKERRQ(ierr); 1647f885a11aSToby Isaac if (mpiComparison != MPI_IDENT && mpiComparison != MPI_CONGRUENT) SETERRQ(dmcomm,PETSC_ERR_SUP,"No support for different communicators yet"); 16484098eed7SToby Isaac } 16495421bac9SToby Isaac ierr = DMGetCoarseDM(dm,&coarseDM);CHKERRQ(ierr); 16505421bac9SToby Isaac if (coarseDM) { 16515421bac9SToby Isaac ierr = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr); 16525421bac9SToby Isaac *dmCoarsened = coarseDM; 16535421bac9SToby Isaac PetscFunctionReturn(0); 16545421bac9SToby Isaac } 16555421bac9SToby Isaac ierr = DMForestTemplate(dm,comm,dmCoarsened);CHKERRQ(ierr); 16564098eed7SToby Isaac ierr = DMForestSetAdaptivityPurpose(coarseDM,DM_FOREST_COARSEN);CHKERRQ(ierr); 16575421bac9SToby Isaac ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr); 16585421bac9SToby Isaac if (!coarsen) { 16595421bac9SToby Isaac ierr = DMCreateLabel(dm,"coarsen");CHKERRQ(ierr); 16605421bac9SToby Isaac ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr); 16615421bac9SToby Isaac ierr = DMLabelSetDefaultValue(coarsen,DM_FOREST_COARSEN);CHKERRQ(ierr); 16625421bac9SToby Isaac } 16635421bac9SToby Isaac ierr = DMForestSetAdaptivityLabel(*dmCoarsened,"coarsen");CHKERRQ(ierr); 16645421bac9SToby Isaac PetscFunctionReturn(0); 16655421bac9SToby Isaac } 16665421bac9SToby Isaac 16675421bac9SToby Isaac #undef __FUNCT__ 1668d222f98bSToby Isaac #define __FUNCT__ "DMInitialize_Forest" 1669d222f98bSToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm) 1670d222f98bSToby Isaac { 1671d222f98bSToby Isaac PetscErrorCode ierr; 1672d222f98bSToby Isaac 1673d222f98bSToby Isaac PetscFunctionBegin; 1674d222f98bSToby Isaac ierr = PetscMemzero(dm->ops,sizeof(*(dm->ops)));CHKERRQ(ierr); 1675d222f98bSToby Isaac 1676d222f98bSToby Isaac dm->ops->clone = DMClone_Forest; 1677d222f98bSToby Isaac dm->ops->setfromoptions = DMSetFromOptions_Forest; 1678d222f98bSToby Isaac dm->ops->destroy = DMDestroy_Forest; 1679d8984e3bSMatthew G. Knepley dm->ops->createsubdm = DMCreateSubDM_Forest; 16805421bac9SToby Isaac dm->ops->refine = DMRefine_Forest; 16815421bac9SToby Isaac dm->ops->coarsen = DMCoarsen_Forest; 1682d222f98bSToby Isaac PetscFunctionReturn(0); 1683d222f98bSToby Isaac } 1684d222f98bSToby Isaac 16859be51f97SToby Isaac /*MC 16869be51f97SToby 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. 16879be51f97SToby Isaac 16889be51f97SToby Isaac To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the previous mesh whose values indicate which cells should be refined (DM_FOREST_REFINE) or coarsened (DM_FOREST_COARSEN) and how (subtypes are free to allow additional values for things like anisotropic refinement). The name of the label should be given to the *new* mesh with DMForestSetAdaptivityLabel(). 16899be51f97SToby Isaac 16909be51f97SToby Isaac Level: advanced 16919be51f97SToby Isaac 16929be51f97SToby Isaac .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel() 16939be51f97SToby Isaac M*/ 16949be51f97SToby Isaac 1695d222f98bSToby Isaac #undef __FUNCT__ 1696db4d5e8cSToby Isaac #define __FUNCT__ "DMCreate_Forest" 1697db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm) 1698db4d5e8cSToby Isaac { 1699db4d5e8cSToby Isaac DM_Forest *forest; 1700db4d5e8cSToby Isaac PetscErrorCode ierr; 1701db4d5e8cSToby Isaac 1702db4d5e8cSToby Isaac PetscFunctionBegin; 1703db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1704db4d5e8cSToby Isaac ierr = PetscNewLog(dm,&forest);CHKERRQ(ierr); 1705db4d5e8cSToby Isaac dm->dim = 0; 1706db4d5e8cSToby Isaac dm->data = forest; 1707db4d5e8cSToby Isaac forest->refct = 1; 1708db4d5e8cSToby Isaac forest->data = NULL; 170958762b62SToby Isaac forest->setfromoptionscalled = PETSC_FALSE; 1710db4d5e8cSToby Isaac forest->topology = NULL; 1711db4d5e8cSToby Isaac forest->base = NULL; 1712db4d5e8cSToby Isaac forest->adjDim = PETSC_DEFAULT; 1713db4d5e8cSToby Isaac forest->overlap = PETSC_DEFAULT; 1714db4d5e8cSToby Isaac forest->minRefinement = PETSC_DEFAULT; 1715db4d5e8cSToby Isaac forest->maxRefinement = PETSC_DEFAULT; 171656ba9f64SToby Isaac forest->initRefinement = PETSC_DEFAULT; 1717c7eeac06SToby Isaac forest->cStart = PETSC_DETERMINE; 1718c7eeac06SToby Isaac forest->cEnd = PETSC_DETERMINE; 1719db4d5e8cSToby Isaac forest->cellSF = 0; 1720ebdf65a2SToby Isaac forest->adaptLabel = NULL; 1721db4d5e8cSToby Isaac forest->gradeFactor = 2; 1722db4d5e8cSToby Isaac forest->cellWeights = NULL; 1723db4d5e8cSToby Isaac forest->cellWeightsCopyMode = PETSC_USE_POINTER; 1724db4d5e8cSToby Isaac forest->weightsFactor = 1.; 1725db4d5e8cSToby Isaac forest->weightCapacity = 1.; 1726a73e2921SToby Isaac ierr = DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);CHKERRQ(ierr); 1727d222f98bSToby Isaac ierr = DMInitialize_Forest(dm);CHKERRQ(ierr); 1728db4d5e8cSToby Isaac PetscFunctionReturn(0); 1729db4d5e8cSToby Isaac } 1730db4d5e8cSToby Isaac 1731