19be51f97SToby Isaac #include <petsc/private/dmforestimpl.h> /*I "petscdmforest.h" I*/ 29be51f97SToby Isaac #include <petsc/private/dmimpl.h> /*I "petscdm.h" I*/ 3a1b0c543SToby Isaac #include <petsc/private/dmlabelimpl.h> /*I "petscdmlabel.h" I*/ 4ef19d27cSToby Isaac #include <petscsf.h> 5db4d5e8cSToby Isaac 627d4645fSToby Isaac PetscBool DMForestPackageInitialized = PETSC_FALSE; 727d4645fSToby Isaac 827d4645fSToby Isaac typedef struct _DMForestTypeLink*DMForestTypeLink; 927d4645fSToby Isaac 1027d4645fSToby Isaac struct _DMForestTypeLink 1127d4645fSToby Isaac { 1227d4645fSToby Isaac char *name; 1327d4645fSToby Isaac DMForestTypeLink next; 1427d4645fSToby Isaac }; 1527d4645fSToby Isaac 1627d4645fSToby Isaac DMForestTypeLink DMForestTypeList; 1727d4645fSToby Isaac 1827d4645fSToby Isaac static PetscErrorCode DMForestPackageFinalize(void) 1927d4645fSToby Isaac { 2027d4645fSToby Isaac DMForestTypeLink oldLink, link = DMForestTypeList; 2127d4645fSToby Isaac PetscErrorCode ierr; 2227d4645fSToby Isaac 2327d4645fSToby Isaac PetscFunctionBegin; 2427d4645fSToby Isaac while (link) { 2527d4645fSToby Isaac oldLink = link; 26f416af30SBarry Smith ierr = PetscFree(oldLink->name);CHKERRQ(ierr); 2727d4645fSToby Isaac link = oldLink->next; 2827d4645fSToby Isaac ierr = PetscFree(oldLink);CHKERRQ(ierr); 2927d4645fSToby Isaac } 3027d4645fSToby Isaac PetscFunctionReturn(0); 3127d4645fSToby Isaac } 3227d4645fSToby Isaac 3327d4645fSToby Isaac static PetscErrorCode DMForestPackageInitialize(void) 3427d4645fSToby Isaac { 3527d4645fSToby Isaac PetscErrorCode ierr; 3627d4645fSToby Isaac 3727d4645fSToby Isaac PetscFunctionBegin; 3827d4645fSToby Isaac if (DMForestPackageInitialized) PetscFunctionReturn(0); 3927d4645fSToby Isaac DMForestPackageInitialized = PETSC_TRUE; 40f885a11aSToby Isaac 4127d4645fSToby Isaac ierr = DMForestRegisterType(DMFOREST);CHKERRQ(ierr); 4227d4645fSToby Isaac ierr = PetscRegisterFinalize(DMForestPackageFinalize);CHKERRQ(ierr); 4327d4645fSToby Isaac PetscFunctionReturn(0); 4427d4645fSToby Isaac } 4527d4645fSToby Isaac 469be51f97SToby Isaac /*@C 479be51f97SToby Isaac DMForestRegisterType - Registers a DMType as a subtype of DMFOREST (so that DMIsForest() will be correct) 489be51f97SToby Isaac 499be51f97SToby Isaac Not Collective 509be51f97SToby Isaac 519be51f97SToby Isaac Input parameter: 529be51f97SToby Isaac . name - the name of the type 539be51f97SToby Isaac 549be51f97SToby Isaac Level: advanced 559be51f97SToby Isaac 569be51f97SToby Isaac .seealso: DMFOREST, DMIsForest() 579be51f97SToby Isaac @*/ 5827d4645fSToby Isaac PetscErrorCode DMForestRegisterType(DMType name) 5927d4645fSToby Isaac { 6027d4645fSToby Isaac DMForestTypeLink link; 6127d4645fSToby Isaac PetscErrorCode ierr; 6227d4645fSToby Isaac 6327d4645fSToby Isaac PetscFunctionBegin; 6427d4645fSToby Isaac ierr = DMForestPackageInitialize();CHKERRQ(ierr); 6527d4645fSToby Isaac ierr = PetscNew(&link);CHKERRQ(ierr); 6627d4645fSToby Isaac ierr = PetscStrallocpy(name,&link->name);CHKERRQ(ierr); 6727d4645fSToby Isaac link->next = DMForestTypeList; 6827d4645fSToby Isaac DMForestTypeList = link; 6927d4645fSToby Isaac PetscFunctionReturn(0); 7027d4645fSToby Isaac } 7127d4645fSToby Isaac 729be51f97SToby Isaac /*@ 739be51f97SToby Isaac DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes 749be51f97SToby Isaac 759be51f97SToby Isaac Not Collective 769be51f97SToby Isaac 779be51f97SToby Isaac Input parameter: 789be51f97SToby Isaac . dm - the DM object 799be51f97SToby Isaac 809be51f97SToby Isaac Output parameter: 819be51f97SToby Isaac . isForest - whether dm is a subtype of DMFOREST 829be51f97SToby Isaac 839be51f97SToby Isaac Level: intermediate 849be51f97SToby Isaac 859be51f97SToby Isaac .seealso: DMFOREST, DMForestRegisterType() 869be51f97SToby Isaac @*/ 8727d4645fSToby Isaac PetscErrorCode DMIsForest(DM dm, PetscBool *isForest) 8827d4645fSToby Isaac { 8927d4645fSToby Isaac DMForestTypeLink link = DMForestTypeList; 9027d4645fSToby Isaac PetscErrorCode ierr; 9127d4645fSToby Isaac 9227d4645fSToby Isaac PetscFunctionBegin; 9327d4645fSToby Isaac while (link) { 9427d4645fSToby Isaac PetscBool sameType; 9527d4645fSToby Isaac ierr = PetscObjectTypeCompare((PetscObject)dm,link->name,&sameType);CHKERRQ(ierr); 9627d4645fSToby Isaac if (sameType) { 9727d4645fSToby Isaac *isForest = PETSC_TRUE; 9827d4645fSToby Isaac PetscFunctionReturn(0); 9927d4645fSToby Isaac } 10027d4645fSToby Isaac link = link->next; 10127d4645fSToby Isaac } 10227d4645fSToby Isaac *isForest = PETSC_FALSE; 10327d4645fSToby Isaac PetscFunctionReturn(0); 10427d4645fSToby Isaac } 10527d4645fSToby Isaac 1069be51f97SToby Isaac /*@ 1079be51f97SToby Isaac DMForestTemplate - Create a new DM that will be adapted from a source DM. The new DM reproduces the configuration 1089be51f97SToby 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 1099be51f97SToby Isaac (by, e.g., refinement or repartitioning). The source DM is also set as the adaptivity source DM of the new DM (see 1109be51f97SToby Isaac DMForestSetAdaptivityForest()). 1119be51f97SToby Isaac 1129be51f97SToby Isaac Collective on dm 1139be51f97SToby Isaac 1149be51f97SToby Isaac Input Parameters: 1159be51f97SToby Isaac + dm - the source DM object 1169be51f97SToby Isaac - comm - the communicator for the new DM (this communicator is currently ignored, but is present so that DMForestTemplate() can be used within DMCoarsen()) 1179be51f97SToby Isaac 1189be51f97SToby Isaac Output Parameter: 1199be51f97SToby Isaac . tdm - the new DM object 1209be51f97SToby Isaac 1219be51f97SToby Isaac Level: intermediate 1229be51f97SToby Isaac 1239be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest() 1249be51f97SToby Isaac @*/ 1259be51f97SToby Isaac PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm) 126a0452a8eSToby Isaac { 127a0452a8eSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 12820e8089bSToby Isaac DMType type; 129a0452a8eSToby Isaac DM base; 130a0452a8eSToby Isaac DMForestTopology topology; 131*05e99e11SStefano Zampini MatType mtype; 132a0452a8eSToby Isaac PetscInt dim, overlap, ref, factor; 133a0452a8eSToby Isaac DMForestAdaptivityStrategy strat; 134795844e7SToby Isaac PetscDS ds; 135795844e7SToby Isaac void *ctx; 13649fc9a2fSToby Isaac PetscErrorCode (*map)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void*); 1373e58adeeSToby Isaac void *mapCtx; 138a0452a8eSToby Isaac PetscErrorCode ierr; 139a0452a8eSToby Isaac 140a0452a8eSToby Isaac PetscFunctionBegin; 141a0452a8eSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 14220e8089bSToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),tdm);CHKERRQ(ierr); 14320e8089bSToby Isaac ierr = DMGetType(dm,&type);CHKERRQ(ierr); 14420e8089bSToby Isaac ierr = DMSetType(*tdm,type);CHKERRQ(ierr); 145a0452a8eSToby Isaac ierr = DMForestGetBaseDM(dm,&base);CHKERRQ(ierr); 14620e8089bSToby Isaac ierr = DMForestSetBaseDM(*tdm,base);CHKERRQ(ierr); 147a0452a8eSToby Isaac ierr = DMForestGetTopology(dm,&topology);CHKERRQ(ierr); 14820e8089bSToby Isaac ierr = DMForestSetTopology(*tdm,topology);CHKERRQ(ierr); 149a0452a8eSToby Isaac ierr = DMForestGetAdjacencyDimension(dm,&dim);CHKERRQ(ierr); 15020e8089bSToby Isaac ierr = DMForestSetAdjacencyDimension(*tdm,dim);CHKERRQ(ierr); 151a0452a8eSToby Isaac ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr); 15220e8089bSToby Isaac ierr = DMForestSetPartitionOverlap(*tdm,overlap);CHKERRQ(ierr); 153a0452a8eSToby Isaac ierr = DMForestGetMinimumRefinement(dm,&ref);CHKERRQ(ierr); 15420e8089bSToby Isaac ierr = DMForestSetMinimumRefinement(*tdm,ref);CHKERRQ(ierr); 155a0452a8eSToby Isaac ierr = DMForestGetMaximumRefinement(dm,&ref);CHKERRQ(ierr); 15620e8089bSToby Isaac ierr = DMForestSetMaximumRefinement(*tdm,ref);CHKERRQ(ierr); 157a0452a8eSToby Isaac ierr = DMForestGetAdaptivityStrategy(dm,&strat);CHKERRQ(ierr); 15820e8089bSToby Isaac ierr = DMForestSetAdaptivityStrategy(*tdm,strat);CHKERRQ(ierr); 159a0452a8eSToby Isaac ierr = DMForestGetGradeFactor(dm,&factor);CHKERRQ(ierr); 16020e8089bSToby Isaac ierr = DMForestSetGradeFactor(*tdm,factor);CHKERRQ(ierr); 1613e58adeeSToby Isaac ierr = DMForestGetBaseCoordinateMapping(dm,&map,&mapCtx);CHKERRQ(ierr); 1625e8c540aSToby Isaac ierr = DMForestSetBaseCoordinateMapping(*tdm,map,mapCtx);CHKERRQ(ierr); 163a0452a8eSToby Isaac if (forest->ftemplate) { 16420e8089bSToby Isaac ierr = (forest->ftemplate)(dm, *tdm);CHKERRQ(ierr); 165a0452a8eSToby Isaac } 16620e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(*tdm,dm);CHKERRQ(ierr); 167795844e7SToby Isaac ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 168795844e7SToby Isaac ierr = DMSetDS(*tdm,ds);CHKERRQ(ierr); 169795844e7SToby Isaac ierr = DMGetApplicationContext(dm,&ctx);CHKERRQ(ierr); 170795844e7SToby Isaac ierr = DMSetApplicationContext(*tdm,&ctx);CHKERRQ(ierr); 17190b157c4SStefano Zampini { 17290b157c4SStefano Zampini PetscBool isper; 173795844e7SToby Isaac const PetscReal *maxCell, *L; 174795844e7SToby Isaac const DMBoundaryType *bd; 175795844e7SToby Isaac 17690b157c4SStefano Zampini ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr); 17790b157c4SStefano Zampini ierr = DMSetPeriodicity(*tdm,isper,maxCell,L,bd);CHKERRQ(ierr); 178795844e7SToby Isaac } 179bff67a9bSToby Isaac ierr = DMCopyBoundary(dm,*tdm);CHKERRQ(ierr); 180*05e99e11SStefano Zampini ierr = DMGetMatType(dm,&mtype);CHKERRQ(ierr); 181*05e99e11SStefano Zampini ierr = DMSetMatType(*tdm,mtype);CHKERRQ(ierr); 182a0452a8eSToby Isaac PetscFunctionReturn(0); 183a0452a8eSToby Isaac } 184a0452a8eSToby Isaac 18501d9d024SToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm); 18601d9d024SToby Isaac 187db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm) 188db4d5e8cSToby Isaac { 189db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 190db4d5e8cSToby Isaac const char *type; 191db4d5e8cSToby Isaac PetscErrorCode ierr; 192db4d5e8cSToby Isaac 193db4d5e8cSToby Isaac PetscFunctionBegin; 194db4d5e8cSToby Isaac forest->refct++; 195db4d5e8cSToby Isaac (*newdm)->data = forest; 196db4d5e8cSToby Isaac ierr = PetscObjectGetType((PetscObject) dm, &type);CHKERRQ(ierr); 197db4d5e8cSToby Isaac ierr = PetscObjectChangeTypeName((PetscObject) *newdm, type);CHKERRQ(ierr); 19801d9d024SToby Isaac ierr = DMInitialize_Forest(*newdm);CHKERRQ(ierr); 199db4d5e8cSToby Isaac PetscFunctionReturn(0); 200db4d5e8cSToby Isaac } 201db4d5e8cSToby Isaac 202d222f98bSToby Isaac static PetscErrorCode DMDestroy_Forest(DM dm) 203db4d5e8cSToby Isaac { 204db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 205db4d5e8cSToby Isaac PetscErrorCode ierr; 206db4d5e8cSToby Isaac 207db4d5e8cSToby Isaac PetscFunctionBegin; 208db4d5e8cSToby Isaac if (--forest->refct > 0) PetscFunctionReturn(0); 209d222f98bSToby Isaac if (forest->destroy) {ierr = forest->destroy(dm);CHKERRQ(ierr);} 210db4d5e8cSToby Isaac ierr = PetscSFDestroy(&forest->cellSF);CHKERRQ(ierr); 2110f17b9e3SToby Isaac ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr); 2120f17b9e3SToby Isaac ierr = PetscSFDestroy(&forest->coarseToPreFine);CHKERRQ(ierr); 213a1b0c543SToby Isaac ierr = DMLabelDestroy(&forest->adaptLabel);CHKERRQ(ierr); 2149a81d013SToby Isaac ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr); 21556ba9f64SToby Isaac ierr = DMDestroy(&forest->base);CHKERRQ(ierr); 216ba936b91SToby Isaac ierr = DMDestroy(&forest->adapt);CHKERRQ(ierr); 21730f902e7SToby Isaac ierr = PetscFree(forest->topology);CHKERRQ(ierr); 21830f902e7SToby Isaac ierr = PetscFree(forest);CHKERRQ(ierr); 219db4d5e8cSToby Isaac PetscFunctionReturn(0); 220db4d5e8cSToby Isaac } 221db4d5e8cSToby Isaac 2229be51f97SToby Isaac /*@C 2239be51f97SToby Isaac DMForestSetTopology - Set the topology of a DMForest during the pre-setup phase. The topology is a string (e.g. 22468d54884SBarry Smith "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest during 2259be51f97SToby Isaac DMSetUp(). 2269be51f97SToby Isaac 2279be51f97SToby Isaac Logically collective on dm 2289be51f97SToby Isaac 2299be51f97SToby Isaac Input parameters: 2309be51f97SToby Isaac + dm - the forest 2319be51f97SToby Isaac - topology - the topology of the forest 2329be51f97SToby Isaac 2339be51f97SToby Isaac Level: intermediate 2349be51f97SToby Isaac 2359be51f97SToby Isaac .seealso(): DMForestGetTopology(), DMForestSetBaseDM() 2369be51f97SToby Isaac @*/ 237dd8e54a2SToby Isaac PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology) 238db4d5e8cSToby Isaac { 239db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 240db4d5e8cSToby Isaac PetscErrorCode ierr; 241db4d5e8cSToby Isaac 242db4d5e8cSToby Isaac PetscFunctionBegin; 243db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 244ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the topology after setup"); 245dd8e54a2SToby Isaac ierr = PetscFree(forest->topology);CHKERRQ(ierr); 246dd8e54a2SToby Isaac ierr = PetscStrallocpy((const char*)topology,(char**) &forest->topology);CHKERRQ(ierr); 247db4d5e8cSToby Isaac PetscFunctionReturn(0); 248db4d5e8cSToby Isaac } 249db4d5e8cSToby Isaac 2509be51f97SToby Isaac /*@C 2519be51f97SToby Isaac DMForestGetTopology - Get a string describing the topology of a DMForest. 2529be51f97SToby Isaac 2539be51f97SToby Isaac Not collective 2549be51f97SToby Isaac 2559be51f97SToby Isaac Input parameter: 2569be51f97SToby Isaac . dm - the forest 2579be51f97SToby Isaac 2589be51f97SToby Isaac Output parameter: 2599be51f97SToby Isaac . topology - the topology of the forest (e.g., 'cube', 'shell') 2609be51f97SToby Isaac 2619be51f97SToby Isaac Level: intermediate 2629be51f97SToby Isaac 2639be51f97SToby Isaac .seealso: DMForestSetTopology() 2649be51f97SToby Isaac @*/ 265dd8e54a2SToby Isaac PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology) 266dd8e54a2SToby Isaac { 267dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 268dd8e54a2SToby Isaac 269dd8e54a2SToby Isaac PetscFunctionBegin; 270dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 271dd8e54a2SToby Isaac PetscValidPointer(topology,2); 272dd8e54a2SToby Isaac *topology = forest->topology; 273dd8e54a2SToby Isaac PetscFunctionReturn(0); 274dd8e54a2SToby Isaac } 275dd8e54a2SToby Isaac 2769be51f97SToby Isaac /*@ 2779be51f97SToby Isaac DMForestSetBaseDM - During the pre-setup phase, set the DM that defines the base mesh of a DMForest forest. The 2789be51f97SToby Isaac forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its 279765b024eSBarry Smith base. In general, two forest must share a base to be comparable, to do things like construct interpolators. 2809be51f97SToby Isaac 2819be51f97SToby Isaac Logically collective on dm 2829be51f97SToby Isaac 2839be51f97SToby Isaac Input Parameters: 2849be51f97SToby Isaac + dm - the forest 2859be51f97SToby Isaac - base - the base DM of the forest 2869be51f97SToby Isaac 28795452b02SPatrick Sanan Notes: 28895452b02SPatrick Sanan Currently the base DM must be a DMPLEX 289765b024eSBarry Smith 2909be51f97SToby Isaac Level: intermediate 2919be51f97SToby Isaac 2929be51f97SToby Isaac .seealso(): DMForestGetBaseDM() 2939be51f97SToby Isaac @*/ 294dd8e54a2SToby Isaac PetscErrorCode DMForestSetBaseDM(DM dm, DM base) 295dd8e54a2SToby Isaac { 296dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 297dd8e54a2SToby Isaac PetscInt dim, dimEmbed; 298dd8e54a2SToby Isaac PetscErrorCode ierr; 299dd8e54a2SToby Isaac 300dd8e54a2SToby Isaac PetscFunctionBegin; 301dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 302ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the base after setup"); 303dd8e54a2SToby Isaac ierr = PetscObjectReference((PetscObject)base);CHKERRQ(ierr); 304dd8e54a2SToby Isaac ierr = DMDestroy(&forest->base);CHKERRQ(ierr); 305dd8e54a2SToby Isaac forest->base = base; 306a0452a8eSToby Isaac if (base) { 307a0452a8eSToby Isaac PetscValidHeaderSpecific(base, DM_CLASSID, 2); 308dd8e54a2SToby Isaac ierr = DMGetDimension(base,&dim);CHKERRQ(ierr); 309dd8e54a2SToby Isaac ierr = DMSetDimension(dm,dim);CHKERRQ(ierr); 310dd8e54a2SToby Isaac ierr = DMGetCoordinateDim(base,&dimEmbed);CHKERRQ(ierr); 311dd8e54a2SToby Isaac ierr = DMSetCoordinateDim(dm,dimEmbed);CHKERRQ(ierr); 312a0452a8eSToby Isaac } 313dd8e54a2SToby Isaac PetscFunctionReturn(0); 314dd8e54a2SToby Isaac } 315dd8e54a2SToby Isaac 3169be51f97SToby Isaac /*@ 3179be51f97SToby Isaac DMForestGetBaseDM - Get the base DM of a DMForest forest. The forest will be hierarchically refined from the base, 31868d54884SBarry Smith and all refinements/coarsenings of the forest will share its base. In general, two forest must share a base to be 3199be51f97SToby Isaac comparable, to do things like construct interpolators. 3209be51f97SToby Isaac 3219be51f97SToby Isaac Not collective 3229be51f97SToby Isaac 3239be51f97SToby Isaac Input Parameter: 3249be51f97SToby Isaac . dm - the forest 3259be51f97SToby Isaac 3269be51f97SToby Isaac Output Parameter: 3279be51f97SToby Isaac . base - the base DM of the forest 3289be51f97SToby Isaac 329367003a6SStefano Zampini Notes: 330367003a6SStefano Zampini After DMSetUp(), the base DM will be redundantly distributed across MPI processes 331367003a6SStefano Zampini 3329be51f97SToby Isaac Level: intermediate 3339be51f97SToby Isaac 3349be51f97SToby Isaac .seealso(); DMForestSetBaseDM() 3359be51f97SToby Isaac @*/ 336dd8e54a2SToby Isaac PetscErrorCode DMForestGetBaseDM(DM dm, DM *base) 337dd8e54a2SToby Isaac { 338dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 339dd8e54a2SToby Isaac 340dd8e54a2SToby Isaac PetscFunctionBegin; 341dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 342dd8e54a2SToby Isaac PetscValidPointer(base, 2); 343dd8e54a2SToby Isaac *base = forest->base; 344dd8e54a2SToby Isaac PetscFunctionReturn(0); 345dd8e54a2SToby Isaac } 346dd8e54a2SToby Isaac 34799478f86SToby Isaac PetscErrorCode DMForestSetBaseCoordinateMapping(DM dm, PetscErrorCode (*func)(DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx) 348cf38a08cSToby Isaac { 349cf38a08cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 350cf38a08cSToby Isaac 351cf38a08cSToby Isaac PetscFunctionBegin; 352cf38a08cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 353cf38a08cSToby Isaac forest->mapcoordinates = func; 354cf38a08cSToby Isaac forest->mapcoordinatesctx = ctx; 355cf38a08cSToby Isaac PetscFunctionReturn(0); 356cf38a08cSToby Isaac } 357cf38a08cSToby Isaac 35899478f86SToby Isaac PetscErrorCode DMForestGetBaseCoordinateMapping(DM dm, PetscErrorCode (**func) (DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx) 359cf38a08cSToby Isaac { 360cf38a08cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 361cf38a08cSToby Isaac 362cf38a08cSToby Isaac PetscFunctionBegin; 363cf38a08cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 364cf38a08cSToby Isaac if (func) *func = forest->mapcoordinates; 365cf38a08cSToby Isaac if (ctx) *((void**) ctx) = forest->mapcoordinatesctx; 366cf38a08cSToby Isaac PetscFunctionReturn(0); 367cf38a08cSToby Isaac } 368cf38a08cSToby Isaac 3699be51f97SToby Isaac /*@ 3709be51f97SToby Isaac DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be 3719be51f97SToby Isaac adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) im DMSetUp(). Usually not needed 3729be51f97SToby Isaac by users directly: DMForestTemplate() constructs a new forest to be adapted from an old forest and calls this 3739be51f97SToby Isaac routine. 3749be51f97SToby Isaac 375dffe73a3SToby Isaac Note that this can be called after setup with adapt = NULL, which will clear all internal data related to the 376dffe73a3SToby Isaac adaptivity forest from dm. This way, repeatedly adapting does not leave stale DM objects in memory. 377dffe73a3SToby Isaac 3789be51f97SToby Isaac Logically collective on dm 3799be51f97SToby Isaac 3809be51f97SToby Isaac Input Parameter: 3819be51f97SToby Isaac + dm - the new forest, which will be constructed from adapt 3829be51f97SToby Isaac - adapt - the old forest 3839be51f97SToby Isaac 3849be51f97SToby Isaac Level: intermediate 3859be51f97SToby Isaac 3869be51f97SToby Isaac .seealso: DMForestGetAdaptivityForest(), DMForestSetAdaptivityPurpose() 3879be51f97SToby Isaac @*/ 388ba936b91SToby Isaac PetscErrorCode DMForestSetAdaptivityForest(DM dm,DM adapt) 389dd8e54a2SToby Isaac { 390dffe73a3SToby Isaac DM_Forest *forest, *adaptForest, *oldAdaptForest; 391dffe73a3SToby Isaac DM oldAdapt; 392456cc5b7SMatthew G. Knepley PetscBool isForest; 393dd8e54a2SToby Isaac PetscErrorCode ierr; 394dd8e54a2SToby Isaac 395dd8e54a2SToby Isaac PetscFunctionBegin; 396dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 397ba936b91SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 398456cc5b7SMatthew G. Knepley ierr = DMIsForest(dm, &isForest);CHKERRQ(ierr); 399456cc5b7SMatthew G. Knepley if (!isForest) PetscFunctionReturn(0); 400ba936b91SToby Isaac forest = (DM_Forest*) dm->data; 401dffe73a3SToby Isaac ierr = DMForestGetAdaptivityForest(dm,&oldAdapt);CHKERRQ(ierr); 402dffe73a3SToby Isaac if (adapt != NULL && dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup"); 403193eb951SToby Isaac adaptForest = (DM_Forest*) (adapt ? adapt->data : NULL); 404193eb951SToby Isaac oldAdaptForest = (DM_Forest*) (oldAdapt ? oldAdapt->data : NULL); 405dffe73a3SToby Isaac if (adaptForest != oldAdaptForest) { 406dffe73a3SToby Isaac ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr); 407dffe73a3SToby Isaac ierr = PetscSFDestroy(&forest->coarseToPreFine);CHKERRQ(ierr); 408dffe73a3SToby Isaac if (forest->clearadaptivityforest) {ierr = (forest->clearadaptivityforest)(dm);CHKERRQ(ierr);} 409dffe73a3SToby Isaac } 41026d9498aSToby Isaac switch (forest->adaptPurpose) { 411cd3c525cSToby Isaac case DM_ADAPT_DETERMINE: 412ba936b91SToby Isaac ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr); 413ba936b91SToby Isaac ierr = DMDestroy(&(forest->adapt));CHKERRQ(ierr); 414ba936b91SToby Isaac forest->adapt = adapt; 41526d9498aSToby Isaac break; 416a1b0c543SToby Isaac case DM_ADAPT_REFINE: 41726d9498aSToby Isaac ierr = DMSetCoarseDM(dm,adapt);CHKERRQ(ierr); 41826d9498aSToby Isaac break; 419a1b0c543SToby Isaac case DM_ADAPT_COARSEN: 42026d9498aSToby Isaac ierr = DMSetFineDM(dm,adapt);CHKERRQ(ierr); 42126d9498aSToby Isaac break; 42226d9498aSToby Isaac default: 42326d9498aSToby Isaac SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose"); 42426d9498aSToby Isaac } 425dd8e54a2SToby Isaac PetscFunctionReturn(0); 426dd8e54a2SToby Isaac } 427dd8e54a2SToby Isaac 4289be51f97SToby Isaac /*@ 4299be51f97SToby Isaac DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted. 4309be51f97SToby Isaac 4319be51f97SToby Isaac Not collective 4329be51f97SToby Isaac 4339be51f97SToby Isaac Input Parameter: 4349be51f97SToby Isaac . dm - the forest 4359be51f97SToby Isaac 4369be51f97SToby Isaac Output Parameter: 4379be51f97SToby Isaac . adapt - the forest from which dm is/was adapted 4389be51f97SToby Isaac 4399be51f97SToby Isaac Level: intermediate 4409be51f97SToby Isaac 4419be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose() 4429be51f97SToby Isaac @*/ 443ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt) 444dd8e54a2SToby Isaac { 445ba936b91SToby Isaac DM_Forest *forest; 44626d9498aSToby Isaac PetscErrorCode ierr; 447dd8e54a2SToby Isaac 448dd8e54a2SToby Isaac PetscFunctionBegin; 449dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 450ba936b91SToby Isaac forest = (DM_Forest*) dm->data; 45126d9498aSToby Isaac switch (forest->adaptPurpose) { 452cd3c525cSToby Isaac case DM_ADAPT_DETERMINE: 453ba936b91SToby Isaac *adapt = forest->adapt; 45426d9498aSToby Isaac break; 455a1b0c543SToby Isaac case DM_ADAPT_REFINE: 45626d9498aSToby Isaac ierr = DMGetCoarseDM(dm,adapt);CHKERRQ(ierr); 45726d9498aSToby Isaac break; 458a1b0c543SToby Isaac case DM_ADAPT_COARSEN: 45926d9498aSToby Isaac ierr = DMGetFineDM(dm,adapt);CHKERRQ(ierr); 46026d9498aSToby Isaac break; 46126d9498aSToby Isaac default: 46226d9498aSToby Isaac SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose"); 46326d9498aSToby Isaac } 46426d9498aSToby Isaac PetscFunctionReturn(0); 46526d9498aSToby Isaac } 46626d9498aSToby Isaac 4679be51f97SToby Isaac /*@ 4689be51f97SToby Isaac DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its 469a1b0c543SToby Isaac source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening 470cd3c525cSToby Isaac (DM_ADAPT_COARSEN), or undefined (DM_ADAPT_DETERMINE). This only matters for the purposes of reference counting: 4719be51f97SToby Isaac during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse 4729be51f97SToby Isaac relationship (see DMSetFineDM()/DMSetCoarseDM()). If the purpose is not refinement or coarsening, and the user does 4739be51f97SToby Isaac not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can 4749be51f97SToby Isaac cause a memory leak. This method is used by subtypes of DMForest when automatically constructing mesh hierarchies. 4759be51f97SToby Isaac 4769be51f97SToby Isaac Logically collective on dm 4779be51f97SToby Isaac 4789be51f97SToby Isaac Input Parameters: 4799be51f97SToby Isaac + dm - the forest 480cd3c525cSToby Isaac - purpose - the adaptivity purpose (DM_ADAPT_DETERMINE/DM_ADAPT_REFINE/DM_ADAPT_COARSEN) 4819be51f97SToby Isaac 4829be51f97SToby Isaac Level: advanced 4839be51f97SToby Isaac 4849be51f97SToby Isaac .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest() 4859be51f97SToby Isaac @*/ 486a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMAdaptFlag purpose) 48726d9498aSToby Isaac { 48826d9498aSToby Isaac DM_Forest *forest; 48926d9498aSToby Isaac PetscErrorCode ierr; 49026d9498aSToby Isaac 49126d9498aSToby Isaac PetscFunctionBegin; 49226d9498aSToby Isaac forest = (DM_Forest*) dm->data; 4939be51f97SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup"); 49426d9498aSToby Isaac if (purpose != forest->adaptPurpose) { 49526d9498aSToby Isaac DM adapt; 49626d9498aSToby Isaac 49726d9498aSToby Isaac ierr = DMForestGetAdaptivityForest(dm,&adapt);CHKERRQ(ierr); 49826d9498aSToby Isaac ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr); 49926d9498aSToby Isaac ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr); 500f885a11aSToby Isaac 50126d9498aSToby Isaac forest->adaptPurpose = purpose; 502f885a11aSToby Isaac 50326d9498aSToby Isaac ierr = DMForestSetAdaptivityForest(dm,adapt);CHKERRQ(ierr); 50426d9498aSToby Isaac ierr = DMDestroy(&adapt);CHKERRQ(ierr); 50526d9498aSToby Isaac } 506dd8e54a2SToby Isaac PetscFunctionReturn(0); 507dd8e54a2SToby Isaac } 508dd8e54a2SToby Isaac 50956c0450aSToby Isaac /*@ 51056c0450aSToby Isaac DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with 511a1b0c543SToby Isaac DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening (DM_ADAPT_COARSEN), or 512cd3c525cSToby Isaac undefined (DM_ADAPT_DETERMINE). This only matters for the purposes of reference counting: during DMDestroy(), cyclic 51356c0450aSToby Isaac references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see 51456c0450aSToby Isaac DMSetFineDM()/DMSetCoarseDM()). If the purpose is not refinement or coarsening, and the user does not maintain a 51556c0450aSToby Isaac reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory 51656c0450aSToby Isaac leak. This method is used by subtypes of DMForest when automatically constructing mesh hierarchies. 51756c0450aSToby Isaac 51856c0450aSToby Isaac Not collective 51956c0450aSToby Isaac 52056c0450aSToby Isaac Input Parameter: 52156c0450aSToby Isaac . dm - the forest 52256c0450aSToby Isaac 52356c0450aSToby Isaac Output Parameter: 524cd3c525cSToby Isaac . purpose - the adaptivity purpose (DM_ADAPT_DETERMINE/DM_ADAPT_REFINE/DM_ADAPT_COARSEN) 52556c0450aSToby Isaac 52656c0450aSToby Isaac Level: advanced 52756c0450aSToby Isaac 52856c0450aSToby Isaac .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest() 52956c0450aSToby Isaac @*/ 530a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMAdaptFlag *purpose) 53156c0450aSToby Isaac { 53256c0450aSToby Isaac DM_Forest *forest; 53356c0450aSToby Isaac 53456c0450aSToby Isaac PetscFunctionBegin; 53556c0450aSToby Isaac forest = (DM_Forest*) dm->data; 53656c0450aSToby Isaac *purpose = forest->adaptPurpose; 53756c0450aSToby Isaac PetscFunctionReturn(0); 53856c0450aSToby Isaac } 53956c0450aSToby Isaac 5409be51f97SToby Isaac /*@ 5419be51f97SToby Isaac DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine 5429be51f97SToby Isaac cell adjacency (for the purposes of partitioning and overlap). 5439be51f97SToby Isaac 5449be51f97SToby Isaac Logically collective on dm 5459be51f97SToby Isaac 5469be51f97SToby Isaac Input Parameters: 5479be51f97SToby Isaac + dm - the forest 5489be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency) 5499be51f97SToby Isaac 5509be51f97SToby Isaac Level: intermediate 5519be51f97SToby Isaac 5529be51f97SToby Isaac .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap() 5539be51f97SToby Isaac @*/ 554dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim) 555dd8e54a2SToby Isaac { 556dd8e54a2SToby Isaac PetscInt dim; 557dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 558dd8e54a2SToby Isaac PetscErrorCode ierr; 559dd8e54a2SToby Isaac 560dd8e54a2SToby Isaac PetscFunctionBegin; 561dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 562ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup"); 563dd8e54a2SToby Isaac if (adjDim < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %d", adjDim); 564dd8e54a2SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 565dd8e54a2SToby Isaac if (adjDim > dim) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %d: %d", dim, adjDim); 566dd8e54a2SToby Isaac forest->adjDim = adjDim; 567dd8e54a2SToby Isaac PetscFunctionReturn(0); 568dd8e54a2SToby Isaac } 569dd8e54a2SToby Isaac 5709be51f97SToby Isaac /*@ 5719be51f97SToby Isaac DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that, 5729be51f97SToby Isaac e.g., adjacency based on facets can be specified by codimension 1 in all cases) 5739be51f97SToby Isaac 5749be51f97SToby Isaac Logically collective on dm 5759be51f97SToby Isaac 5769be51f97SToby Isaac Input Parameters: 5779be51f97SToby Isaac + dm - the forest 5789be51f97SToby Isaac - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices 5799be51f97SToby Isaac 5809be51f97SToby Isaac Level: intermediate 5819be51f97SToby Isaac 5829be51f97SToby Isaac .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension() 5839be51f97SToby Isaac @*/ 584dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim) 585dd8e54a2SToby Isaac { 586dd8e54a2SToby Isaac PetscInt dim; 587dd8e54a2SToby Isaac PetscErrorCode ierr; 588dd8e54a2SToby Isaac 589dd8e54a2SToby Isaac PetscFunctionBegin; 590dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 591dd8e54a2SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 592dd8e54a2SToby Isaac ierr = DMForestSetAdjacencyDimension(dm,dim-adjCodim);CHKERRQ(ierr); 593dd8e54a2SToby Isaac PetscFunctionReturn(0); 594dd8e54a2SToby Isaac } 595dd8e54a2SToby Isaac 5969be51f97SToby Isaac /*@ 5979be51f97SToby Isaac DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the 5989be51f97SToby Isaac purposes of partitioning and overlap). 5999be51f97SToby Isaac 6009be51f97SToby Isaac Not collective 6019be51f97SToby Isaac 6029be51f97SToby Isaac Input Parameter: 6039be51f97SToby Isaac . dm - the forest 6049be51f97SToby Isaac 6059be51f97SToby Isaac Output Parameter: 6069be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency) 6079be51f97SToby Isaac 6089be51f97SToby Isaac Level: intermediate 6099be51f97SToby Isaac 6109be51f97SToby Isaac .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap() 6119be51f97SToby Isaac @*/ 612dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim) 613dd8e54a2SToby Isaac { 614dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 615dd8e54a2SToby Isaac 616dd8e54a2SToby Isaac PetscFunctionBegin; 617dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 618dd8e54a2SToby Isaac PetscValidIntPointer(adjDim,2); 619dd8e54a2SToby Isaac *adjDim = forest->adjDim; 620dd8e54a2SToby Isaac PetscFunctionReturn(0); 621dd8e54a2SToby Isaac } 622dd8e54a2SToby Isaac 6239be51f97SToby Isaac /*@ 6249be51f97SToby Isaac DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that, 6259be51f97SToby Isaac e.g., adjacency based on facets can be specified by codimension 1 in all cases) 6269be51f97SToby Isaac 6279be51f97SToby Isaac Not collective 6289be51f97SToby Isaac 6299be51f97SToby Isaac Input Parameter: 6309be51f97SToby Isaac . dm - the forest 6319be51f97SToby Isaac 6329be51f97SToby Isaac Output Parameter: 6339be51f97SToby Isaac . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices 6349be51f97SToby Isaac 6359be51f97SToby Isaac Level: intermediate 6369be51f97SToby Isaac 6379be51f97SToby Isaac .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension() 6389be51f97SToby Isaac @*/ 639dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim) 640dd8e54a2SToby Isaac { 641dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 642dd8e54a2SToby Isaac PetscInt dim; 643dd8e54a2SToby Isaac PetscErrorCode ierr; 644dd8e54a2SToby Isaac 645dd8e54a2SToby Isaac PetscFunctionBegin; 646dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 647dd8e54a2SToby Isaac PetscValidIntPointer(adjCodim,2); 648dd8e54a2SToby Isaac ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 649dd8e54a2SToby Isaac *adjCodim = dim - forest->adjDim; 650dd8e54a2SToby Isaac PetscFunctionReturn(0); 651dd8e54a2SToby Isaac } 652dd8e54a2SToby Isaac 6539be51f97SToby Isaac /*@ 6549be51f97SToby Isaac DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel 6559be51f97SToby Isaac partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding 6569be51f97SToby Isaac adjacent cells 6579be51f97SToby Isaac 6589be51f97SToby Isaac Logically collective on dm 6599be51f97SToby Isaac 6609be51f97SToby Isaac Input Parameters: 6619be51f97SToby Isaac + dm - the forest 6629be51f97SToby Isaac - overlap - default 0 6639be51f97SToby Isaac 6649be51f97SToby Isaac Level: intermediate 6659be51f97SToby Isaac 6669be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension() 6679be51f97SToby Isaac @*/ 668dd8e54a2SToby Isaac PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap) 669dd8e54a2SToby Isaac { 670dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 671dd8e54a2SToby Isaac 672dd8e54a2SToby Isaac PetscFunctionBegin; 673dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 674ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup"); 675dd8e54a2SToby Isaac if (overlap < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %d", overlap); 676dd8e54a2SToby Isaac forest->overlap = overlap; 677dd8e54a2SToby Isaac PetscFunctionReturn(0); 678dd8e54a2SToby Isaac } 679dd8e54a2SToby Isaac 6809be51f97SToby Isaac /*@ 6819be51f97SToby Isaac DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values 6829be51f97SToby Isaac > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells 6839be51f97SToby Isaac 6849be51f97SToby Isaac Not collective 6859be51f97SToby Isaac 6869be51f97SToby Isaac Input Parameter: 6879be51f97SToby Isaac . dm - the forest 6889be51f97SToby Isaac 6899be51f97SToby Isaac Output Parameter: 6909be51f97SToby Isaac . overlap - default 0 6919be51f97SToby Isaac 6929be51f97SToby Isaac Level: intermediate 6939be51f97SToby Isaac 6949be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension() 6959be51f97SToby Isaac @*/ 696dd8e54a2SToby Isaac PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap) 697dd8e54a2SToby Isaac { 698dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 699dd8e54a2SToby Isaac 700dd8e54a2SToby Isaac PetscFunctionBegin; 701dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 702dd8e54a2SToby Isaac PetscValidIntPointer(overlap,2); 703dd8e54a2SToby Isaac *overlap = forest->overlap; 704dd8e54a2SToby Isaac PetscFunctionReturn(0); 705dd8e54a2SToby Isaac } 706dd8e54a2SToby Isaac 7079be51f97SToby Isaac /*@ 7089be51f97SToby Isaac DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base 7099be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. If the forest is being created by coarsening a previous forest 7109be51f97SToby Isaac (see DMForestGetAdaptivityForest()) this limits the amount of coarsening. 7119be51f97SToby Isaac 7129be51f97SToby Isaac Logically collective on dm 7139be51f97SToby Isaac 7149be51f97SToby Isaac Input Parameters: 7159be51f97SToby Isaac + dm - the forest 7169be51f97SToby Isaac - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7179be51f97SToby Isaac 7189be51f97SToby Isaac Level: intermediate 7199be51f97SToby Isaac 7209be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 7219be51f97SToby Isaac @*/ 722dd8e54a2SToby Isaac PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement) 723dd8e54a2SToby Isaac { 724dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 725dd8e54a2SToby Isaac 726dd8e54a2SToby Isaac PetscFunctionBegin; 727dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 728ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup"); 729dd8e54a2SToby Isaac forest->minRefinement = minRefinement; 730dd8e54a2SToby Isaac PetscFunctionReturn(0); 731dd8e54a2SToby Isaac } 732dd8e54a2SToby Isaac 7339be51f97SToby Isaac /*@ 7349be51f97SToby Isaac DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see 7359be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. If the forest is being created by coarsening a previous forest (see 7369be51f97SToby Isaac DMForestGetAdaptivityForest()), this limits the amount of coarsening. 7379be51f97SToby Isaac 7389be51f97SToby Isaac Not collective 7399be51f97SToby Isaac 7409be51f97SToby Isaac Input Parameter: 7419be51f97SToby Isaac . dm - the forest 7429be51f97SToby Isaac 7439be51f97SToby Isaac Output Parameter: 7449be51f97SToby Isaac . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7459be51f97SToby Isaac 7469be51f97SToby Isaac Level: intermediate 7479be51f97SToby Isaac 7489be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestGetMaximumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 7499be51f97SToby Isaac @*/ 750dd8e54a2SToby Isaac PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement) 751dd8e54a2SToby Isaac { 752dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 753dd8e54a2SToby Isaac 754dd8e54a2SToby Isaac PetscFunctionBegin; 755dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 756dd8e54a2SToby Isaac PetscValidIntPointer(minRefinement,2); 757dd8e54a2SToby Isaac *minRefinement = forest->minRefinement; 758dd8e54a2SToby Isaac PetscFunctionReturn(0); 759dd8e54a2SToby Isaac } 760dd8e54a2SToby Isaac 7619be51f97SToby Isaac /*@ 7629be51f97SToby Isaac DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base 7639be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. 7649be51f97SToby Isaac 7659be51f97SToby Isaac Logically collective on dm 7669be51f97SToby Isaac 7679be51f97SToby Isaac Input Parameters: 7689be51f97SToby Isaac + dm - the forest 7699be51f97SToby Isaac - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7709be51f97SToby Isaac 7719be51f97SToby Isaac Level: intermediate 7729be51f97SToby Isaac 7739be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM() 7749be51f97SToby Isaac @*/ 77556ba9f64SToby Isaac PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement) 77656ba9f64SToby Isaac { 77756ba9f64SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 77856ba9f64SToby Isaac 77956ba9f64SToby Isaac PetscFunctionBegin; 78056ba9f64SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 78156ba9f64SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup"); 78256ba9f64SToby Isaac forest->initRefinement = initRefinement; 78356ba9f64SToby Isaac PetscFunctionReturn(0); 78456ba9f64SToby Isaac } 78556ba9f64SToby Isaac 7869be51f97SToby Isaac /*@ 7879be51f97SToby Isaac DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see 7889be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. 7899be51f97SToby Isaac 7909be51f97SToby Isaac Not collective 7919be51f97SToby Isaac 7929be51f97SToby Isaac Input Parameter: 7939be51f97SToby Isaac . dm - the forest 7949be51f97SToby Isaac 7959be51f97SToby Isaac Output Paramater: 7969be51f97SToby Isaac . initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 7979be51f97SToby Isaac 7989be51f97SToby Isaac Level: intermediate 7999be51f97SToby Isaac 8009be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM() 8019be51f97SToby Isaac @*/ 80256ba9f64SToby Isaac PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement) 80356ba9f64SToby Isaac { 80456ba9f64SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 80556ba9f64SToby Isaac 80656ba9f64SToby Isaac PetscFunctionBegin; 80756ba9f64SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 80856ba9f64SToby Isaac PetscValidIntPointer(initRefinement,2); 80956ba9f64SToby Isaac *initRefinement = forest->initRefinement; 81056ba9f64SToby Isaac PetscFunctionReturn(0); 81156ba9f64SToby Isaac } 81256ba9f64SToby Isaac 8139be51f97SToby Isaac /*@ 8149be51f97SToby Isaac DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base 8159be51f97SToby Isaac DM, see DMForestGetBaseDM()) allowed in the forest. If the forest is being created by refining a previous forest 8169be51f97SToby Isaac (see DMForestGetAdaptivityForest()), this limits the amount of refinement. 8179be51f97SToby Isaac 8189be51f97SToby Isaac Logically collective on dm 8199be51f97SToby Isaac 8209be51f97SToby Isaac Input Parameters: 8219be51f97SToby Isaac + dm - the forest 8229be51f97SToby Isaac - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 8239be51f97SToby Isaac 8249be51f97SToby Isaac Level: intermediate 8259be51f97SToby Isaac 8269be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityDM() 8279be51f97SToby Isaac @*/ 828c7eeac06SToby Isaac PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement) 829dd8e54a2SToby Isaac { 830dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 831dd8e54a2SToby Isaac 832dd8e54a2SToby Isaac PetscFunctionBegin; 833dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 834ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup"); 835c7eeac06SToby Isaac forest->maxRefinement = maxRefinement; 836dd8e54a2SToby Isaac PetscFunctionReturn(0); 837dd8e54a2SToby Isaac } 838dd8e54a2SToby Isaac 8399be51f97SToby Isaac /*@ 8409be51f97SToby Isaac DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see 8419be51f97SToby Isaac DMForestGetBaseDM()) allowed in the forest. If the forest is being created by refining a previous forest (see 8429be51f97SToby Isaac DMForestGetAdaptivityForest()), this limits the amount of refinement. 8439be51f97SToby Isaac 8449be51f97SToby Isaac Not collective 8459be51f97SToby Isaac 8469be51f97SToby Isaac Input Parameter: 8479be51f97SToby Isaac . dm - the forest 8489be51f97SToby Isaac 8499be51f97SToby Isaac Output Parameter: 8509be51f97SToby Isaac . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest) 8519be51f97SToby Isaac 8529be51f97SToby Isaac Level: intermediate 8539be51f97SToby Isaac 8549be51f97SToby Isaac .seealso: DMForestSetMaximumRefinement(), DMForestGetMinimumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest() 8559be51f97SToby Isaac @*/ 856c7eeac06SToby Isaac PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement) 857dd8e54a2SToby Isaac { 858dd8e54a2SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 859dd8e54a2SToby Isaac 860dd8e54a2SToby Isaac PetscFunctionBegin; 861dd8e54a2SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 862c7eeac06SToby Isaac PetscValidIntPointer(maxRefinement,2); 863c7eeac06SToby Isaac *maxRefinement = forest->maxRefinement; 864dd8e54a2SToby Isaac PetscFunctionReturn(0); 865dd8e54a2SToby Isaac } 866c7eeac06SToby Isaac 8679be51f97SToby Isaac /*@C 8689be51f97SToby Isaac DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes. 8699be51f97SToby Isaac Subtypes of DMForest may define their own strategies. Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree 8709be51f97SToby Isaac for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to 8719be51f97SToby Isaac specify refinement/coarsening. 8729be51f97SToby Isaac 8739be51f97SToby Isaac Logically collective on dm 8749be51f97SToby Isaac 8759be51f97SToby Isaac Input Parameters: 8769be51f97SToby Isaac + dm - the forest 8779be51f97SToby Isaac - adaptStrategy - default DMFORESTADAPTALL 8789be51f97SToby Isaac 8799be51f97SToby Isaac Level: advanced 8809be51f97SToby Isaac 8819be51f97SToby Isaac .seealso: DMForestGetAdaptivityStrategy() 8829be51f97SToby Isaac @*/ 883c7eeac06SToby Isaac PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy) 884c7eeac06SToby Isaac { 885c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 886c7eeac06SToby Isaac PetscErrorCode ierr; 887c7eeac06SToby Isaac 888c7eeac06SToby Isaac PetscFunctionBegin; 889c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 890c7eeac06SToby Isaac ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr); 891a73e2921SToby Isaac ierr = PetscStrallocpy((const char*) adaptStrategy,(char**)&forest->adaptStrategy);CHKERRQ(ierr); 892c7eeac06SToby Isaac PetscFunctionReturn(0); 893c7eeac06SToby Isaac } 894c7eeac06SToby Isaac 8959be51f97SToby Isaac /*@C 8969be51f97SToby Isaac DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes. Subtypes 8979be51f97SToby Isaac of DMForest may define their own strategies. Two default strategies are DMFORESTADAPTALL, which indicates that all 8989be51f97SToby Isaac processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only 8999be51f97SToby Isaac one process needs to specify refinement/coarsening. 9009be51f97SToby Isaac 9019be51f97SToby Isaac Not collective 9029be51f97SToby Isaac 9039be51f97SToby Isaac Input Parameter: 9049be51f97SToby Isaac . dm - the forest 9059be51f97SToby Isaac 9069be51f97SToby Isaac Output Parameter: 9079be51f97SToby Isaac . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL) 9089be51f97SToby Isaac 9099be51f97SToby Isaac Level: advanced 9109be51f97SToby Isaac 9119be51f97SToby Isaac .seealso: DMForestSetAdaptivityStrategy() 9129be51f97SToby Isaac @*/ 913c7eeac06SToby Isaac PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy) 914c7eeac06SToby Isaac { 915c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 916c7eeac06SToby Isaac 917c7eeac06SToby Isaac PetscFunctionBegin; 918c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 919c7eeac06SToby Isaac PetscValidPointer(adaptStrategy,2); 920c7eeac06SToby Isaac *adaptStrategy = forest->adaptStrategy; 921c7eeac06SToby Isaac PetscFunctionReturn(0); 922c7eeac06SToby Isaac } 923c7eeac06SToby Isaac 9242a133e43SToby Isaac /*@ 9252a133e43SToby Isaac DMForestGetAdaptivitySuccess - Return whether the requested adaptation (refinement, coarsening, repartitioning, 9262a133e43SToby Isaac etc.) was successful. PETSC_FALSE indicates that the post-adaptation forest is the same as the pre-adpatation 9272a133e43SToby Isaac forest. A requested adaptation may have been unsuccessful if, for example, the requested refinement would have 9282a133e43SToby Isaac exceeded the maximum refinement level. 9292a133e43SToby Isaac 9302a133e43SToby Isaac Collective on dm 9312a133e43SToby Isaac 9322a133e43SToby Isaac Input Parameter: 9332a133e43SToby Isaac 9342a133e43SToby Isaac . dm - the post-adaptation forest 9352a133e43SToby Isaac 9362a133e43SToby Isaac Output Parameter: 9372a133e43SToby Isaac 9382a133e43SToby Isaac . success - PETSC_TRUE if the post-adaptation forest is different from the pre-adaptation forest. 9392a133e43SToby Isaac 9402a133e43SToby Isaac Level: intermediate 9412a133e43SToby Isaac 9422a133e43SToby Isaac .see 9432a133e43SToby Isaac @*/ 9442a133e43SToby Isaac PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success) 9452a133e43SToby Isaac { 9462a133e43SToby Isaac DM_Forest *forest; 9472a133e43SToby Isaac PetscErrorCode ierr; 9482a133e43SToby Isaac 9492a133e43SToby Isaac PetscFunctionBegin; 9502a133e43SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 9512a133e43SToby Isaac if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() has not been called yet."); 9522a133e43SToby Isaac forest = (DM_Forest *) dm->data; 9532a133e43SToby Isaac ierr = (forest->getadaptivitysuccess)(dm,success);CHKERRQ(ierr); 9542a133e43SToby Isaac PetscFunctionReturn(0); 9552a133e43SToby Isaac } 9562a133e43SToby Isaac 957bf9b5d84SToby Isaac /*@ 958bf9b5d84SToby Isaac DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed 959bf9b5d84SToby 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(). 960bf9b5d84SToby Isaac 961bf9b5d84SToby Isaac Logically collective on dm 962bf9b5d84SToby Isaac 963bf9b5d84SToby Isaac Input Parameters: 964bf9b5d84SToby Isaac + dm - the post-adaptation forest 965bf9b5d84SToby Isaac - computeSF - default PETSC_TRUE 966bf9b5d84SToby Isaac 967bf9b5d84SToby Isaac Level: advanced 968bf9b5d84SToby Isaac 969bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF() 970bf9b5d84SToby Isaac @*/ 971bf9b5d84SToby Isaac PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF) 972bf9b5d84SToby Isaac { 973bf9b5d84SToby Isaac DM_Forest *forest; 974bf9b5d84SToby Isaac 975bf9b5d84SToby Isaac PetscFunctionBegin; 976bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 977bf9b5d84SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called"); 978bf9b5d84SToby Isaac forest = (DM_Forest*) dm->data; 979bf9b5d84SToby Isaac forest->computeAdaptSF = computeSF; 980bf9b5d84SToby Isaac PetscFunctionReturn(0); 981bf9b5d84SToby Isaac } 982bf9b5d84SToby Isaac 9830eb7e1eaSToby Isaac PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time) 98480b27e07SToby Isaac { 98580b27e07SToby Isaac DM_Forest *forest; 98680b27e07SToby Isaac PetscErrorCode ierr; 98780b27e07SToby Isaac 98880b27e07SToby Isaac PetscFunctionBegin; 98980b27e07SToby Isaac PetscValidHeaderSpecific(dmIn ,DM_CLASSID ,1); 99080b27e07SToby Isaac PetscValidHeaderSpecific(vecIn ,VEC_CLASSID ,2); 99180b27e07SToby Isaac PetscValidHeaderSpecific(dmOut ,DM_CLASSID ,3); 99280b27e07SToby Isaac PetscValidHeaderSpecific(vecOut ,VEC_CLASSID ,4); 99380b27e07SToby Isaac forest = (DM_Forest *) dmIn->data; 99480b27e07SToby Isaac if (!forest->transfervec) SETERRQ(PetscObjectComm((PetscObject)dmIn),PETSC_ERR_SUP,"DMForestTransferVec() not implemented"); 9950eb7e1eaSToby Isaac ierr = (forest->transfervec)(dmIn,vecIn,dmOut,vecOut,useBCs,time);CHKERRQ(ierr); 99680b27e07SToby Isaac PetscFunctionReturn(0); 99780b27e07SToby Isaac } 99880b27e07SToby Isaac 999bf9b5d84SToby Isaac /*@ 1000bf9b5d84SToby Isaac DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the 1001bf9b5d84SToby Isaac pre-adaptation forest to the post-adaptiation forest. After DMSetUp() is called, these transfer PetscSFs can be 1002bf9b5d84SToby Isaac accessed with DMForestGetAdaptivitySF(). 1003bf9b5d84SToby Isaac 1004bf9b5d84SToby Isaac Not collective 1005bf9b5d84SToby Isaac 1006bf9b5d84SToby Isaac Input Parameter: 1007bf9b5d84SToby Isaac . dm - the post-adaptation forest 1008bf9b5d84SToby Isaac 1009bf9b5d84SToby Isaac Output Parameter: 1010bf9b5d84SToby Isaac . computeSF - default PETSC_TRUE 1011bf9b5d84SToby Isaac 1012bf9b5d84SToby Isaac Level: advanced 1013bf9b5d84SToby Isaac 1014bf9b5d84SToby Isaac .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF() 1015bf9b5d84SToby Isaac @*/ 1016bf9b5d84SToby Isaac PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF) 1017bf9b5d84SToby Isaac { 1018bf9b5d84SToby Isaac DM_Forest *forest; 1019bf9b5d84SToby Isaac 1020bf9b5d84SToby Isaac PetscFunctionBegin; 1021bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1022bf9b5d84SToby Isaac forest = (DM_Forest*) dm->data; 1023bf9b5d84SToby Isaac *computeSF = forest->computeAdaptSF; 1024bf9b5d84SToby Isaac PetscFunctionReturn(0); 1025bf9b5d84SToby Isaac } 1026bf9b5d84SToby Isaac 1027bf9b5d84SToby Isaac /*@ 1028bf9b5d84SToby Isaac DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest. 1029bf9b5d84SToby Isaac Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be 1030bf9b5d84SToby Isaac some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa. Therefore there are two 1031bf9b5d84SToby Isaac PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates 1032bf9b5d84SToby Isaac pre-adaptation fine cells to post-adaptation coarse cells. 1033bf9b5d84SToby Isaac 1034bf9b5d84SToby Isaac Not collective 1035bf9b5d84SToby Isaac 1036bf9b5d84SToby Isaac Input Parameter: 1037bf9b5d84SToby Isaac dm - the post-adaptation forest 1038bf9b5d84SToby Isaac 1039bf9b5d84SToby Isaac Output Parameter: 10400f17b9e3SToby Isaac preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post- 10410f17b9e3SToby Isaac coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre- 1042bf9b5d84SToby Isaac 1043bf9b5d84SToby Isaac Level: advanced 1044bf9b5d84SToby Isaac 1045bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF() 1046bf9b5d84SToby Isaac @*/ 10470f17b9e3SToby Isaac PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine) 1048bf9b5d84SToby Isaac { 1049bf9b5d84SToby Isaac DM_Forest *forest; 1050bf9b5d84SToby Isaac PetscErrorCode ierr; 1051bf9b5d84SToby Isaac 1052bf9b5d84SToby Isaac PetscFunctionBegin; 1053bf9b5d84SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1054bf9b5d84SToby Isaac ierr = DMSetUp(dm);CHKERRQ(ierr); 1055bf9b5d84SToby Isaac forest = (DM_Forest*) dm->data; 1056f885a11aSToby Isaac if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine; 1057f885a11aSToby Isaac if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine; 1058bf9b5d84SToby Isaac PetscFunctionReturn(0); 1059bf9b5d84SToby Isaac } 1060bf9b5d84SToby Isaac 10619be51f97SToby Isaac /*@ 10629be51f97SToby Isaac DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to 10639be51f97SToby Isaac indicate that the diameter of neighboring cells should differ by at most a factor of 2. Subtypes of DMForest may 10649be51f97SToby Isaac only support one particular choice of grading factor. 10659be51f97SToby Isaac 10669be51f97SToby Isaac Logically collective on dm 10679be51f97SToby Isaac 10689be51f97SToby Isaac Input Parameters: 10699be51f97SToby Isaac + dm - the forest 10709be51f97SToby Isaac - grade - the grading factor 10719be51f97SToby Isaac 10729be51f97SToby Isaac Level: advanced 10739be51f97SToby Isaac 10749be51f97SToby Isaac .seealso: DMForestGetGradeFactor() 10759be51f97SToby Isaac @*/ 1076c7eeac06SToby Isaac PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade) 1077c7eeac06SToby Isaac { 1078c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1079c7eeac06SToby Isaac 1080c7eeac06SToby Isaac PetscFunctionBegin; 1081c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1082ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup"); 1083c7eeac06SToby Isaac forest->gradeFactor = grade; 1084c7eeac06SToby Isaac PetscFunctionReturn(0); 1085c7eeac06SToby Isaac } 1086c7eeac06SToby Isaac 10879be51f97SToby Isaac /*@ 10889be51f97SToby Isaac DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of 10899be51f97SToby Isaac neighboring cells should differ by at most a factor of 2. Subtypes of DMForest may only support one particular 10909be51f97SToby Isaac choice of grading factor. 10919be51f97SToby Isaac 10929be51f97SToby Isaac Not collective 10939be51f97SToby Isaac 10949be51f97SToby Isaac Input Parameter: 10959be51f97SToby Isaac . dm - the forest 10969be51f97SToby Isaac 10979be51f97SToby Isaac Output Parameter: 10989be51f97SToby Isaac . grade - the grading factor 10999be51f97SToby Isaac 11009be51f97SToby Isaac Level: advanced 11019be51f97SToby Isaac 11029be51f97SToby Isaac .seealso: DMForestSetGradeFactor() 11039be51f97SToby Isaac @*/ 1104c7eeac06SToby Isaac PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade) 1105c7eeac06SToby Isaac { 1106c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1107c7eeac06SToby Isaac 1108c7eeac06SToby Isaac PetscFunctionBegin; 1109c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1110c7eeac06SToby Isaac PetscValidIntPointer(grade,2); 1111c7eeac06SToby Isaac *grade = forest->gradeFactor; 1112c7eeac06SToby Isaac PetscFunctionReturn(0); 1113c7eeac06SToby Isaac } 1114c7eeac06SToby Isaac 11159be51f97SToby Isaac /*@ 11169be51f97SToby Isaac DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes 11179be51f97SToby Isaac the cell weight (see DMForestSetCellWeights()) when calculating partitions. The final weight of a cell will be 11189be51f97SToby Isaac (cellWeight) * (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on 11199be51f97SToby Isaac its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the 11209be51f97SToby Isaac computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement. 11219be51f97SToby Isaac 11229be51f97SToby Isaac Logically collective on dm 11239be51f97SToby Isaac 11249be51f97SToby Isaac Input Parameters: 11259be51f97SToby Isaac + dm - the forest 11269be51f97SToby Isaac - weightsFactors - default 1. 11279be51f97SToby Isaac 11289be51f97SToby Isaac Level: advanced 11299be51f97SToby Isaac 11309be51f97SToby Isaac .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights() 11319be51f97SToby Isaac @*/ 1132ef51cf95SToby Isaac PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor) 1133c7eeac06SToby Isaac { 1134c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1135c7eeac06SToby Isaac 1136c7eeac06SToby Isaac PetscFunctionBegin; 1137c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1138ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup"); 1139c7eeac06SToby Isaac forest->weightsFactor = weightsFactor; 1140c7eeac06SToby Isaac PetscFunctionReturn(0); 1141c7eeac06SToby Isaac } 1142c7eeac06SToby Isaac 11439be51f97SToby Isaac /*@ 11449be51f97SToby Isaac DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see 11459be51f97SToby Isaac DMForestSetCellWeights()) when calculating partitions. The final weight of a cell will be (cellWeight) * 11469be51f97SToby Isaac (weightFactor^refinementLevel). A factor of 1 indicates that the weight of a cell does not depend on its level; a 11479be51f97SToby Isaac factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation 11489be51f97SToby Isaac associated with a cell is multiplied by a factor of 2 for each additional level of refinement. 11499be51f97SToby Isaac 11509be51f97SToby Isaac Not collective 11519be51f97SToby Isaac 11529be51f97SToby Isaac Input Parameter: 11539be51f97SToby Isaac . dm - the forest 11549be51f97SToby Isaac 11559be51f97SToby Isaac Output Parameter: 11569be51f97SToby Isaac . weightsFactors - default 1. 11579be51f97SToby Isaac 11589be51f97SToby Isaac Level: advanced 11599be51f97SToby Isaac 11609be51f97SToby Isaac .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights() 11619be51f97SToby Isaac @*/ 1162ef51cf95SToby Isaac PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor) 1163c7eeac06SToby Isaac { 1164c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1165c7eeac06SToby Isaac 1166c7eeac06SToby Isaac PetscFunctionBegin; 1167c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1168c7eeac06SToby Isaac PetscValidRealPointer(weightsFactor,2); 1169c7eeac06SToby Isaac *weightsFactor = forest->weightsFactor; 1170c7eeac06SToby Isaac PetscFunctionReturn(0); 1171c7eeac06SToby Isaac } 1172c7eeac06SToby Isaac 11739be51f97SToby Isaac /*@ 11749be51f97SToby Isaac DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process 11759be51f97SToby Isaac 11769be51f97SToby Isaac Not collective 11779be51f97SToby Isaac 11789be51f97SToby Isaac Input Parameter: 11799be51f97SToby Isaac . dm - the forest 11809be51f97SToby Isaac 11819be51f97SToby Isaac Output Parameters: 11829be51f97SToby Isaac + cStart - the first cell on this process 11839be51f97SToby Isaac - cEnd - one after the final cell on this process 11849be51f97SToby Isaac 11851a244344SSatish Balay Level: intermediate 11869be51f97SToby Isaac 11879be51f97SToby Isaac .seealso: DMForestGetCellSF() 11889be51f97SToby Isaac @*/ 1189c7eeac06SToby Isaac PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd) 1190c7eeac06SToby Isaac { 1191c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1192c7eeac06SToby Isaac PetscErrorCode ierr; 1193c7eeac06SToby Isaac 1194c7eeac06SToby Isaac PetscFunctionBegin; 1195c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1196c7eeac06SToby Isaac PetscValidIntPointer(cStart,2); 1197c7eeac06SToby Isaac PetscValidIntPointer(cEnd,2); 1198c7eeac06SToby Isaac if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) { 1199c7eeac06SToby Isaac ierr = forest->createcellchart(dm,&forest->cStart,&forest->cEnd);CHKERRQ(ierr); 1200c7eeac06SToby Isaac } 1201c7eeac06SToby Isaac *cStart = forest->cStart; 1202c7eeac06SToby Isaac *cEnd = forest->cEnd; 1203c7eeac06SToby Isaac PetscFunctionReturn(0); 1204c7eeac06SToby Isaac } 1205c7eeac06SToby Isaac 12069be51f97SToby Isaac /*@ 12079be51f97SToby Isaac DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes 12089be51f97SToby Isaac 12099be51f97SToby Isaac Not collective 12109be51f97SToby Isaac 12119be51f97SToby Isaac Input Parameter: 12129be51f97SToby Isaac . dm - the forest 12139be51f97SToby Isaac 12149be51f97SToby Isaac Output Parameter: 12159be51f97SToby Isaac . cellSF - the PetscSF 12169be51f97SToby Isaac 12171a244344SSatish Balay Level: intermediate 12189be51f97SToby Isaac 12199be51f97SToby Isaac .seealso: DMForestGetCellChart() 12209be51f97SToby Isaac @*/ 1221c7eeac06SToby Isaac PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF) 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 PetscValidPointer(cellSF,2); 1229c7eeac06SToby Isaac if ((!forest->cellSF) && forest->createcellsf) { 1230c7eeac06SToby Isaac ierr = forest->createcellsf(dm,&forest->cellSF);CHKERRQ(ierr); 1231c7eeac06SToby Isaac } 1232c7eeac06SToby Isaac *cellSF = forest->cellSF; 1233c7eeac06SToby Isaac PetscFunctionReturn(0); 1234c7eeac06SToby Isaac } 1235c7eeac06SToby Isaac 12369be51f97SToby Isaac /*@C 12379be51f97SToby Isaac DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see 12389be51f97SToby Isaac DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination). The 1239cd3c525cSToby Isaac interpretation of the label values is up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, 1240cd3c525cSToby Isaac DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have been reserved as choices that should be accepted by all subtypes. 12419be51f97SToby Isaac 12429be51f97SToby Isaac Logically collective on dm 12439be51f97SToby Isaac 12449be51f97SToby Isaac Input Parameters: 12459be51f97SToby Isaac - dm - the forest 1246a1b0c543SToby Isaac + adaptLabel - the label in the pre-adaptation forest 12479be51f97SToby Isaac 12489be51f97SToby Isaac Level: intermediate 12499be51f97SToby Isaac 12509be51f97SToby Isaac .seealso DMForestGetAdaptivityLabel() 12519be51f97SToby Isaac @*/ 1252a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityLabel(DM dm, DMLabel adaptLabel) 1253c7eeac06SToby Isaac { 1254c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1255c7eeac06SToby Isaac PetscErrorCode ierr; 1256c7eeac06SToby Isaac 1257c7eeac06SToby Isaac PetscFunctionBegin; 1258c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1259a1b0c543SToby Isaac adaptLabel->refct++; 1260a1b0c543SToby Isaac if (forest->adaptLabel) {ierr = DMLabelDestroy(&forest->adaptLabel);CHKERRQ(ierr);} 1261a1b0c543SToby Isaac forest->adaptLabel = adaptLabel; 1262c7eeac06SToby Isaac PetscFunctionReturn(0); 1263c7eeac06SToby Isaac } 1264c7eeac06SToby Isaac 12659be51f97SToby Isaac /*@C 12669be51f97SToby Isaac DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that 12679be51f97SToby Isaac holds the adaptation flags (refinement, coarsening, or some combination). The interpretation of the label values is 1268cd3c525cSToby Isaac up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have 1269cd3c525cSToby Isaac been reserved as choices that should be accepted by all subtypes. 12709be51f97SToby Isaac 12719be51f97SToby Isaac Not collective 12729be51f97SToby Isaac 12739be51f97SToby Isaac Input Parameter: 12749be51f97SToby Isaac . dm - the forest 12759be51f97SToby Isaac 12769be51f97SToby Isaac Output Parameter: 12779be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest 12789be51f97SToby Isaac 12799be51f97SToby Isaac Level: intermediate 12809be51f97SToby Isaac 12819be51f97SToby Isaac .seealso DMForestSetAdaptivityLabel() 12829be51f97SToby Isaac @*/ 1283a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityLabel(DM dm, DMLabel *adaptLabel) 1284c7eeac06SToby Isaac { 1285c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1286c7eeac06SToby Isaac 1287c7eeac06SToby Isaac PetscFunctionBegin; 1288c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1289ba936b91SToby Isaac *adaptLabel = forest->adaptLabel; 1290c7eeac06SToby Isaac PetscFunctionReturn(0); 1291c7eeac06SToby Isaac } 1292c7eeac06SToby Isaac 12939be51f97SToby Isaac /*@ 12949be51f97SToby Isaac DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current 12959be51f97SToby Isaac process: weights are used to determine parallel partitioning. Partitions will be created so that each process's 12969be51f97SToby Isaac ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight 12979be51f97SToby Isaac of 1. 12989be51f97SToby Isaac 12999be51f97SToby Isaac Logically collective on dm 13009be51f97SToby Isaac 13019be51f97SToby Isaac Input Parameters: 13029be51f97SToby Isaac + dm - the forest 13039be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1. 13049be51f97SToby Isaac - copyMode - how weights should reference weights 13059be51f97SToby Isaac 13069be51f97SToby Isaac Level: advanced 13079be51f97SToby Isaac 13089be51f97SToby Isaac .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity() 13099be51f97SToby Isaac @*/ 1310c7eeac06SToby Isaac PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode) 1311c7eeac06SToby Isaac { 1312c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1313c7eeac06SToby Isaac PetscInt cStart, cEnd; 1314c7eeac06SToby Isaac PetscErrorCode ierr; 1315c7eeac06SToby Isaac 1316c7eeac06SToby Isaac PetscFunctionBegin; 1317c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1318c7eeac06SToby Isaac ierr = DMForestGetCellChart(dm,&cStart,&cEnd);CHKERRQ(ierr); 1319c7eeac06SToby Isaac if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd); 1320c7eeac06SToby Isaac if (copyMode == PETSC_COPY_VALUES) { 1321c7eeac06SToby Isaac if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) { 1322c7eeac06SToby Isaac ierr = PetscMalloc1(cEnd-cStart,&forest->cellWeights);CHKERRQ(ierr); 1323c7eeac06SToby Isaac } 1324c7eeac06SToby Isaac ierr = PetscMemcpy(forest->cellWeights,weights,(cEnd-cStart)*sizeof(*weights));CHKERRQ(ierr); 1325c7eeac06SToby Isaac forest->cellWeightsCopyMode = PETSC_OWN_POINTER; 1326c7eeac06SToby Isaac PetscFunctionReturn(0); 1327c7eeac06SToby Isaac } 1328c7eeac06SToby Isaac if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) { 1329c7eeac06SToby Isaac ierr = PetscFree(forest->cellWeights);CHKERRQ(ierr); 1330c7eeac06SToby Isaac } 1331c7eeac06SToby Isaac forest->cellWeights = weights; 1332c7eeac06SToby Isaac forest->cellWeightsCopyMode = copyMode; 1333c7eeac06SToby Isaac PetscFunctionReturn(0); 1334c7eeac06SToby Isaac } 1335c7eeac06SToby Isaac 13369be51f97SToby Isaac /*@ 13379be51f97SToby Isaac DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current 13389be51f97SToby Isaac process: weights are used to determine parallel partitioning. Partitions will be created so that each process's 13399be51f97SToby Isaac ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight 13409be51f97SToby Isaac of 1. 13419be51f97SToby Isaac 13429be51f97SToby Isaac Not collective 13439be51f97SToby Isaac 13449be51f97SToby Isaac Input Parameter: 13459be51f97SToby Isaac . dm - the forest 13469be51f97SToby Isaac 13479be51f97SToby Isaac Output Parameter: 13489be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1. 13499be51f97SToby Isaac 13509be51f97SToby Isaac Level: advanced 13519be51f97SToby Isaac 13529be51f97SToby Isaac .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity() 13539be51f97SToby Isaac @*/ 1354c7eeac06SToby Isaac PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights) 1355c7eeac06SToby Isaac { 1356c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1357c7eeac06SToby Isaac 1358c7eeac06SToby Isaac PetscFunctionBegin; 1359c7eeac06SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1360c7eeac06SToby Isaac PetscValidPointer(weights,2); 1361c7eeac06SToby Isaac *weights = forest->cellWeights; 1362c7eeac06SToby Isaac PetscFunctionReturn(0); 1363c7eeac06SToby Isaac } 1364c7eeac06SToby Isaac 13659be51f97SToby Isaac /*@ 13669be51f97SToby Isaac DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning 13679be51f97SToby Isaac a pre-adaptation forest (see DMForestGetAdaptivityForest()). After partitioning, the ratio of the weight of each 13689be51f97SToby Isaac process's cells to the process's capacity will be roughly equal for all processes. A capacity of 0 indicates that 13699be51f97SToby Isaac the current process should not have any cells after repartitioning. 13709be51f97SToby Isaac 13719be51f97SToby Isaac Logically Collective on dm 13729be51f97SToby Isaac 13739be51f97SToby Isaac Input parameters: 13749be51f97SToby Isaac + dm - the forest 13759be51f97SToby Isaac - capacity - this process's capacity 13769be51f97SToby Isaac 13779be51f97SToby Isaac Level: advanced 13789be51f97SToby Isaac 13799be51f97SToby Isaac .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor() 13809be51f97SToby Isaac @*/ 1381c7eeac06SToby Isaac PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity) 1382c7eeac06SToby Isaac { 1383c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1384c7eeac06SToby Isaac 1385c7eeac06SToby Isaac PetscFunctionBegin; 1386c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1387ef51cf95SToby Isaac if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup"); 1388c7eeac06SToby Isaac if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity); 1389c7eeac06SToby Isaac forest->weightCapacity = capacity; 1390c7eeac06SToby Isaac PetscFunctionReturn(0); 1391c7eeac06SToby Isaac } 1392c7eeac06SToby Isaac 13939be51f97SToby Isaac /*@ 13949be51f97SToby Isaac DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see 13959be51f97SToby Isaac DMForestGetAdaptivityForest()). After partitioning, the ratio of the weight of each process's cells to the 13969be51f97SToby Isaac process's capacity will be roughly equal for all processes. A capacity of 0 indicates that the current process 13979be51f97SToby Isaac should not have any cells after repartitioning. 13989be51f97SToby Isaac 13999be51f97SToby Isaac Not collective 14009be51f97SToby Isaac 14019be51f97SToby Isaac Input parameter: 14029be51f97SToby Isaac . dm - the forest 14039be51f97SToby Isaac 14049be51f97SToby Isaac Output parameter: 14059be51f97SToby Isaac . capacity - this process's capacity 14069be51f97SToby Isaac 14079be51f97SToby Isaac Level: advanced 14089be51f97SToby Isaac 14099be51f97SToby Isaac .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor() 14109be51f97SToby Isaac @*/ 1411c7eeac06SToby Isaac PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity) 1412c7eeac06SToby Isaac { 1413c7eeac06SToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 1414c7eeac06SToby Isaac 1415c7eeac06SToby Isaac PetscFunctionBegin; 1416c7eeac06SToby Isaac PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1417c7eeac06SToby Isaac PetscValidRealPointer(capacity,2); 1418c7eeac06SToby Isaac *capacity = forest->weightCapacity; 1419c7eeac06SToby Isaac PetscFunctionReturn(0); 1420c7eeac06SToby Isaac } 1421c7eeac06SToby Isaac 14220709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm) 1423db4d5e8cSToby Isaac { 1424db4d5e8cSToby Isaac DM_Forest *forest = (DM_Forest*) dm->data; 142556ba9f64SToby Isaac PetscBool flg, flg1, flg2, flg3, flg4; 1426dd8e54a2SToby Isaac DMForestTopology oldTopo; 1427c7eeac06SToby Isaac char stringBuffer[256]; 1428dd8e54a2SToby Isaac PetscViewer viewer; 1429dd8e54a2SToby Isaac PetscViewerFormat format; 143056ba9f64SToby Isaac PetscInt adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade; 1431c7eeac06SToby Isaac PetscReal weightsFactor; 1432c7eeac06SToby Isaac DMForestAdaptivityStrategy adaptStrategy; 1433db4d5e8cSToby Isaac PetscErrorCode ierr; 1434db4d5e8cSToby Isaac 1435db4d5e8cSToby Isaac PetscFunctionBegin; 1436db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 143758762b62SToby Isaac forest->setfromoptionscalled = PETSC_TRUE; 1438dd8e54a2SToby Isaac ierr = DMForestGetTopology(dm, &oldTopo);CHKERRQ(ierr); 1439a3eda1eaSToby Isaac ierr = PetscOptionsHead(PetscOptionsObject,"DMForest Options");CHKERRQ(ierr); 144056ba9f64SToby Isaac ierr = PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);CHKERRQ(ierr); 144156ba9f64SToby Isaac ierr = PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);CHKERRQ(ierr); 144256ba9f64SToby Isaac ierr = PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);CHKERRQ(ierr); 144356ba9f64SToby Isaac ierr = PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);CHKERRQ(ierr); 1444f885a11aSToby 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}"); 144556ba9f64SToby Isaac if (flg1) { 144656ba9f64SToby Isaac ierr = DMForestSetTopology(dm,(DMForestTopology)stringBuffer);CHKERRQ(ierr); 144756ba9f64SToby Isaac ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr); 144820e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr); 144956ba9f64SToby Isaac } 145056ba9f64SToby Isaac if (flg2) { 1451dd8e54a2SToby Isaac DM base; 1452dd8e54a2SToby Isaac 1453dd8e54a2SToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),&base);CHKERRQ(ierr); 1454dd8e54a2SToby Isaac ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 1455dd8e54a2SToby Isaac ierr = DMLoad(base,viewer);CHKERRQ(ierr); 1456dd8e54a2SToby Isaac ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1457dd8e54a2SToby Isaac ierr = DMForestSetBaseDM(dm,base);CHKERRQ(ierr); 1458dd8e54a2SToby Isaac ierr = DMDestroy(&base);CHKERRQ(ierr); 145956ba9f64SToby Isaac ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr); 146020e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr); 1461dd8e54a2SToby Isaac } 146256ba9f64SToby Isaac if (flg3) { 1463dd8e54a2SToby Isaac DM coarse; 1464dd8e54a2SToby Isaac 1465dd8e54a2SToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),&coarse);CHKERRQ(ierr); 1466dd8e54a2SToby Isaac ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 1467dd8e54a2SToby Isaac ierr = DMLoad(coarse,viewer);CHKERRQ(ierr); 1468dd8e54a2SToby Isaac ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 146920e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,coarse);CHKERRQ(ierr); 1470dd8e54a2SToby Isaac ierr = DMDestroy(&coarse);CHKERRQ(ierr); 147156ba9f64SToby Isaac ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr); 147256ba9f64SToby Isaac ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr); 1473dd8e54a2SToby Isaac } 147456ba9f64SToby Isaac if (flg4) { 1475dd8e54a2SToby Isaac DM fine; 1476dd8e54a2SToby Isaac 1477dd8e54a2SToby Isaac ierr = DMCreate(PetscObjectComm((PetscObject)dm),&fine);CHKERRQ(ierr); 1478dd8e54a2SToby Isaac ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr); 1479dd8e54a2SToby Isaac ierr = DMLoad(fine,viewer);CHKERRQ(ierr); 1480dd8e54a2SToby Isaac ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 148120e8089bSToby Isaac ierr = DMForestSetAdaptivityForest(dm,fine);CHKERRQ(ierr); 1482dd8e54a2SToby Isaac ierr = DMDestroy(&fine);CHKERRQ(ierr); 148356ba9f64SToby Isaac ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr); 148456ba9f64SToby Isaac ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr); 1485dd8e54a2SToby Isaac } 1486dd8e54a2SToby Isaac ierr = DMForestGetAdjacencyDimension(dm,&adjDim);CHKERRQ(ierr); 1487dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg);CHKERRQ(ierr); 1488dd8e54a2SToby Isaac if (flg) { 1489dd8e54a2SToby Isaac ierr = DMForestSetAdjacencyDimension(dm,adjDim);CHKERRQ(ierr); 1490f885a11aSToby Isaac } else { 1491dd8e54a2SToby Isaac ierr = DMForestGetAdjacencyCodimension(dm,&adjCodim);CHKERRQ(ierr); 1492dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg);CHKERRQ(ierr); 1493dd8e54a2SToby Isaac if (flg) { 1494dd8e54a2SToby Isaac ierr = DMForestSetAdjacencyCodimension(dm,adjCodim);CHKERRQ(ierr); 1495dd8e54a2SToby Isaac } 1496dd8e54a2SToby Isaac } 1497dd8e54a2SToby Isaac ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr); 1498dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg);CHKERRQ(ierr); 1499dd8e54a2SToby Isaac if (flg) { 1500dd8e54a2SToby Isaac ierr = DMForestSetPartitionOverlap(dm,overlap);CHKERRQ(ierr); 1501dd8e54a2SToby Isaac } 1502a6121fbdSMatthew G. Knepley #if 0 1503a6121fbdSMatthew 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); 1504a6121fbdSMatthew G. Knepley if (flg) { 1505a6121fbdSMatthew G. Knepley ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr); 1506a6121fbdSMatthew G. Knepley ierr = DMForestSetInitialRefinement(dm,minRefinement);CHKERRQ(ierr); 1507a6121fbdSMatthew G. Knepley } 1508a6121fbdSMatthew 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); 1509a6121fbdSMatthew G. Knepley if (flg) { 1510a6121fbdSMatthew G. Knepley ierr = DMForestSetMinimumRefinement(dm,0);CHKERRQ(ierr); 1511a6121fbdSMatthew G. Knepley ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr); 1512a6121fbdSMatthew G. Knepley } 1513a6121fbdSMatthew G. Knepley #endif 1514dd8e54a2SToby Isaac ierr = DMForestGetMinimumRefinement(dm,&minRefinement);CHKERRQ(ierr); 1515dd8e54a2SToby Isaac ierr = PetscOptionsInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg);CHKERRQ(ierr); 1516dd8e54a2SToby Isaac if (flg) { 1517dd8e54a2SToby Isaac ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr); 1518db4d5e8cSToby Isaac } 151956ba9f64SToby Isaac ierr = DMForestGetInitialRefinement(dm,&initRefinement);CHKERRQ(ierr); 152056ba9f64SToby Isaac ierr = PetscOptionsInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg);CHKERRQ(ierr); 152156ba9f64SToby Isaac if (flg) { 152256ba9f64SToby Isaac ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr); 152356ba9f64SToby Isaac } 1524c7eeac06SToby Isaac ierr = DMForestGetMaximumRefinement(dm,&maxRefinement);CHKERRQ(ierr); 1525c7eeac06SToby Isaac ierr = PetscOptionsInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg);CHKERRQ(ierr); 1526c7eeac06SToby Isaac if (flg) { 1527c7eeac06SToby Isaac ierr = DMForestSetMaximumRefinement(dm,maxRefinement);CHKERRQ(ierr); 1528c7eeac06SToby Isaac } 1529c7eeac06SToby Isaac ierr = DMForestGetAdaptivityStrategy(dm,&adaptStrategy);CHKERRQ(ierr); 1530c7eeac06SToby Isaac ierr = PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,256,&flg);CHKERRQ(ierr); 1531c7eeac06SToby Isaac if (flg) { 1532c7eeac06SToby Isaac ierr = DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);CHKERRQ(ierr); 1533c7eeac06SToby Isaac } 1534c7eeac06SToby Isaac ierr = DMForestGetGradeFactor(dm,&grade);CHKERRQ(ierr); 1535c7eeac06SToby Isaac ierr = PetscOptionsInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg);CHKERRQ(ierr); 1536c7eeac06SToby Isaac if (flg) { 1537c7eeac06SToby Isaac ierr = DMForestSetGradeFactor(dm,grade);CHKERRQ(ierr); 1538c7eeac06SToby Isaac } 1539c7eeac06SToby Isaac ierr = DMForestGetCellWeightFactor(dm,&weightsFactor);CHKERRQ(ierr); 1540c7eeac06SToby Isaac ierr = PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);CHKERRQ(ierr); 1541c7eeac06SToby Isaac if (flg) { 1542c7eeac06SToby Isaac ierr = DMForestSetCellWeightFactor(dm,weightsFactor);CHKERRQ(ierr); 1543c7eeac06SToby Isaac } 1544db4d5e8cSToby Isaac ierr = PetscOptionsTail();CHKERRQ(ierr); 1545db4d5e8cSToby Isaac PetscFunctionReturn(0); 1546db4d5e8cSToby Isaac } 1547db4d5e8cSToby Isaac 1548276c5506SMatthew G. Knepley PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm) 1549d8984e3bSMatthew G. Knepley { 1550d8984e3bSMatthew G. Knepley PetscErrorCode ierr; 1551d8984e3bSMatthew G. Knepley 1552d8984e3bSMatthew G. Knepley PetscFunctionBegin; 1553d8984e3bSMatthew G. Knepley if (subdm) {ierr = DMClone(dm, subdm);CHKERRQ(ierr);} 1554d8984e3bSMatthew G. Knepley ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr); 1555d8984e3bSMatthew G. Knepley PetscFunctionReturn(0); 1556d8984e3bSMatthew G. Knepley } 1557d8984e3bSMatthew G. Knepley 15585421bac9SToby Isaac PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined) 15595421bac9SToby Isaac { 15605421bac9SToby Isaac DMLabel refine; 15615421bac9SToby Isaac DM fineDM; 15625421bac9SToby Isaac PetscErrorCode ierr; 15635421bac9SToby Isaac 15645421bac9SToby Isaac PetscFunctionBegin; 15655421bac9SToby Isaac ierr = DMGetFineDM(dm,&fineDM);CHKERRQ(ierr); 15665421bac9SToby Isaac if (fineDM) { 15675421bac9SToby Isaac ierr = PetscObjectReference((PetscObject)fineDM);CHKERRQ(ierr); 15685421bac9SToby Isaac *dmRefined = fineDM; 15695421bac9SToby Isaac PetscFunctionReturn(0); 15705421bac9SToby Isaac } 15715421bac9SToby Isaac ierr = DMForestTemplate(dm,comm,dmRefined);CHKERRQ(ierr); 15725421bac9SToby Isaac ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr); 15735421bac9SToby Isaac if (!refine) { 1574a1b0c543SToby Isaac ierr = DMLabelCreate("refine",&refine);CHKERRQ(ierr); 1575a1b0c543SToby Isaac ierr = DMLabelSetDefaultValue(refine,DM_ADAPT_REFINE);CHKERRQ(ierr); 15765421bac9SToby Isaac } 1577a1b0c543SToby Isaac else { 1578a1b0c543SToby Isaac refine->refct++; 1579a1b0c543SToby Isaac } 1580a1b0c543SToby Isaac ierr = DMForestSetAdaptivityLabel(*dmRefined,refine);CHKERRQ(ierr); 1581a1b0c543SToby Isaac ierr = DMLabelDestroy(&refine);CHKERRQ(ierr); 15825421bac9SToby Isaac PetscFunctionReturn(0); 15835421bac9SToby Isaac } 15845421bac9SToby Isaac 15855421bac9SToby Isaac PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened) 15865421bac9SToby Isaac { 15875421bac9SToby Isaac DMLabel coarsen; 15885421bac9SToby Isaac DM coarseDM; 15895421bac9SToby Isaac PetscErrorCode ierr; 15905421bac9SToby Isaac 15915421bac9SToby Isaac PetscFunctionBegin; 15924098eed7SToby Isaac { 15934098eed7SToby Isaac PetscMPIInt mpiComparison; 15944098eed7SToby Isaac MPI_Comm dmcomm = PetscObjectComm((PetscObject)dm); 15954098eed7SToby Isaac 15964098eed7SToby Isaac ierr = MPI_Comm_compare(comm, dmcomm, &mpiComparison);CHKERRQ(ierr); 1597f885a11aSToby Isaac if (mpiComparison != MPI_IDENT && mpiComparison != MPI_CONGRUENT) SETERRQ(dmcomm,PETSC_ERR_SUP,"No support for different communicators yet"); 15984098eed7SToby Isaac } 15995421bac9SToby Isaac ierr = DMGetCoarseDM(dm,&coarseDM);CHKERRQ(ierr); 16005421bac9SToby Isaac if (coarseDM) { 16015421bac9SToby Isaac ierr = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr); 16025421bac9SToby Isaac *dmCoarsened = coarseDM; 16035421bac9SToby Isaac PetscFunctionReturn(0); 16045421bac9SToby Isaac } 16055421bac9SToby Isaac ierr = DMForestTemplate(dm,comm,dmCoarsened);CHKERRQ(ierr); 16067ee90a76SStefano Zampini ierr = DMForestSetAdaptivityPurpose(*dmCoarsened,DM_ADAPT_COARSEN);CHKERRQ(ierr); 16075421bac9SToby Isaac ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr); 16085421bac9SToby Isaac if (!coarsen) { 1609a1b0c543SToby Isaac ierr = DMLabelCreate("coarsen",&coarsen);CHKERRQ(ierr); 1610a1b0c543SToby Isaac ierr = DMLabelSetDefaultValue(coarsen,DM_ADAPT_COARSEN);CHKERRQ(ierr); 1611a1b0c543SToby Isaac } else { 1612a1b0c543SToby Isaac coarsen->refct++; 16135421bac9SToby Isaac } 1614a1b0c543SToby Isaac ierr = DMForestSetAdaptivityLabel(*dmCoarsened,coarsen);CHKERRQ(ierr); 1615a1b0c543SToby Isaac ierr = DMLabelDestroy(&coarsen);CHKERRQ(ierr); 16165421bac9SToby Isaac PetscFunctionReturn(0); 16175421bac9SToby Isaac } 16185421bac9SToby Isaac 1619a1b0c543SToby Isaac static PetscErrorCode DMAdaptLabel_Forest(DM dm, DMLabel label, DM *adaptedDM) 162009350103SToby Isaac { 162109350103SToby Isaac PetscBool success; 162209350103SToby Isaac PetscErrorCode ierr; 162309350103SToby Isaac 162409350103SToby Isaac PetscFunctionBegin; 162509350103SToby Isaac ierr = DMForestTemplate(dm,PetscObjectComm((PetscObject)dm),adaptedDM);CHKERRQ(ierr); 1626a1b0c543SToby Isaac ierr = DMForestSetAdaptivityLabel(*adaptedDM,label);CHKERRQ(ierr); 162709350103SToby Isaac ierr = DMSetUp(*adaptedDM);CHKERRQ(ierr); 162809350103SToby Isaac ierr = DMForestGetAdaptivitySuccess(*adaptedDM,&success);CHKERRQ(ierr); 162909350103SToby Isaac if (!success) { 163009350103SToby Isaac ierr = DMDestroy(adaptedDM);CHKERRQ(ierr); 163109350103SToby Isaac *adaptedDM = NULL; 163209350103SToby Isaac } 163309350103SToby Isaac PetscFunctionReturn(0); 163409350103SToby Isaac } 163509350103SToby Isaac 1636d222f98bSToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm) 1637d222f98bSToby Isaac { 1638d222f98bSToby Isaac PetscErrorCode ierr; 1639d222f98bSToby Isaac 1640d222f98bSToby Isaac PetscFunctionBegin; 1641d222f98bSToby Isaac ierr = PetscMemzero(dm->ops,sizeof(*(dm->ops)));CHKERRQ(ierr); 1642d222f98bSToby Isaac 1643d222f98bSToby Isaac dm->ops->clone = DMClone_Forest; 1644d222f98bSToby Isaac dm->ops->setfromoptions = DMSetFromOptions_Forest; 1645d222f98bSToby Isaac dm->ops->destroy = DMDestroy_Forest; 1646d8984e3bSMatthew G. Knepley dm->ops->createsubdm = DMCreateSubDM_Forest; 16475421bac9SToby Isaac dm->ops->refine = DMRefine_Forest; 16485421bac9SToby Isaac dm->ops->coarsen = DMCoarsen_Forest; 16490d1cd5e0SMatthew G. Knepley dm->ops->adaptlabel = DMAdaptLabel_Forest; 1650d222f98bSToby Isaac PetscFunctionReturn(0); 1651d222f98bSToby Isaac } 1652d222f98bSToby Isaac 16539be51f97SToby Isaac /*MC 16549be51f97SToby Isaac 1655bae1f979SBarry Smith DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh. Forests usually have a base DM 1656bae1f979SBarry Smith (see DMForestGetBaseDM()), from which it is refined. The refinement and partitioning of forests is considered 1657bae1f979SBarry Smith immutable after DMSetUp() is called. To adapt a mesh, one should call DMForestTemplate() to create a new mesh that 1658bae1f979SBarry Smith will default to being identical to it, specify how that mesh should differ, and then calling DMSetUp() on the new 1659bae1f979SBarry Smith mesh. 1660bae1f979SBarry Smith 1661bae1f979SBarry Smith To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the 1662bae1f979SBarry Smith previous mesh whose values indicate which cells should be refined (DM_ADAPT_REFINE) or coarsened (DM_ADAPT_COARSEN) 1663bae1f979SBarry Smith and how (subtypes are free to allow additional values for things like anisotropic refinement). The label should be 1664bae1f979SBarry Smith given to the *new* mesh with DMForestSetAdaptivityLabel(). 16659be51f97SToby Isaac 16669be51f97SToby Isaac Level: advanced 16679be51f97SToby Isaac 16689be51f97SToby Isaac .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel() 16699be51f97SToby Isaac M*/ 16709be51f97SToby Isaac 1671db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm) 1672db4d5e8cSToby Isaac { 1673db4d5e8cSToby Isaac DM_Forest *forest; 1674db4d5e8cSToby Isaac PetscErrorCode ierr; 1675db4d5e8cSToby Isaac 1676db4d5e8cSToby Isaac PetscFunctionBegin; 1677db4d5e8cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1678db4d5e8cSToby Isaac ierr = PetscNewLog(dm,&forest);CHKERRQ(ierr); 1679db4d5e8cSToby Isaac dm->dim = 0; 1680db4d5e8cSToby Isaac dm->data = forest; 1681db4d5e8cSToby Isaac forest->refct = 1; 1682db4d5e8cSToby Isaac forest->data = NULL; 168358762b62SToby Isaac forest->setfromoptionscalled = PETSC_FALSE; 1684db4d5e8cSToby Isaac forest->topology = NULL; 1685cd3c525cSToby Isaac forest->adapt = NULL; 1686db4d5e8cSToby Isaac forest->base = NULL; 16876a87ffbfSToby Isaac forest->adaptPurpose = DM_ADAPT_DETERMINE; 1688db4d5e8cSToby Isaac forest->adjDim = PETSC_DEFAULT; 1689db4d5e8cSToby Isaac forest->overlap = PETSC_DEFAULT; 1690db4d5e8cSToby Isaac forest->minRefinement = PETSC_DEFAULT; 1691db4d5e8cSToby Isaac forest->maxRefinement = PETSC_DEFAULT; 169256ba9f64SToby Isaac forest->initRefinement = PETSC_DEFAULT; 1693c7eeac06SToby Isaac forest->cStart = PETSC_DETERMINE; 1694c7eeac06SToby Isaac forest->cEnd = PETSC_DETERMINE; 1695cd3c525cSToby Isaac forest->cellSF = NULL; 1696ebdf65a2SToby Isaac forest->adaptLabel = NULL; 1697db4d5e8cSToby Isaac forest->gradeFactor = 2; 1698db4d5e8cSToby Isaac forest->cellWeights = NULL; 1699db4d5e8cSToby Isaac forest->cellWeightsCopyMode = PETSC_USE_POINTER; 1700db4d5e8cSToby Isaac forest->weightsFactor = 1.; 1701db4d5e8cSToby Isaac forest->weightCapacity = 1.; 1702a73e2921SToby Isaac ierr = DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);CHKERRQ(ierr); 1703d222f98bSToby Isaac ierr = DMInitialize_Forest(dm);CHKERRQ(ierr); 1704db4d5e8cSToby Isaac PetscFunctionReturn(0); 1705db4d5e8cSToby Isaac } 1706