xref: /petsc/src/dm/impls/forest/forest.c (revision f885a11a8c1ae8666db3d7394d8efe7138068c68)
19be51f97SToby Isaac #include <petsc/private/dmforestimpl.h> /*I "petscdmforest.h" I*/
29be51f97SToby Isaac #include <petsc/private/dmimpl.h>       /*I "petscdm.h" I*/
3ef19d27cSToby Isaac #include <petscsf.h>
4db4d5e8cSToby Isaac 
527d4645fSToby Isaac PetscBool DMForestPackageInitialized = PETSC_FALSE;
627d4645fSToby Isaac 
727d4645fSToby Isaac typedef struct _DMForestTypeLink*DMForestTypeLink;
827d4645fSToby Isaac 
927d4645fSToby Isaac struct _DMForestTypeLink
1027d4645fSToby Isaac {
1127d4645fSToby Isaac   char             *name;
1227d4645fSToby Isaac   DMForestTypeLink next;
1327d4645fSToby Isaac };
1427d4645fSToby Isaac 
1527d4645fSToby Isaac DMForestTypeLink DMForestTypeList;
1627d4645fSToby Isaac 
1727d4645fSToby Isaac #undef __FUNCT__
1827d4645fSToby Isaac #define __FUNCT__ "DMForestPackageFinalize"
1927d4645fSToby Isaac static PetscErrorCode DMForestPackageFinalize(void)
2027d4645fSToby Isaac {
2127d4645fSToby Isaac   DMForestTypeLink oldLink, link = DMForestTypeList;
2227d4645fSToby Isaac   PetscErrorCode   ierr;
2327d4645fSToby Isaac 
2427d4645fSToby Isaac   PetscFunctionBegin;
2527d4645fSToby Isaac   while (link) {
2627d4645fSToby Isaac     oldLink = link;
27f416af30SBarry Smith     ierr    = PetscFree(oldLink->name);CHKERRQ(ierr);
2827d4645fSToby Isaac     link    = oldLink->next;
2927d4645fSToby Isaac     ierr    = PetscFree(oldLink);CHKERRQ(ierr);
3027d4645fSToby Isaac   }
3127d4645fSToby Isaac   PetscFunctionReturn(0);
3227d4645fSToby Isaac }
3327d4645fSToby Isaac 
3427d4645fSToby Isaac #undef __FUNCT__
3527d4645fSToby Isaac #define __FUNCT__ "DMForestPackageInitialize"
3627d4645fSToby Isaac static PetscErrorCode DMForestPackageInitialize(void)
3727d4645fSToby Isaac {
3827d4645fSToby Isaac   PetscErrorCode ierr;
3927d4645fSToby Isaac 
4027d4645fSToby Isaac   PetscFunctionBegin;
4127d4645fSToby Isaac   if (DMForestPackageInitialized) PetscFunctionReturn(0);
4227d4645fSToby Isaac   DMForestPackageInitialized = PETSC_TRUE;
43*f885a11aSToby Isaac 
4427d4645fSToby Isaac   ierr = DMForestRegisterType(DMFOREST);CHKERRQ(ierr);
4527d4645fSToby Isaac   ierr = PetscRegisterFinalize(DMForestPackageFinalize);CHKERRQ(ierr);
4627d4645fSToby Isaac   PetscFunctionReturn(0);
4727d4645fSToby Isaac }
4827d4645fSToby Isaac 
4927d4645fSToby Isaac #undef __FUNCT__
5027d4645fSToby Isaac #define __FUNCT__ "DMForestRegisterType"
519be51f97SToby Isaac /*@C
529be51f97SToby Isaac   DMForestRegisterType - Registers a DMType as a subtype of DMFOREST (so that DMIsForest() will be correct)
539be51f97SToby Isaac 
549be51f97SToby Isaac   Not Collective
559be51f97SToby Isaac 
569be51f97SToby Isaac   Input parameter:
579be51f97SToby Isaac . name - the name of the type
589be51f97SToby Isaac 
599be51f97SToby Isaac   Level: advanced
609be51f97SToby Isaac 
619be51f97SToby Isaac .seealso: DMFOREST, DMIsForest()
629be51f97SToby Isaac @*/
6327d4645fSToby Isaac PetscErrorCode DMForestRegisterType(DMType name)
6427d4645fSToby Isaac {
6527d4645fSToby Isaac   DMForestTypeLink link;
6627d4645fSToby Isaac   PetscErrorCode   ierr;
6727d4645fSToby Isaac 
6827d4645fSToby Isaac   PetscFunctionBegin;
6927d4645fSToby Isaac   ierr             = DMForestPackageInitialize();CHKERRQ(ierr);
7027d4645fSToby Isaac   ierr             = PetscNew(&link);CHKERRQ(ierr);
7127d4645fSToby Isaac   ierr             = PetscStrallocpy(name,&link->name);CHKERRQ(ierr);
7227d4645fSToby Isaac   link->next       = DMForestTypeList;
7327d4645fSToby Isaac   DMForestTypeList = link;
7427d4645fSToby Isaac   PetscFunctionReturn(0);
7527d4645fSToby Isaac }
7627d4645fSToby Isaac 
7727d4645fSToby Isaac #undef __FUNCT__
7827d4645fSToby Isaac #define __FUNCT__ "DMIsForest"
799be51f97SToby Isaac /*@
809be51f97SToby Isaac   DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes
819be51f97SToby Isaac 
829be51f97SToby Isaac   Not Collective
839be51f97SToby Isaac 
849be51f97SToby Isaac   Input parameter:
859be51f97SToby Isaac . dm - the DM object
869be51f97SToby Isaac 
879be51f97SToby Isaac   Output parameter:
889be51f97SToby Isaac . isForest - whether dm is a subtype of DMFOREST
899be51f97SToby Isaac 
909be51f97SToby Isaac   Level: intermediate
919be51f97SToby Isaac 
929be51f97SToby Isaac .seealso: DMFOREST, DMForestRegisterType()
939be51f97SToby Isaac @*/
9427d4645fSToby Isaac PetscErrorCode DMIsForest(DM dm, PetscBool *isForest)
9527d4645fSToby Isaac {
9627d4645fSToby Isaac   DMForestTypeLink link = DMForestTypeList;
9727d4645fSToby Isaac   PetscErrorCode   ierr;
9827d4645fSToby Isaac 
9927d4645fSToby Isaac   PetscFunctionBegin;
10027d4645fSToby Isaac   while (link) {
10127d4645fSToby Isaac     PetscBool sameType;
10227d4645fSToby Isaac     ierr = PetscObjectTypeCompare((PetscObject)dm,link->name,&sameType);CHKERRQ(ierr);
10327d4645fSToby Isaac     if (sameType) {
10427d4645fSToby Isaac       *isForest = PETSC_TRUE;
10527d4645fSToby Isaac       PetscFunctionReturn(0);
10627d4645fSToby Isaac     }
10727d4645fSToby Isaac     link = link->next;
10827d4645fSToby Isaac   }
10927d4645fSToby Isaac   *isForest = PETSC_FALSE;
11027d4645fSToby Isaac   PetscFunctionReturn(0);
11127d4645fSToby Isaac }
11227d4645fSToby Isaac 
113db4d5e8cSToby Isaac #undef __FUNCT__
114a0452a8eSToby Isaac #define __FUNCT__ "DMForestTemplate"
1159be51f97SToby Isaac /*@
1169be51f97SToby Isaac   DMForestTemplate - Create a new DM that will be adapted from a source DM.  The new DM reproduces the configuration
1179be51f97SToby Isaac   of the source, but is not yet setup, so that the user can then define only the ways that the new DM should differ
1189be51f97SToby Isaac   (by, e.g., refinement or repartitioning).  The source DM is also set as the adaptivity source DM of the new DM (see
1199be51f97SToby Isaac   DMForestSetAdaptivityForest()).
1209be51f97SToby Isaac 
1219be51f97SToby Isaac   Collective on dm
1229be51f97SToby Isaac 
1239be51f97SToby Isaac   Input Parameters:
1249be51f97SToby Isaac + dm - the source DM object
1259be51f97SToby Isaac - comm - the communicator for the new DM (this communicator is currently ignored, but is present so that DMForestTemplate() can be used within DMCoarsen())
1269be51f97SToby Isaac 
1279be51f97SToby Isaac   Output Parameter:
1289be51f97SToby Isaac . tdm - the new DM object
1299be51f97SToby Isaac 
1309be51f97SToby Isaac   Level: intermediate
1319be51f97SToby Isaac 
1329be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest()
1339be51f97SToby Isaac @*/
1349be51f97SToby Isaac PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm)
135a0452a8eSToby Isaac {
136a0452a8eSToby Isaac   DM_Forest                  *forest = (DM_Forest*) dm->data;
13720e8089bSToby Isaac   DMType                     type;
138a0452a8eSToby Isaac   DM                         base;
139a0452a8eSToby Isaac   DMForestTopology           topology;
140a0452a8eSToby Isaac   PetscInt                   dim, overlap, ref, factor;
141a0452a8eSToby Isaac   DMForestAdaptivityStrategy strat;
142795844e7SToby Isaac   PetscDS                    ds;
143795844e7SToby Isaac   void                       *ctx;
14449fc9a2fSToby Isaac   PetscErrorCode             (*map)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void*);
1453e58adeeSToby Isaac   void                       *mapCtx;
146a0452a8eSToby Isaac   PetscErrorCode             ierr;
147a0452a8eSToby Isaac 
148a0452a8eSToby Isaac   PetscFunctionBegin;
149a0452a8eSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
15020e8089bSToby Isaac   ierr = DMCreate(PetscObjectComm((PetscObject)dm),tdm);CHKERRQ(ierr);
15120e8089bSToby Isaac   ierr = DMGetType(dm,&type);CHKERRQ(ierr);
15220e8089bSToby Isaac   ierr = DMSetType(*tdm,type);CHKERRQ(ierr);
153a0452a8eSToby Isaac   ierr = DMForestGetBaseDM(dm,&base);CHKERRQ(ierr);
15420e8089bSToby Isaac   ierr = DMForestSetBaseDM(*tdm,base);CHKERRQ(ierr);
155a0452a8eSToby Isaac   ierr = DMForestGetTopology(dm,&topology);CHKERRQ(ierr);
15620e8089bSToby Isaac   ierr = DMForestSetTopology(*tdm,topology);CHKERRQ(ierr);
157a0452a8eSToby Isaac   ierr = DMForestGetAdjacencyDimension(dm,&dim);CHKERRQ(ierr);
15820e8089bSToby Isaac   ierr = DMForestSetAdjacencyDimension(*tdm,dim);CHKERRQ(ierr);
159a0452a8eSToby Isaac   ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr);
16020e8089bSToby Isaac   ierr = DMForestSetPartitionOverlap(*tdm,overlap);CHKERRQ(ierr);
161a0452a8eSToby Isaac   ierr = DMForestGetMinimumRefinement(dm,&ref);CHKERRQ(ierr);
16220e8089bSToby Isaac   ierr = DMForestSetMinimumRefinement(*tdm,ref);CHKERRQ(ierr);
163a0452a8eSToby Isaac   ierr = DMForestGetMaximumRefinement(dm,&ref);CHKERRQ(ierr);
16420e8089bSToby Isaac   ierr = DMForestSetMaximumRefinement(*tdm,ref);CHKERRQ(ierr);
165a0452a8eSToby Isaac   ierr = DMForestGetAdaptivityStrategy(dm,&strat);CHKERRQ(ierr);
16620e8089bSToby Isaac   ierr = DMForestSetAdaptivityStrategy(*tdm,strat);CHKERRQ(ierr);
167a0452a8eSToby Isaac   ierr = DMForestGetGradeFactor(dm,&factor);CHKERRQ(ierr);
16820e8089bSToby Isaac   ierr = DMForestSetGradeFactor(*tdm,factor);CHKERRQ(ierr);
1693e58adeeSToby Isaac   ierr = DMForestGetBaseCoordinateMapping(dm,&map,&mapCtx);CHKERRQ(ierr);
1705e8c540aSToby Isaac   ierr = DMForestSetBaseCoordinateMapping(*tdm,map,mapCtx);CHKERRQ(ierr);
171a0452a8eSToby Isaac   if (forest->ftemplate) {
17220e8089bSToby Isaac     ierr = (forest->ftemplate)(dm, *tdm);CHKERRQ(ierr);
173a0452a8eSToby Isaac   }
17420e8089bSToby Isaac   ierr = DMForestSetAdaptivityForest(*tdm,dm);CHKERRQ(ierr);
175795844e7SToby Isaac   ierr = DMGetDS(dm,&ds);CHKERRQ(ierr);
176795844e7SToby Isaac   ierr = DMSetDS(*tdm,ds);CHKERRQ(ierr);
177795844e7SToby Isaac   ierr = DMGetApplicationContext(dm,&ctx);CHKERRQ(ierr);
178795844e7SToby Isaac   ierr = DMSetApplicationContext(*tdm,&ctx);CHKERRQ(ierr);
179795844e7SToby Isaac   if (dm->maxCell) {
180795844e7SToby Isaac     const PetscReal      *maxCell, *L;
181795844e7SToby Isaac     const DMBoundaryType *bd;
182795844e7SToby Isaac 
183795844e7SToby Isaac     ierr = DMGetPeriodicity(dm,&maxCell,&L,&bd);CHKERRQ(ierr);
184795844e7SToby Isaac     ierr = DMSetPeriodicity(*tdm,maxCell,L,bd);CHKERRQ(ierr);
185795844e7SToby Isaac   }
186a0452a8eSToby Isaac   PetscFunctionReturn(0);
187a0452a8eSToby Isaac }
188a0452a8eSToby Isaac 
18901d9d024SToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm);
19001d9d024SToby Isaac 
191a0452a8eSToby Isaac #undef __FUNCT__
192db4d5e8cSToby Isaac #define __FUNCT__ "DMClone_Forest"
193db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm)
194db4d5e8cSToby Isaac {
195db4d5e8cSToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
196db4d5e8cSToby Isaac   const char     *type;
197db4d5e8cSToby Isaac   PetscErrorCode ierr;
198db4d5e8cSToby Isaac 
199db4d5e8cSToby Isaac   PetscFunctionBegin;
200db4d5e8cSToby Isaac   forest->refct++;
201db4d5e8cSToby Isaac   (*newdm)->data = forest;
202db4d5e8cSToby Isaac   ierr           = PetscObjectGetType((PetscObject) dm, &type);CHKERRQ(ierr);
203db4d5e8cSToby Isaac   ierr           = PetscObjectChangeTypeName((PetscObject) *newdm, type);CHKERRQ(ierr);
20401d9d024SToby Isaac   ierr           = DMInitialize_Forest(*newdm);CHKERRQ(ierr);
205db4d5e8cSToby Isaac   PetscFunctionReturn(0);
206db4d5e8cSToby Isaac }
207db4d5e8cSToby Isaac 
208db4d5e8cSToby Isaac #undef __FUNCT__
209db4d5e8cSToby Isaac #define __FUNCT__ "DMDestroy_Forest"
210d222f98bSToby Isaac static PetscErrorCode DMDestroy_Forest(DM dm)
211db4d5e8cSToby Isaac {
212db4d5e8cSToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
213db4d5e8cSToby Isaac   PetscErrorCode ierr;
214db4d5e8cSToby Isaac 
215db4d5e8cSToby Isaac   PetscFunctionBegin;
216db4d5e8cSToby Isaac   if (--forest->refct > 0) PetscFunctionReturn(0);
217d222f98bSToby Isaac   if (forest->destroy) {ierr = forest->destroy(dm);CHKERRQ(ierr);}
218db4d5e8cSToby Isaac   ierr = PetscSFDestroy(&forest->cellSF);CHKERRQ(ierr);
2190f17b9e3SToby Isaac   ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr);
2200f17b9e3SToby Isaac   ierr = PetscSFDestroy(&forest->coarseToPreFine);CHKERRQ(ierr);
221ebdf65a2SToby Isaac   ierr = PetscFree(forest->adaptLabel);CHKERRQ(ierr);
2229a81d013SToby Isaac   ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr);
22356ba9f64SToby Isaac   ierr = DMDestroy(&forest->base);CHKERRQ(ierr);
224ba936b91SToby Isaac   ierr = DMDestroy(&forest->adapt);CHKERRQ(ierr);
22530f902e7SToby Isaac   ierr = PetscFree(forest->topology);CHKERRQ(ierr);
22630f902e7SToby Isaac   ierr = PetscFree(forest);CHKERRQ(ierr);
227db4d5e8cSToby Isaac   PetscFunctionReturn(0);
228db4d5e8cSToby Isaac }
229db4d5e8cSToby Isaac 
230db4d5e8cSToby Isaac #undef __FUNCT__
231dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetTopology"
2329be51f97SToby Isaac /*@C
2339be51f97SToby Isaac   DMForestSetTopology - Set the topology of a DMForest during the pre-setup phase.  The topology is a string (e.g.
2349be51f97SToby Isaac   "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest durint
2359be51f97SToby Isaac   DMSetUp().
2369be51f97SToby Isaac 
2379be51f97SToby Isaac   Logically collective on dm
2389be51f97SToby Isaac 
2399be51f97SToby Isaac   Input parameters:
2409be51f97SToby Isaac + dm - the forest
2419be51f97SToby Isaac - topology - the topology of the forest
2429be51f97SToby Isaac 
2439be51f97SToby Isaac   Level: intermediate
2449be51f97SToby Isaac 
2459be51f97SToby Isaac .seealso(): DMForestGetTopology(), DMForestSetBaseDM()
2469be51f97SToby Isaac @*/
247dd8e54a2SToby Isaac PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology)
248db4d5e8cSToby Isaac {
249db4d5e8cSToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
250db4d5e8cSToby Isaac   PetscErrorCode ierr;
251db4d5e8cSToby Isaac 
252db4d5e8cSToby Isaac   PetscFunctionBegin;
253db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
254ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the topology after setup");
255dd8e54a2SToby Isaac   ierr = PetscFree(forest->topology);CHKERRQ(ierr);
256dd8e54a2SToby Isaac   ierr = PetscStrallocpy((const char*)topology,(char**) &forest->topology);CHKERRQ(ierr);
257db4d5e8cSToby Isaac   PetscFunctionReturn(0);
258db4d5e8cSToby Isaac }
259db4d5e8cSToby Isaac 
260db4d5e8cSToby Isaac #undef __FUNCT__
261dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetTopology"
2629be51f97SToby Isaac /*@C
2639be51f97SToby Isaac   DMForestGetTopology - Get a string describing the topology of a DMForest.
2649be51f97SToby Isaac 
2659be51f97SToby Isaac   Not collective
2669be51f97SToby Isaac 
2679be51f97SToby Isaac   Input parameter:
2689be51f97SToby Isaac . dm - the forest
2699be51f97SToby Isaac 
2709be51f97SToby Isaac   Output parameter:
2719be51f97SToby Isaac . topology - the topology of the forest (e.g., 'cube', 'shell')
2729be51f97SToby Isaac 
2739be51f97SToby Isaac   Level: intermediate
2749be51f97SToby Isaac 
2759be51f97SToby Isaac .seealso: DMForestSetTopology()
2769be51f97SToby Isaac @*/
277dd8e54a2SToby Isaac PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology)
278dd8e54a2SToby Isaac {
279dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
280dd8e54a2SToby Isaac 
281dd8e54a2SToby Isaac   PetscFunctionBegin;
282dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
283dd8e54a2SToby Isaac   PetscValidPointer(topology,2);
284dd8e54a2SToby Isaac   *topology = forest->topology;
285dd8e54a2SToby Isaac   PetscFunctionReturn(0);
286dd8e54a2SToby Isaac }
287dd8e54a2SToby Isaac 
288dd8e54a2SToby Isaac #undef __FUNCT__
289dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetBaseDM"
2909be51f97SToby Isaac /*@
2919be51f97SToby Isaac   DMForestSetBaseDM - During the pre-setup phase, set the DM that defines the base mesh of a DMForest forest.  The
2929be51f97SToby Isaac   forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its
2939be51f97SToby Isaac   base.  In general, two forest must share a bse to be comparable, to do things like construct interpolators.
2949be51f97SToby Isaac 
2959be51f97SToby Isaac   Logically collective on dm
2969be51f97SToby Isaac 
2979be51f97SToby Isaac   Input Parameters:
2989be51f97SToby Isaac + dm - the forest
2999be51f97SToby Isaac - base - the base DM of the forest
3009be51f97SToby Isaac 
3019be51f97SToby Isaac   Level: intermediate
3029be51f97SToby Isaac 
3039be51f97SToby Isaac .seealso(): DMForestGetBaseDM()
3049be51f97SToby Isaac @*/
305dd8e54a2SToby Isaac PetscErrorCode DMForestSetBaseDM(DM dm, DM base)
306dd8e54a2SToby Isaac {
307dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
308dd8e54a2SToby Isaac   PetscInt       dim, dimEmbed;
309dd8e54a2SToby Isaac   PetscErrorCode ierr;
310dd8e54a2SToby Isaac 
311dd8e54a2SToby Isaac   PetscFunctionBegin;
312dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
313ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the base after setup");
314dd8e54a2SToby Isaac   ierr         = PetscObjectReference((PetscObject)base);CHKERRQ(ierr);
315dd8e54a2SToby Isaac   ierr         = DMDestroy(&forest->base);CHKERRQ(ierr);
316dd8e54a2SToby Isaac   forest->base = base;
317a0452a8eSToby Isaac   if (base) {
318a0452a8eSToby Isaac     PetscValidHeaderSpecific(base, DM_CLASSID, 2);
319dd8e54a2SToby Isaac     ierr = DMGetDimension(base,&dim);CHKERRQ(ierr);
320dd8e54a2SToby Isaac     ierr = DMSetDimension(dm,dim);CHKERRQ(ierr);
321dd8e54a2SToby Isaac     ierr = DMGetCoordinateDim(base,&dimEmbed);CHKERRQ(ierr);
322dd8e54a2SToby Isaac     ierr = DMSetCoordinateDim(dm,dimEmbed);CHKERRQ(ierr);
323a0452a8eSToby Isaac   }
324dd8e54a2SToby Isaac   PetscFunctionReturn(0);
325dd8e54a2SToby Isaac }
326dd8e54a2SToby Isaac 
327dd8e54a2SToby Isaac #undef __FUNCT__
328dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetBaseDM"
3299be51f97SToby Isaac /*@
3309be51f97SToby Isaac   DMForestGetBaseDM - Get the base DM of a DMForest forest.  The forest will be hierarchically refined from the base,
3319be51f97SToby Isaac   and all refinements/coarsenings of the forest will share its base.  In general, two forest must share a bse to be
3329be51f97SToby Isaac   comparable, to do things like construct interpolators.
3339be51f97SToby Isaac 
3349be51f97SToby Isaac   Not collective
3359be51f97SToby Isaac 
3369be51f97SToby Isaac   Input Parameter:
3379be51f97SToby Isaac . dm - the forest
3389be51f97SToby Isaac 
3399be51f97SToby Isaac   Output Parameter:
3409be51f97SToby Isaac . base - the base DM of the forest
3419be51f97SToby Isaac 
3429be51f97SToby Isaac   Level: intermediate
3439be51f97SToby Isaac 
3449be51f97SToby Isaac .seealso(); DMForestSetBaseDM()
3459be51f97SToby Isaac @*/
346dd8e54a2SToby Isaac PetscErrorCode DMForestGetBaseDM(DM dm, DM *base)
347dd8e54a2SToby Isaac {
348dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
349dd8e54a2SToby Isaac 
350dd8e54a2SToby Isaac   PetscFunctionBegin;
351dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
352dd8e54a2SToby Isaac   PetscValidPointer(base, 2);
353dd8e54a2SToby Isaac   *base = forest->base;
354dd8e54a2SToby Isaac   PetscFunctionReturn(0);
355dd8e54a2SToby Isaac }
356dd8e54a2SToby Isaac 
357dd8e54a2SToby Isaac #undef __FUNCT__
358cf38a08cSToby Isaac #define __FUNCT__ "DMForestSetBaseCoordinateMapping"
35999478f86SToby Isaac PetscErrorCode DMForestSetBaseCoordinateMapping(DM dm, PetscErrorCode (*func)(DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx)
360cf38a08cSToby Isaac {
361cf38a08cSToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
362cf38a08cSToby Isaac 
363cf38a08cSToby Isaac   PetscFunctionBegin;
364cf38a08cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
365cf38a08cSToby Isaac   forest->mapcoordinates    = func;
366cf38a08cSToby Isaac   forest->mapcoordinatesctx = ctx;
367cf38a08cSToby Isaac   PetscFunctionReturn(0);
368cf38a08cSToby Isaac }
369cf38a08cSToby Isaac 
370cf38a08cSToby Isaac #undef __FUNCT__
371cf38a08cSToby Isaac #define __FUNCT__ "DMForestGetBaseCoordinateMapping"
37299478f86SToby Isaac PetscErrorCode DMForestGetBaseCoordinateMapping(DM dm, PetscErrorCode (**func) (DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx)
373cf38a08cSToby Isaac {
374cf38a08cSToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
375cf38a08cSToby Isaac 
376cf38a08cSToby Isaac   PetscFunctionBegin;
377cf38a08cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
378cf38a08cSToby Isaac   if (func) *func = forest->mapcoordinates;
379cf38a08cSToby Isaac   if (ctx) *((void**) ctx) = forest->mapcoordinatesctx;
380cf38a08cSToby Isaac   PetscFunctionReturn(0);
381cf38a08cSToby Isaac }
382cf38a08cSToby Isaac 
383cf38a08cSToby Isaac #undef __FUNCT__
384ba936b91SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityForest"
3859be51f97SToby Isaac /*@
3869be51f97SToby Isaac   DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be
3879be51f97SToby Isaac   adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) im DMSetUp().  Usually not needed
3889be51f97SToby Isaac   by users directly: DMForestTemplate() constructs a new forest to be adapted from an old forest and calls this
3899be51f97SToby Isaac   routine.
3909be51f97SToby Isaac 
3919be51f97SToby Isaac   Logically collective on dm
3929be51f97SToby Isaac 
3939be51f97SToby Isaac   Input Parameter:
3949be51f97SToby Isaac + dm - the new forest, which will be constructed from adapt
3959be51f97SToby Isaac - adapt - the old forest
3969be51f97SToby Isaac 
3979be51f97SToby Isaac   Level: intermediate
3989be51f97SToby Isaac 
3999be51f97SToby Isaac .seealso: DMForestGetAdaptivityForest(), DMForestSetAdaptivityPurpose()
4009be51f97SToby Isaac @*/
401ba936b91SToby Isaac PetscErrorCode DMForestSetAdaptivityForest(DM dm,DM adapt)
402dd8e54a2SToby Isaac {
403ba936b91SToby Isaac   DM_Forest      *forest;
404dd8e54a2SToby Isaac   PetscErrorCode ierr;
405dd8e54a2SToby Isaac 
406dd8e54a2SToby Isaac   PetscFunctionBegin;
407dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
408ba936b91SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
409ba936b91SToby Isaac   forest = (DM_Forest*) dm->data;
410ba936b91SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
41126d9498aSToby Isaac   switch (forest->adaptPurpose) {
41226d9498aSToby Isaac   case DM_FOREST_KEEP:
413ba936b91SToby Isaac     ierr          = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr);
414ba936b91SToby Isaac     ierr          = DMDestroy(&(forest->adapt));CHKERRQ(ierr);
415ba936b91SToby Isaac     forest->adapt = adapt;
41626d9498aSToby Isaac     break;
41726d9498aSToby Isaac   case DM_FOREST_REFINE:
41826d9498aSToby Isaac     ierr = DMSetCoarseDM(dm,adapt);CHKERRQ(ierr);
41926d9498aSToby Isaac     break;
42026d9498aSToby Isaac   case DM_FOREST_COARSEN:
42126d9498aSToby Isaac     ierr = DMSetFineDM(dm,adapt);CHKERRQ(ierr);
42226d9498aSToby Isaac     break;
42326d9498aSToby Isaac   default:
42426d9498aSToby Isaac     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
42526d9498aSToby Isaac   }
426dd8e54a2SToby Isaac   PetscFunctionReturn(0);
427dd8e54a2SToby Isaac }
428dd8e54a2SToby Isaac 
429dd8e54a2SToby Isaac #undef __FUNCT__
430ba936b91SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityForest"
4319be51f97SToby Isaac /*@
4329be51f97SToby Isaac   DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted.
4339be51f97SToby Isaac 
4349be51f97SToby Isaac   Not collective
4359be51f97SToby Isaac 
4369be51f97SToby Isaac   Input Parameter:
4379be51f97SToby Isaac . dm - the forest
4389be51f97SToby Isaac 
4399be51f97SToby Isaac   Output Parameter:
4409be51f97SToby Isaac . adapt - the forest from which dm is/was adapted
4419be51f97SToby Isaac 
4429be51f97SToby Isaac   Level: intermediate
4439be51f97SToby Isaac 
4449be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose()
4459be51f97SToby Isaac @*/
446ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt)
447dd8e54a2SToby Isaac {
448ba936b91SToby Isaac   DM_Forest      *forest;
44926d9498aSToby Isaac   PetscErrorCode ierr;
450dd8e54a2SToby Isaac 
451dd8e54a2SToby Isaac   PetscFunctionBegin;
452dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
453ba936b91SToby Isaac   forest = (DM_Forest*) dm->data;
45426d9498aSToby Isaac   switch (forest->adaptPurpose) {
45526d9498aSToby Isaac   case DM_FOREST_KEEP:
456ba936b91SToby Isaac     *adapt = forest->adapt;
45726d9498aSToby Isaac     break;
45826d9498aSToby Isaac   case DM_FOREST_REFINE:
45926d9498aSToby Isaac     ierr = DMGetCoarseDM(dm,adapt);CHKERRQ(ierr);
46026d9498aSToby Isaac     break;
46126d9498aSToby Isaac   case DM_FOREST_COARSEN:
46226d9498aSToby Isaac     ierr = DMGetFineDM(dm,adapt);CHKERRQ(ierr);
46326d9498aSToby Isaac     break;
46426d9498aSToby Isaac   default:
46526d9498aSToby Isaac     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
46626d9498aSToby Isaac   }
46726d9498aSToby Isaac   PetscFunctionReturn(0);
46826d9498aSToby Isaac }
46926d9498aSToby Isaac 
47026d9498aSToby Isaac #undef __FUNCT__
47126d9498aSToby Isaac #define __FUNCT__ "DMForestSetAdaptivityPurpose"
4729be51f97SToby Isaac /*@
4739be51f97SToby Isaac   DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its
4749be51f97SToby Isaac   source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_FOREST_REFINE), coarsening
4759be51f97SToby Isaac   (DM_FOREST_COARSEN), or undefined (DM_FOREST_NONE).  This only matters for the purposes of reference counting:
4769be51f97SToby Isaac   during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse
4779be51f97SToby Isaac   relationship (see DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does
4789be51f97SToby Isaac   not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can
4799be51f97SToby Isaac   cause a memory leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
4809be51f97SToby Isaac 
4819be51f97SToby Isaac   Logically collective on dm
4829be51f97SToby Isaac 
4839be51f97SToby Isaac   Input Parameters:
4849be51f97SToby Isaac + dm - the forest
4859be51f97SToby Isaac - purpose - the adaptivity purpose (DM_FOREST_NONE/DM_FOREST_REFINE/DM_FOREST_COARSEN)
4869be51f97SToby Isaac 
4879be51f97SToby Isaac   Level: advanced
4889be51f97SToby Isaac 
4899be51f97SToby Isaac .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest()
4909be51f97SToby Isaac @*/
49126d9498aSToby Isaac PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMForestAdaptivityPurpose purpose)
49226d9498aSToby Isaac {
49326d9498aSToby Isaac   DM_Forest      *forest;
49426d9498aSToby Isaac   PetscErrorCode ierr;
49526d9498aSToby Isaac 
49626d9498aSToby Isaac   PetscFunctionBegin;
49726d9498aSToby Isaac   forest = (DM_Forest*) dm->data;
4989be51f97SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
49926d9498aSToby Isaac   if (purpose != forest->adaptPurpose) {
50026d9498aSToby Isaac     DM adapt;
50126d9498aSToby Isaac 
50226d9498aSToby Isaac     ierr = DMForestGetAdaptivityForest(dm,&adapt);CHKERRQ(ierr);
50326d9498aSToby Isaac     ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr);
50426d9498aSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
505*f885a11aSToby Isaac 
50626d9498aSToby Isaac     forest->adaptPurpose = purpose;
507*f885a11aSToby Isaac 
50826d9498aSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,adapt);CHKERRQ(ierr);
50926d9498aSToby Isaac     ierr = DMDestroy(&adapt);CHKERRQ(ierr);
51026d9498aSToby Isaac   }
511dd8e54a2SToby Isaac   PetscFunctionReturn(0);
512dd8e54a2SToby Isaac }
513dd8e54a2SToby Isaac 
514dd8e54a2SToby Isaac #undef __FUNCT__
51556c0450aSToby Isaac #define __FUNCT__ "DMForestGetAdaptivityPurpose"
51656c0450aSToby Isaac /*@
51756c0450aSToby Isaac   DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with
51856c0450aSToby Isaac   DMForestSetAdaptivityForest()) for the purpose of refinement (DM_FOREST_REFINE), coarsening (DM_FOREST_COARSEN), or
51956c0450aSToby Isaac   undefined (DM_FOREST_NONE).  This only matters for the purposes of reference counting: during DMDestroy(), cyclic
52056c0450aSToby Isaac   references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see
52156c0450aSToby Isaac   DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does not maintain a
52256c0450aSToby Isaac   reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory
52356c0450aSToby Isaac   leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
52456c0450aSToby Isaac 
52556c0450aSToby Isaac   Not collective
52656c0450aSToby Isaac 
52756c0450aSToby Isaac   Input Parameter:
52856c0450aSToby Isaac . dm - the forest
52956c0450aSToby Isaac 
53056c0450aSToby Isaac   Output Parameter:
53156c0450aSToby Isaac . purpose - the adaptivity purpose (DM_FOREST_NONE/DM_FOREST_REFINE/DM_FOREST_COARSEN)
53256c0450aSToby Isaac 
53356c0450aSToby Isaac   Level: advanced
53456c0450aSToby Isaac 
53556c0450aSToby Isaac .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest()
53656c0450aSToby Isaac @*/
53756c0450aSToby Isaac PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMForestAdaptivityPurpose *purpose)
53856c0450aSToby Isaac {
53956c0450aSToby Isaac   DM_Forest *forest;
54056c0450aSToby Isaac 
54156c0450aSToby Isaac   PetscFunctionBegin;
54256c0450aSToby Isaac   forest   = (DM_Forest*) dm->data;
54356c0450aSToby Isaac   *purpose = forest->adaptPurpose;
54456c0450aSToby Isaac   PetscFunctionReturn(0);
54556c0450aSToby Isaac }
54656c0450aSToby Isaac 
54756c0450aSToby Isaac #undef __FUNCT__
548dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetAdjacencyDimension"
5499be51f97SToby Isaac /*@
5509be51f97SToby Isaac   DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine
5519be51f97SToby Isaac   cell adjacency (for the purposes of partitioning and overlap).
5529be51f97SToby Isaac 
5539be51f97SToby Isaac   Logically collective on dm
5549be51f97SToby Isaac 
5559be51f97SToby Isaac   Input Parameters:
5569be51f97SToby Isaac + dm - the forest
5579be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency)
5589be51f97SToby Isaac 
5599be51f97SToby Isaac   Level: intermediate
5609be51f97SToby Isaac 
5619be51f97SToby Isaac .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap()
5629be51f97SToby Isaac @*/
563dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim)
564dd8e54a2SToby Isaac {
565dd8e54a2SToby Isaac   PetscInt       dim;
566dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
567dd8e54a2SToby Isaac   PetscErrorCode ierr;
568dd8e54a2SToby Isaac 
569dd8e54a2SToby Isaac   PetscFunctionBegin;
570dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
571ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup");
572dd8e54a2SToby Isaac   if (adjDim < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %d", adjDim);
573dd8e54a2SToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
574dd8e54a2SToby Isaac   if (adjDim > dim) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %d: %d", dim, adjDim);
575dd8e54a2SToby Isaac   forest->adjDim = adjDim;
576dd8e54a2SToby Isaac   PetscFunctionReturn(0);
577dd8e54a2SToby Isaac }
578dd8e54a2SToby Isaac 
579dd8e54a2SToby Isaac #undef __FUNCT__
580dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetAdjacencyCodimension"
5819be51f97SToby Isaac /*@
5829be51f97SToby Isaac   DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that,
5839be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
5849be51f97SToby Isaac 
5859be51f97SToby Isaac   Logically collective on dm
5869be51f97SToby Isaac 
5879be51f97SToby Isaac   Input Parameters:
5889be51f97SToby Isaac + dm - the forest
5899be51f97SToby Isaac - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
5909be51f97SToby Isaac 
5919be51f97SToby Isaac   Level: intermediate
5929be51f97SToby Isaac 
5939be51f97SToby Isaac .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension()
5949be51f97SToby Isaac @*/
595dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim)
596dd8e54a2SToby Isaac {
597dd8e54a2SToby Isaac   PetscInt       dim;
598dd8e54a2SToby Isaac   PetscErrorCode ierr;
599dd8e54a2SToby Isaac 
600dd8e54a2SToby Isaac   PetscFunctionBegin;
601dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
602dd8e54a2SToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
603dd8e54a2SToby Isaac   ierr = DMForestSetAdjacencyDimension(dm,dim-adjCodim);CHKERRQ(ierr);
604dd8e54a2SToby Isaac   PetscFunctionReturn(0);
605dd8e54a2SToby Isaac }
606dd8e54a2SToby Isaac 
607dd8e54a2SToby Isaac #undef __FUNCT__
608dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetAdjacencyDimension"
6099be51f97SToby Isaac /*@
6109be51f97SToby Isaac   DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the
6119be51f97SToby Isaac   purposes of partitioning and overlap).
6129be51f97SToby Isaac 
6139be51f97SToby Isaac   Not collective
6149be51f97SToby Isaac 
6159be51f97SToby Isaac   Input Parameter:
6169be51f97SToby Isaac . dm - the forest
6179be51f97SToby Isaac 
6189be51f97SToby Isaac   Output Parameter:
6199be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency)
6209be51f97SToby Isaac 
6219be51f97SToby Isaac   Level: intermediate
6229be51f97SToby Isaac 
6239be51f97SToby Isaac .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap()
6249be51f97SToby Isaac @*/
625dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim)
626dd8e54a2SToby Isaac {
627dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
628dd8e54a2SToby Isaac 
629dd8e54a2SToby Isaac   PetscFunctionBegin;
630dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
631dd8e54a2SToby Isaac   PetscValidIntPointer(adjDim,2);
632dd8e54a2SToby Isaac   *adjDim = forest->adjDim;
633dd8e54a2SToby Isaac   PetscFunctionReturn(0);
634dd8e54a2SToby Isaac }
635dd8e54a2SToby Isaac 
636dd8e54a2SToby Isaac #undef __FUNCT__
637dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetAdjacencyCodimension"
6389be51f97SToby Isaac /*@
6399be51f97SToby Isaac   DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that,
6409be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
6419be51f97SToby Isaac 
6429be51f97SToby Isaac   Not collective
6439be51f97SToby Isaac 
6449be51f97SToby Isaac   Input Parameter:
6459be51f97SToby Isaac . dm - the forest
6469be51f97SToby Isaac 
6479be51f97SToby Isaac   Output Parameter:
6489be51f97SToby Isaac . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
6499be51f97SToby Isaac 
6509be51f97SToby Isaac   Level: intermediate
6519be51f97SToby Isaac 
6529be51f97SToby Isaac .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension()
6539be51f97SToby Isaac @*/
654dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim)
655dd8e54a2SToby Isaac {
656dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
657dd8e54a2SToby Isaac   PetscInt       dim;
658dd8e54a2SToby Isaac   PetscErrorCode ierr;
659dd8e54a2SToby Isaac 
660dd8e54a2SToby Isaac   PetscFunctionBegin;
661dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
662dd8e54a2SToby Isaac   PetscValidIntPointer(adjCodim,2);
663dd8e54a2SToby Isaac   ierr      = DMGetDimension(dm,&dim);CHKERRQ(ierr);
664dd8e54a2SToby Isaac   *adjCodim = dim - forest->adjDim;
665dd8e54a2SToby Isaac   PetscFunctionReturn(0);
666dd8e54a2SToby Isaac }
667dd8e54a2SToby Isaac 
668dd8e54a2SToby Isaac #undef __FUNCT__
669ef51cf95SToby Isaac #define __FUNCT__ "DMForestSetPartitionOverlap"
6709be51f97SToby Isaac /*@
6719be51f97SToby Isaac   DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel
6729be51f97SToby Isaac   partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding
6739be51f97SToby Isaac   adjacent cells
6749be51f97SToby Isaac 
6759be51f97SToby Isaac   Logically collective on dm
6769be51f97SToby Isaac 
6779be51f97SToby Isaac   Input Parameters:
6789be51f97SToby Isaac + dm - the forest
6799be51f97SToby Isaac - overlap - default 0
6809be51f97SToby Isaac 
6819be51f97SToby Isaac   Level: intermediate
6829be51f97SToby Isaac 
6839be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
6849be51f97SToby Isaac @*/
685dd8e54a2SToby Isaac PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap)
686dd8e54a2SToby Isaac {
687dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
688dd8e54a2SToby Isaac 
689dd8e54a2SToby Isaac   PetscFunctionBegin;
690dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
691ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup");
692dd8e54a2SToby Isaac   if (overlap < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %d", overlap);
693dd8e54a2SToby Isaac   forest->overlap = overlap;
694dd8e54a2SToby Isaac   PetscFunctionReturn(0);
695dd8e54a2SToby Isaac }
696dd8e54a2SToby Isaac 
697dd8e54a2SToby Isaac #undef __FUNCT__
698dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetPartitionOverlap"
6999be51f97SToby Isaac /*@
7009be51f97SToby Isaac   DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values
7019be51f97SToby Isaac   > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells
7029be51f97SToby Isaac 
7039be51f97SToby Isaac   Not collective
7049be51f97SToby Isaac 
7059be51f97SToby Isaac   Input Parameter:
7069be51f97SToby Isaac . dm - the forest
7079be51f97SToby Isaac 
7089be51f97SToby Isaac   Output Parameter:
7099be51f97SToby Isaac . overlap - default 0
7109be51f97SToby Isaac 
7119be51f97SToby Isaac   Level: intermediate
7129be51f97SToby Isaac 
7139be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
7149be51f97SToby Isaac @*/
715dd8e54a2SToby Isaac PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap)
716dd8e54a2SToby Isaac {
717dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
718dd8e54a2SToby Isaac 
719dd8e54a2SToby Isaac   PetscFunctionBegin;
720dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
721dd8e54a2SToby Isaac   PetscValidIntPointer(overlap,2);
722dd8e54a2SToby Isaac   *overlap = forest->overlap;
723dd8e54a2SToby Isaac   PetscFunctionReturn(0);
724dd8e54a2SToby Isaac }
725dd8e54a2SToby Isaac 
726dd8e54a2SToby Isaac #undef __FUNCT__
727dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetMinimumRefinement"
7289be51f97SToby Isaac /*@
7299be51f97SToby Isaac   DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base
7309be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest
7319be51f97SToby Isaac   (see DMForestGetAdaptivityForest()) this limits the amount of coarsening.
7329be51f97SToby Isaac 
7339be51f97SToby Isaac   Logically collective on dm
7349be51f97SToby Isaac 
7359be51f97SToby Isaac   Input Parameters:
7369be51f97SToby Isaac + dm - the forest
7379be51f97SToby Isaac - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7389be51f97SToby Isaac 
7399be51f97SToby Isaac   Level: intermediate
7409be51f97SToby Isaac 
7419be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
7429be51f97SToby Isaac @*/
743dd8e54a2SToby Isaac PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement)
744dd8e54a2SToby Isaac {
745dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
746dd8e54a2SToby Isaac 
747dd8e54a2SToby Isaac   PetscFunctionBegin;
748dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
749ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup");
750dd8e54a2SToby Isaac   forest->minRefinement = minRefinement;
751dd8e54a2SToby Isaac   PetscFunctionReturn(0);
752dd8e54a2SToby Isaac }
753dd8e54a2SToby Isaac 
754dd8e54a2SToby Isaac #undef __FUNCT__
755dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetMinimumRefinement"
7569be51f97SToby Isaac /*@
7579be51f97SToby Isaac   DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see
7589be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest (see
7599be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of coarsening.
7609be51f97SToby Isaac 
7619be51f97SToby Isaac   Not collective
7629be51f97SToby Isaac 
7639be51f97SToby Isaac   Input Parameter:
7649be51f97SToby Isaac . dm - the forest
7659be51f97SToby Isaac 
7669be51f97SToby Isaac   Output Parameter:
7679be51f97SToby Isaac . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7689be51f97SToby Isaac 
7699be51f97SToby Isaac   Level: intermediate
7709be51f97SToby Isaac 
7719be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestGetMaximumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
7729be51f97SToby Isaac @*/
773dd8e54a2SToby Isaac PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement)
774dd8e54a2SToby Isaac {
775dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
776dd8e54a2SToby Isaac 
777dd8e54a2SToby Isaac   PetscFunctionBegin;
778dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
779dd8e54a2SToby Isaac   PetscValidIntPointer(minRefinement,2);
780dd8e54a2SToby Isaac   *minRefinement = forest->minRefinement;
781dd8e54a2SToby Isaac   PetscFunctionReturn(0);
782dd8e54a2SToby Isaac }
783dd8e54a2SToby Isaac 
784dd8e54a2SToby Isaac #undef __FUNCT__
78556ba9f64SToby Isaac #define __FUNCT__ "DMForestSetInitialRefinement"
7869be51f97SToby Isaac /*@
7879be51f97SToby Isaac   DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base
7889be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.
7899be51f97SToby Isaac 
7909be51f97SToby Isaac   Logically collective on dm
7919be51f97SToby Isaac 
7929be51f97SToby Isaac   Input Parameters:
7939be51f97SToby Isaac + dm - the forest
7949be51f97SToby Isaac - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7959be51f97SToby Isaac 
7969be51f97SToby Isaac   Level: intermediate
7979be51f97SToby Isaac 
7989be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
7999be51f97SToby Isaac @*/
80056ba9f64SToby Isaac PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement)
80156ba9f64SToby Isaac {
80256ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
80356ba9f64SToby Isaac 
80456ba9f64SToby Isaac   PetscFunctionBegin;
80556ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
80656ba9f64SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup");
80756ba9f64SToby Isaac   forest->initRefinement = initRefinement;
80856ba9f64SToby Isaac   PetscFunctionReturn(0);
80956ba9f64SToby Isaac }
81056ba9f64SToby Isaac 
81156ba9f64SToby Isaac #undef __FUNCT__
81256ba9f64SToby Isaac #define __FUNCT__ "DMForestGetInitialRefinement"
8139be51f97SToby Isaac /*@
8149be51f97SToby Isaac   DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see
8159be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.
8169be51f97SToby Isaac 
8179be51f97SToby Isaac   Not collective
8189be51f97SToby Isaac 
8199be51f97SToby Isaac   Input Parameter:
8209be51f97SToby Isaac . dm - the forest
8219be51f97SToby Isaac 
8229be51f97SToby Isaac   Output Paramater:
8239be51f97SToby Isaac . initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8249be51f97SToby Isaac 
8259be51f97SToby Isaac   Level: intermediate
8269be51f97SToby Isaac 
8279be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
8289be51f97SToby Isaac @*/
82956ba9f64SToby Isaac PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement)
83056ba9f64SToby Isaac {
83156ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
83256ba9f64SToby Isaac 
83356ba9f64SToby Isaac   PetscFunctionBegin;
83456ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
83556ba9f64SToby Isaac   PetscValidIntPointer(initRefinement,2);
83656ba9f64SToby Isaac   *initRefinement = forest->initRefinement;
83756ba9f64SToby Isaac   PetscFunctionReturn(0);
83856ba9f64SToby Isaac }
83956ba9f64SToby Isaac 
84056ba9f64SToby Isaac #undef __FUNCT__
841c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetMaximumRefinement"
8429be51f97SToby Isaac /*@
8439be51f97SToby Isaac   DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base
8449be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest
8459be51f97SToby Isaac   (see DMForestGetAdaptivityForest()), this limits the amount of refinement.
8469be51f97SToby Isaac 
8479be51f97SToby Isaac   Logically collective on dm
8489be51f97SToby Isaac 
8499be51f97SToby Isaac   Input Parameters:
8509be51f97SToby Isaac + dm - the forest
8519be51f97SToby Isaac - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8529be51f97SToby Isaac 
8539be51f97SToby Isaac   Level: intermediate
8549be51f97SToby Isaac 
8559be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityDM()
8569be51f97SToby Isaac @*/
857c7eeac06SToby Isaac PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement)
858dd8e54a2SToby Isaac {
859dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
860dd8e54a2SToby Isaac 
861dd8e54a2SToby Isaac   PetscFunctionBegin;
862dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
863ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup");
864c7eeac06SToby Isaac   forest->maxRefinement = maxRefinement;
865dd8e54a2SToby Isaac   PetscFunctionReturn(0);
866dd8e54a2SToby Isaac }
867dd8e54a2SToby Isaac 
868dd8e54a2SToby Isaac #undef __FUNCT__
869c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetMaximumRefinement"
8709be51f97SToby Isaac /*@
8719be51f97SToby Isaac   DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see
8729be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest (see
8739be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of refinement.
8749be51f97SToby Isaac 
8759be51f97SToby Isaac   Not collective
8769be51f97SToby Isaac 
8779be51f97SToby Isaac   Input Parameter:
8789be51f97SToby Isaac . dm - the forest
8799be51f97SToby Isaac 
8809be51f97SToby Isaac   Output Parameter:
8819be51f97SToby Isaac . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8829be51f97SToby Isaac 
8839be51f97SToby Isaac   Level: intermediate
8849be51f97SToby Isaac 
8859be51f97SToby Isaac .seealso: DMForestSetMaximumRefinement(), DMForestGetMinimumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
8869be51f97SToby Isaac @*/
887c7eeac06SToby Isaac PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement)
888dd8e54a2SToby Isaac {
889dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
890dd8e54a2SToby Isaac 
891dd8e54a2SToby Isaac   PetscFunctionBegin;
892dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
893c7eeac06SToby Isaac   PetscValidIntPointer(maxRefinement,2);
894c7eeac06SToby Isaac   *maxRefinement = forest->maxRefinement;
895dd8e54a2SToby Isaac   PetscFunctionReturn(0);
896dd8e54a2SToby Isaac }
897c7eeac06SToby Isaac 
898c7eeac06SToby Isaac #undef __FUNCT__
899c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityStrategy"
9009be51f97SToby Isaac /*@C
9019be51f97SToby Isaac   DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes.
9029be51f97SToby Isaac   Subtypes of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree
9039be51f97SToby Isaac   for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to
9049be51f97SToby Isaac   specify refinement/coarsening.
9059be51f97SToby Isaac 
9069be51f97SToby Isaac   Logically collective on dm
9079be51f97SToby Isaac 
9089be51f97SToby Isaac   Input Parameters:
9099be51f97SToby Isaac + dm - the forest
9109be51f97SToby Isaac - adaptStrategy - default DMFORESTADAPTALL
9119be51f97SToby Isaac 
9129be51f97SToby Isaac   Level: advanced
9139be51f97SToby Isaac 
9149be51f97SToby Isaac .seealso: DMForestGetAdaptivityStrategy()
9159be51f97SToby Isaac @*/
916c7eeac06SToby Isaac PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy)
917c7eeac06SToby Isaac {
918c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
919c7eeac06SToby Isaac   PetscErrorCode ierr;
920c7eeac06SToby Isaac 
921c7eeac06SToby Isaac   PetscFunctionBegin;
922c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
923c7eeac06SToby Isaac   ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr);
924a73e2921SToby Isaac   ierr = PetscStrallocpy((const char*) adaptStrategy,(char**)&forest->adaptStrategy);CHKERRQ(ierr);
925c7eeac06SToby Isaac   PetscFunctionReturn(0);
926c7eeac06SToby Isaac }
927c7eeac06SToby Isaac 
928c7eeac06SToby Isaac #undef __FUNCT__
929c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityStrategy"
9309be51f97SToby Isaac /*@C
9319be51f97SToby Isaac   DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes.  Subtypes
9329be51f97SToby Isaac   of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all
9339be51f97SToby Isaac   processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only
9349be51f97SToby Isaac   one process needs to specify refinement/coarsening.
9359be51f97SToby Isaac 
9369be51f97SToby Isaac   Not collective
9379be51f97SToby Isaac 
9389be51f97SToby Isaac   Input Parameter:
9399be51f97SToby Isaac . dm - the forest
9409be51f97SToby Isaac 
9419be51f97SToby Isaac   Output Parameter:
9429be51f97SToby Isaac . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL)
9439be51f97SToby Isaac 
9449be51f97SToby Isaac   Level: advanced
9459be51f97SToby Isaac 
9469be51f97SToby Isaac .seealso: DMForestSetAdaptivityStrategy()
9479be51f97SToby Isaac @*/
948c7eeac06SToby Isaac PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy)
949c7eeac06SToby Isaac {
950c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
951c7eeac06SToby Isaac 
952c7eeac06SToby Isaac   PetscFunctionBegin;
953c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
954c7eeac06SToby Isaac   PetscValidPointer(adaptStrategy,2);
955c7eeac06SToby Isaac   *adaptStrategy = forest->adaptStrategy;
956c7eeac06SToby Isaac   PetscFunctionReturn(0);
957c7eeac06SToby Isaac }
958c7eeac06SToby Isaac 
959c7eeac06SToby Isaac #undef __FUNCT__
960bf9b5d84SToby Isaac #define __FUNCT__ "DMForestSetComputeAdaptivitySF"
961bf9b5d84SToby Isaac /*@
962bf9b5d84SToby Isaac   DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed
963bf9b5d84SToby 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().
964bf9b5d84SToby Isaac 
965bf9b5d84SToby Isaac   Logically collective on dm
966bf9b5d84SToby Isaac 
967bf9b5d84SToby Isaac   Input Parameters:
968bf9b5d84SToby Isaac + dm - the post-adaptation forest
969bf9b5d84SToby Isaac - computeSF - default PETSC_TRUE
970bf9b5d84SToby Isaac 
971bf9b5d84SToby Isaac   Level: advanced
972bf9b5d84SToby Isaac 
973bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
974bf9b5d84SToby Isaac @*/
975bf9b5d84SToby Isaac PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF)
976bf9b5d84SToby Isaac {
977bf9b5d84SToby Isaac   DM_Forest *forest;
978bf9b5d84SToby Isaac 
979bf9b5d84SToby Isaac   PetscFunctionBegin;
980bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
981bf9b5d84SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called");
982bf9b5d84SToby Isaac   forest                 = (DM_Forest*) dm->data;
983bf9b5d84SToby Isaac   forest->computeAdaptSF = computeSF;
984bf9b5d84SToby Isaac   PetscFunctionReturn(0);
985bf9b5d84SToby Isaac }
986bf9b5d84SToby Isaac 
987bf9b5d84SToby Isaac #undef __FUNCT__
988bf9b5d84SToby Isaac #define __FUNCT__ "DMForestGetComputeAdaptivitySF"
989bf9b5d84SToby Isaac /*@
990bf9b5d84SToby Isaac   DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the
991bf9b5d84SToby Isaac   pre-adaptation forest to the post-adaptiation forest.  After DMSetUp() is called, these transfer PetscSFs can be
992bf9b5d84SToby Isaac   accessed with DMForestGetAdaptivitySF().
993bf9b5d84SToby Isaac 
994bf9b5d84SToby Isaac   Not collective
995bf9b5d84SToby Isaac 
996bf9b5d84SToby Isaac   Input Parameter:
997bf9b5d84SToby Isaac . dm - the post-adaptation forest
998bf9b5d84SToby Isaac 
999bf9b5d84SToby Isaac   Output Parameter:
1000bf9b5d84SToby Isaac . computeSF - default PETSC_TRUE
1001bf9b5d84SToby Isaac 
1002bf9b5d84SToby Isaac   Level: advanced
1003bf9b5d84SToby Isaac 
1004bf9b5d84SToby Isaac .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
1005bf9b5d84SToby Isaac @*/
1006bf9b5d84SToby Isaac PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF)
1007bf9b5d84SToby Isaac {
1008bf9b5d84SToby Isaac   DM_Forest *forest;
1009bf9b5d84SToby Isaac 
1010bf9b5d84SToby Isaac   PetscFunctionBegin;
1011bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1012bf9b5d84SToby Isaac   forest     = (DM_Forest*) dm->data;
1013bf9b5d84SToby Isaac   *computeSF = forest->computeAdaptSF;
1014bf9b5d84SToby Isaac   PetscFunctionReturn(0);
1015bf9b5d84SToby Isaac }
1016bf9b5d84SToby Isaac 
1017bf9b5d84SToby Isaac #undef __FUNCT__
1018bf9b5d84SToby Isaac #define __FUNCT__ "DMForestGetAdaptivitySF"
1019bf9b5d84SToby Isaac /*@
1020bf9b5d84SToby Isaac   DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest.
1021bf9b5d84SToby Isaac   Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be
1022bf9b5d84SToby Isaac   some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa.  Therefore there are two
1023bf9b5d84SToby Isaac   PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates
1024bf9b5d84SToby Isaac   pre-adaptation fine cells to post-adaptation coarse cells.
1025bf9b5d84SToby Isaac 
1026bf9b5d84SToby Isaac   Not collective
1027bf9b5d84SToby Isaac 
1028bf9b5d84SToby Isaac   Input Parameter:
1029bf9b5d84SToby Isaac   dm - the post-adaptation forest
1030bf9b5d84SToby Isaac 
1031bf9b5d84SToby Isaac   Output Parameter:
10320f17b9e3SToby Isaac   preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post-
10330f17b9e3SToby Isaac   coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre-
1034bf9b5d84SToby Isaac 
1035bf9b5d84SToby Isaac   Level: advanced
1036bf9b5d84SToby Isaac 
1037bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF()
1038bf9b5d84SToby Isaac @*/
10390f17b9e3SToby Isaac PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine)
1040bf9b5d84SToby Isaac {
1041bf9b5d84SToby Isaac   DM_Forest      *forest;
1042bf9b5d84SToby Isaac   PetscErrorCode ierr;
1043bf9b5d84SToby Isaac 
1044bf9b5d84SToby Isaac   PetscFunctionBegin;
1045bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1046bf9b5d84SToby Isaac   ierr   = DMSetUp(dm);CHKERRQ(ierr);
1047bf9b5d84SToby Isaac   forest = (DM_Forest*) dm->data;
1048*f885a11aSToby Isaac   if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine;
1049*f885a11aSToby Isaac   if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine;
1050bf9b5d84SToby Isaac   PetscFunctionReturn(0);
1051bf9b5d84SToby Isaac }
1052bf9b5d84SToby Isaac 
1053bf9b5d84SToby Isaac #undef __FUNCT__
1054c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetGradeFactor"
10559be51f97SToby Isaac /*@
10569be51f97SToby Isaac   DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to
10579be51f97SToby Isaac   indicate that the diameter of neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may
10589be51f97SToby Isaac   only support one particular choice of grading factor.
10599be51f97SToby Isaac 
10609be51f97SToby Isaac   Logically collective on dm
10619be51f97SToby Isaac 
10629be51f97SToby Isaac   Input Parameters:
10639be51f97SToby Isaac + dm - the forest
10649be51f97SToby Isaac - grade - the grading factor
10659be51f97SToby Isaac 
10669be51f97SToby Isaac   Level: advanced
10679be51f97SToby Isaac 
10689be51f97SToby Isaac .seealso: DMForestGetGradeFactor()
10699be51f97SToby Isaac @*/
1070c7eeac06SToby Isaac PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade)
1071c7eeac06SToby Isaac {
1072c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1073c7eeac06SToby Isaac 
1074c7eeac06SToby Isaac   PetscFunctionBegin;
1075c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1076ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup");
1077c7eeac06SToby Isaac   forest->gradeFactor = grade;
1078c7eeac06SToby Isaac   PetscFunctionReturn(0);
1079c7eeac06SToby Isaac }
1080c7eeac06SToby Isaac 
1081c7eeac06SToby Isaac #undef __FUNCT__
1082c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetGradeFactor"
10839be51f97SToby Isaac /*@
10849be51f97SToby Isaac   DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of
10859be51f97SToby Isaac   neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may only support one particular
10869be51f97SToby Isaac   choice of grading factor.
10879be51f97SToby Isaac 
10889be51f97SToby Isaac   Not collective
10899be51f97SToby Isaac 
10909be51f97SToby Isaac   Input Parameter:
10919be51f97SToby Isaac . dm - the forest
10929be51f97SToby Isaac 
10939be51f97SToby Isaac   Output Parameter:
10949be51f97SToby Isaac . grade - the grading factor
10959be51f97SToby Isaac 
10969be51f97SToby Isaac   Level: advanced
10979be51f97SToby Isaac 
10989be51f97SToby Isaac .seealso: DMForestSetGradeFactor()
10999be51f97SToby Isaac @*/
1100c7eeac06SToby Isaac PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade)
1101c7eeac06SToby Isaac {
1102c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1103c7eeac06SToby Isaac 
1104c7eeac06SToby Isaac   PetscFunctionBegin;
1105c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1106c7eeac06SToby Isaac   PetscValidIntPointer(grade,2);
1107c7eeac06SToby Isaac   *grade = forest->gradeFactor;
1108c7eeac06SToby Isaac   PetscFunctionReturn(0);
1109c7eeac06SToby Isaac }
1110c7eeac06SToby Isaac 
1111c7eeac06SToby Isaac #undef __FUNCT__
1112ef51cf95SToby Isaac #define __FUNCT__ "DMForestSetCellWeightFactor"
11139be51f97SToby Isaac /*@
11149be51f97SToby Isaac   DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes
11159be51f97SToby Isaac   the cell weight (see DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be
11169be51f97SToby Isaac   (cellWeight) * (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on
11179be51f97SToby Isaac   its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the
11189be51f97SToby Isaac   computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
11199be51f97SToby Isaac 
11209be51f97SToby Isaac   Logically collective on dm
11219be51f97SToby Isaac 
11229be51f97SToby Isaac   Input Parameters:
11239be51f97SToby Isaac + dm - the forest
11249be51f97SToby Isaac - weightsFactors - default 1.
11259be51f97SToby Isaac 
11269be51f97SToby Isaac   Level: advanced
11279be51f97SToby Isaac 
11289be51f97SToby Isaac .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights()
11299be51f97SToby Isaac @*/
1130ef51cf95SToby Isaac PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor)
1131c7eeac06SToby Isaac {
1132c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1133c7eeac06SToby Isaac 
1134c7eeac06SToby Isaac   PetscFunctionBegin;
1135c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1136ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup");
1137c7eeac06SToby Isaac   forest->weightsFactor = weightsFactor;
1138c7eeac06SToby Isaac   PetscFunctionReturn(0);
1139c7eeac06SToby Isaac }
1140c7eeac06SToby Isaac 
1141c7eeac06SToby Isaac #undef __FUNCT__
1142ef51cf95SToby Isaac #define __FUNCT__ "DMForestGetCellWeightFactor"
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 
1173c7eeac06SToby Isaac #undef __FUNCT__
1174c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellChart"
11759be51f97SToby Isaac /*@
11769be51f97SToby Isaac   DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process
11779be51f97SToby Isaac 
11789be51f97SToby Isaac   Not collective
11799be51f97SToby Isaac 
11809be51f97SToby Isaac   Input Parameter:
11819be51f97SToby Isaac . dm - the forest
11829be51f97SToby Isaac 
11839be51f97SToby Isaac   Output Parameters:
11849be51f97SToby Isaac + cStart - the first cell on this process
11859be51f97SToby Isaac - cEnd - one after the final cell on this process
11869be51f97SToby Isaac 
11871a244344SSatish Balay   Level: intermediate
11889be51f97SToby Isaac 
11899be51f97SToby Isaac .seealso: DMForestGetCellSF()
11909be51f97SToby Isaac @*/
1191c7eeac06SToby Isaac PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd)
1192c7eeac06SToby Isaac {
1193c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1194c7eeac06SToby Isaac   PetscErrorCode ierr;
1195c7eeac06SToby Isaac 
1196c7eeac06SToby Isaac   PetscFunctionBegin;
1197c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1198c7eeac06SToby Isaac   PetscValidIntPointer(cStart,2);
1199c7eeac06SToby Isaac   PetscValidIntPointer(cEnd,2);
1200c7eeac06SToby Isaac   if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) {
1201c7eeac06SToby Isaac     ierr = forest->createcellchart(dm,&forest->cStart,&forest->cEnd);CHKERRQ(ierr);
1202c7eeac06SToby Isaac   }
1203c7eeac06SToby Isaac   *cStart =  forest->cStart;
1204c7eeac06SToby Isaac   *cEnd   =  forest->cEnd;
1205c7eeac06SToby Isaac   PetscFunctionReturn(0);
1206c7eeac06SToby Isaac }
1207c7eeac06SToby Isaac 
1208c7eeac06SToby Isaac #undef __FUNCT__
1209c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellSF"
12109be51f97SToby Isaac /*@
12119be51f97SToby Isaac   DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes
12129be51f97SToby Isaac 
12139be51f97SToby Isaac   Not collective
12149be51f97SToby Isaac 
12159be51f97SToby Isaac   Input Parameter:
12169be51f97SToby Isaac . dm - the forest
12179be51f97SToby Isaac 
12189be51f97SToby Isaac   Output Parameter:
12199be51f97SToby Isaac . cellSF - the PetscSF
12209be51f97SToby Isaac 
12211a244344SSatish Balay   Level: intermediate
12229be51f97SToby Isaac 
12239be51f97SToby Isaac .seealso: DMForestGetCellChart()
12249be51f97SToby Isaac @*/
1225c7eeac06SToby Isaac PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF)
1226c7eeac06SToby Isaac {
1227c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1228c7eeac06SToby Isaac   PetscErrorCode ierr;
1229c7eeac06SToby Isaac 
1230c7eeac06SToby Isaac   PetscFunctionBegin;
1231c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1232c7eeac06SToby Isaac   PetscValidPointer(cellSF,2);
1233c7eeac06SToby Isaac   if ((!forest->cellSF) && forest->createcellsf) {
1234c7eeac06SToby Isaac     ierr = forest->createcellsf(dm,&forest->cellSF);CHKERRQ(ierr);
1235c7eeac06SToby Isaac   }
1236c7eeac06SToby Isaac   *cellSF = forest->cellSF;
1237c7eeac06SToby Isaac   PetscFunctionReturn(0);
1238c7eeac06SToby Isaac }
1239c7eeac06SToby Isaac 
1240c7eeac06SToby Isaac #undef __FUNCT__
1241ebdf65a2SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityLabel"
12429be51f97SToby Isaac /*@C
12439be51f97SToby Isaac   DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see
12449be51f97SToby Isaac   DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination).  The
12459be51f97SToby Isaac   interpretation of the label values is up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and
12469be51f97SToby Isaac   DM_FOREST_COARSEN have been reserved as choices that should be accepted by all subtypes.
12479be51f97SToby Isaac 
12489be51f97SToby Isaac   Logically collective on dm
12499be51f97SToby Isaac 
12509be51f97SToby Isaac   Input Parameters:
12519be51f97SToby Isaac - dm - the forest
12529be51f97SToby Isaac + adaptLabel - the name of the label in the pre-adaptation forest
12539be51f97SToby Isaac 
12549be51f97SToby Isaac   Level: intermediate
12559be51f97SToby Isaac 
12569be51f97SToby Isaac .seealso DMForestGetAdaptivityLabel()
12579be51f97SToby Isaac @*/
1258ebdf65a2SToby Isaac PetscErrorCode DMForestSetAdaptivityLabel(DM dm, const char * adaptLabel)
1259c7eeac06SToby Isaac {
1260c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1261c7eeac06SToby Isaac   PetscErrorCode ierr;
1262c7eeac06SToby Isaac 
1263c7eeac06SToby Isaac   PetscFunctionBegin;
1264c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1265ebdf65a2SToby Isaac   ierr = PetscFree(forest->adaptLabel);CHKERRQ(ierr);
1266ebdf65a2SToby Isaac   ierr = PetscStrallocpy(adaptLabel,&forest->adaptLabel);CHKERRQ(ierr);
1267c7eeac06SToby Isaac   PetscFunctionReturn(0);
1268c7eeac06SToby Isaac }
1269c7eeac06SToby Isaac 
1270c7eeac06SToby Isaac #undef __FUNCT__
1271ebdf65a2SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityLabel"
12729be51f97SToby Isaac /*@C
12739be51f97SToby Isaac   DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that
12749be51f97SToby Isaac   holds the adaptation flags (refinement, coarsening, or some combination).  The interpretation of the label values is
12759be51f97SToby Isaac   up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and DM_FOREST_COARSEN have been reserved as
12769be51f97SToby Isaac   choices that should be accepted by all subtypes.
12779be51f97SToby Isaac 
12789be51f97SToby Isaac   Not collective
12799be51f97SToby Isaac 
12809be51f97SToby Isaac   Input Parameter:
12819be51f97SToby Isaac . dm - the forest
12829be51f97SToby Isaac 
12839be51f97SToby Isaac   Output Parameter:
12849be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest
12859be51f97SToby Isaac 
12869be51f97SToby Isaac   Level: intermediate
12879be51f97SToby Isaac 
12889be51f97SToby Isaac .seealso DMForestSetAdaptivityLabel()
12899be51f97SToby Isaac @*/
1290ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityLabel(DM dm, const char ** adaptLabel)
1291c7eeac06SToby Isaac {
1292c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1293c7eeac06SToby Isaac 
1294c7eeac06SToby Isaac   PetscFunctionBegin;
1295c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1296ba936b91SToby Isaac   *adaptLabel = forest->adaptLabel;
1297c7eeac06SToby Isaac   PetscFunctionReturn(0);
1298c7eeac06SToby Isaac }
1299c7eeac06SToby Isaac 
1300c7eeac06SToby Isaac #undef __FUNCT__
1301c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetCellWeights"
13029be51f97SToby Isaac /*@
13039be51f97SToby Isaac   DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
13049be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
13059be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
13069be51f97SToby Isaac   of 1.
13079be51f97SToby Isaac 
13089be51f97SToby Isaac   Logically collective on dm
13099be51f97SToby Isaac 
13109be51f97SToby Isaac   Input Parameters:
13119be51f97SToby Isaac + dm - the forest
13129be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
13139be51f97SToby Isaac - copyMode - how weights should reference weights
13149be51f97SToby Isaac 
13159be51f97SToby Isaac   Level: advanced
13169be51f97SToby Isaac 
13179be51f97SToby Isaac .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity()
13189be51f97SToby Isaac @*/
1319c7eeac06SToby Isaac PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode)
1320c7eeac06SToby Isaac {
1321c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1322c7eeac06SToby Isaac   PetscInt       cStart, cEnd;
1323c7eeac06SToby Isaac   PetscErrorCode ierr;
1324c7eeac06SToby Isaac 
1325c7eeac06SToby Isaac   PetscFunctionBegin;
1326c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1327c7eeac06SToby Isaac   ierr = DMForestGetCellChart(dm,&cStart,&cEnd);CHKERRQ(ierr);
1328c7eeac06SToby Isaac   if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd);
1329c7eeac06SToby Isaac   if (copyMode == PETSC_COPY_VALUES) {
1330c7eeac06SToby Isaac     if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) {
1331c7eeac06SToby Isaac       ierr = PetscMalloc1(cEnd-cStart,&forest->cellWeights);CHKERRQ(ierr);
1332c7eeac06SToby Isaac     }
1333c7eeac06SToby Isaac     ierr                        = PetscMemcpy(forest->cellWeights,weights,(cEnd-cStart)*sizeof(*weights));CHKERRQ(ierr);
1334c7eeac06SToby Isaac     forest->cellWeightsCopyMode = PETSC_OWN_POINTER;
1335c7eeac06SToby Isaac     PetscFunctionReturn(0);
1336c7eeac06SToby Isaac   }
1337c7eeac06SToby Isaac   if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) {
1338c7eeac06SToby Isaac     ierr = PetscFree(forest->cellWeights);CHKERRQ(ierr);
1339c7eeac06SToby Isaac   }
1340c7eeac06SToby Isaac   forest->cellWeights         = weights;
1341c7eeac06SToby Isaac   forest->cellWeightsCopyMode = copyMode;
1342c7eeac06SToby Isaac   PetscFunctionReturn(0);
1343c7eeac06SToby Isaac }
1344c7eeac06SToby Isaac 
1345c7eeac06SToby Isaac #undef __FUNCT__
1346c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellWeights"
13479be51f97SToby Isaac /*@
13489be51f97SToby Isaac   DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
13499be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
13509be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
13519be51f97SToby Isaac   of 1.
13529be51f97SToby Isaac 
13539be51f97SToby Isaac   Not collective
13549be51f97SToby Isaac 
13559be51f97SToby Isaac   Input Parameter:
13569be51f97SToby Isaac . dm - the forest
13579be51f97SToby Isaac 
13589be51f97SToby Isaac   Output Parameter:
13599be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
13609be51f97SToby Isaac 
13619be51f97SToby Isaac   Level: advanced
13629be51f97SToby Isaac 
13639be51f97SToby Isaac .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity()
13649be51f97SToby Isaac @*/
1365c7eeac06SToby Isaac PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights)
1366c7eeac06SToby Isaac {
1367c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1368c7eeac06SToby Isaac 
1369c7eeac06SToby Isaac   PetscFunctionBegin;
1370c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1371c7eeac06SToby Isaac   PetscValidPointer(weights,2);
1372c7eeac06SToby Isaac   *weights = forest->cellWeights;
1373c7eeac06SToby Isaac   PetscFunctionReturn(0);
1374c7eeac06SToby Isaac }
1375c7eeac06SToby Isaac 
1376c7eeac06SToby Isaac #undef __FUNCT__
1377c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetWeightCapacity"
13789be51f97SToby Isaac /*@
13799be51f97SToby Isaac   DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning
13809be51f97SToby Isaac   a pre-adaptation forest (see DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each
13819be51f97SToby Isaac   process's cells to the process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that
13829be51f97SToby Isaac   the current process should not have any cells after repartitioning.
13839be51f97SToby Isaac 
13849be51f97SToby Isaac   Logically Collective on dm
13859be51f97SToby Isaac 
13869be51f97SToby Isaac   Input parameters:
13879be51f97SToby Isaac + dm - the forest
13889be51f97SToby Isaac - capacity - this process's capacity
13899be51f97SToby Isaac 
13909be51f97SToby Isaac   Level: advanced
13919be51f97SToby Isaac 
13929be51f97SToby Isaac .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
13939be51f97SToby Isaac @*/
1394c7eeac06SToby Isaac PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity)
1395c7eeac06SToby Isaac {
1396c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1397c7eeac06SToby Isaac 
1398c7eeac06SToby Isaac   PetscFunctionBegin;
1399c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1400ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup");
1401c7eeac06SToby Isaac   if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity);
1402c7eeac06SToby Isaac   forest->weightCapacity = capacity;
1403c7eeac06SToby Isaac   PetscFunctionReturn(0);
1404c7eeac06SToby Isaac }
1405c7eeac06SToby Isaac 
1406c7eeac06SToby Isaac #undef __FUNCT__
1407c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetWeightCapacity"
14089be51f97SToby Isaac /*@
14099be51f97SToby Isaac   DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see
14109be51f97SToby Isaac   DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each process's cells to the
14119be51f97SToby Isaac   process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that the current process
14129be51f97SToby Isaac   should not have any cells after repartitioning.
14139be51f97SToby Isaac 
14149be51f97SToby Isaac   Not collective
14159be51f97SToby Isaac 
14169be51f97SToby Isaac   Input parameter:
14179be51f97SToby Isaac . dm - the forest
14189be51f97SToby Isaac 
14199be51f97SToby Isaac   Output parameter:
14209be51f97SToby Isaac . capacity - this process's capacity
14219be51f97SToby Isaac 
14229be51f97SToby Isaac   Level: advanced
14239be51f97SToby Isaac 
14249be51f97SToby Isaac .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
14259be51f97SToby Isaac @*/
1426c7eeac06SToby Isaac PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity)
1427c7eeac06SToby Isaac {
1428c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1429c7eeac06SToby Isaac 
1430c7eeac06SToby Isaac   PetscFunctionBegin;
1431c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1432c7eeac06SToby Isaac   PetscValidRealPointer(capacity,2);
1433c7eeac06SToby Isaac   *capacity = forest->weightCapacity;
1434c7eeac06SToby Isaac   PetscFunctionReturn(0);
1435c7eeac06SToby Isaac }
1436c7eeac06SToby Isaac 
1437dd8e54a2SToby Isaac #undef __FUNCT__
1438db4d5e8cSToby Isaac #define __FUNCT__ "DMSetFromOptions_Forest"
14390709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm)
1440db4d5e8cSToby Isaac {
1441db4d5e8cSToby Isaac   DM_Forest                  *forest = (DM_Forest*) dm->data;
144256ba9f64SToby Isaac   PetscBool                  flg, flg1, flg2, flg3, flg4;
1443dd8e54a2SToby Isaac   DMForestTopology           oldTopo;
1444c7eeac06SToby Isaac   char                       stringBuffer[256];
1445dd8e54a2SToby Isaac   PetscViewer                viewer;
1446dd8e54a2SToby Isaac   PetscViewerFormat          format;
144756ba9f64SToby Isaac   PetscInt                   adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade;
1448c7eeac06SToby Isaac   PetscReal                  weightsFactor;
1449c7eeac06SToby Isaac   DMForestAdaptivityStrategy adaptStrategy;
1450db4d5e8cSToby Isaac   PetscErrorCode             ierr;
1451db4d5e8cSToby Isaac 
1452db4d5e8cSToby Isaac   PetscFunctionBegin;
1453db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
145458762b62SToby Isaac   forest->setfromoptionscalled = PETSC_TRUE;
1455dd8e54a2SToby Isaac   ierr                         = DMForestGetTopology(dm, &oldTopo);CHKERRQ(ierr);
1456a3eda1eaSToby Isaac   ierr                         = PetscOptionsHead(PetscOptionsObject,"DMForest Options");CHKERRQ(ierr);
145756ba9f64SToby Isaac   ierr                         = PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);CHKERRQ(ierr);
145856ba9f64SToby Isaac   ierr                         = PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);CHKERRQ(ierr);
145956ba9f64SToby Isaac   ierr                         = PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);CHKERRQ(ierr);
146056ba9f64SToby Isaac   ierr                         = PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);CHKERRQ(ierr);
1461*f885a11aSToby 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}");
146256ba9f64SToby Isaac   if (flg1) {
146356ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,(DMForestTopology)stringBuffer);CHKERRQ(ierr);
146456ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
146520e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
146656ba9f64SToby Isaac   }
146756ba9f64SToby Isaac   if (flg2) {
1468dd8e54a2SToby Isaac     DM base;
1469dd8e54a2SToby Isaac 
1470dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&base);CHKERRQ(ierr);
1471dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1472dd8e54a2SToby Isaac     ierr = DMLoad(base,viewer);CHKERRQ(ierr);
1473dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1474dd8e54a2SToby Isaac     ierr = DMForestSetBaseDM(dm,base);CHKERRQ(ierr);
1475dd8e54a2SToby Isaac     ierr = DMDestroy(&base);CHKERRQ(ierr);
147656ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
147720e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
1478dd8e54a2SToby Isaac   }
147956ba9f64SToby Isaac   if (flg3) {
1480dd8e54a2SToby Isaac     DM coarse;
1481dd8e54a2SToby Isaac 
1482dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&coarse);CHKERRQ(ierr);
1483dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1484dd8e54a2SToby Isaac     ierr = DMLoad(coarse,viewer);CHKERRQ(ierr);
1485dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
148620e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,coarse);CHKERRQ(ierr);
1487dd8e54a2SToby Isaac     ierr = DMDestroy(&coarse);CHKERRQ(ierr);
148856ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
148956ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1490dd8e54a2SToby Isaac   }
149156ba9f64SToby Isaac   if (flg4) {
1492dd8e54a2SToby Isaac     DM fine;
1493dd8e54a2SToby Isaac 
1494dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&fine);CHKERRQ(ierr);
1495dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1496dd8e54a2SToby Isaac     ierr = DMLoad(fine,viewer);CHKERRQ(ierr);
1497dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
149820e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,fine);CHKERRQ(ierr);
1499dd8e54a2SToby Isaac     ierr = DMDestroy(&fine);CHKERRQ(ierr);
150056ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
150156ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1502dd8e54a2SToby Isaac   }
1503dd8e54a2SToby Isaac   ierr = DMForestGetAdjacencyDimension(dm,&adjDim);CHKERRQ(ierr);
1504dd8e54a2SToby Isaac   ierr = PetscOptionsInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg);CHKERRQ(ierr);
1505dd8e54a2SToby Isaac   if (flg) {
1506dd8e54a2SToby Isaac     ierr = DMForestSetAdjacencyDimension(dm,adjDim);CHKERRQ(ierr);
1507*f885a11aSToby Isaac   } else {
1508dd8e54a2SToby Isaac     ierr = DMForestGetAdjacencyCodimension(dm,&adjCodim);CHKERRQ(ierr);
1509dd8e54a2SToby Isaac     ierr = PetscOptionsInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg);CHKERRQ(ierr);
1510dd8e54a2SToby Isaac     if (flg) {
1511dd8e54a2SToby Isaac       ierr = DMForestSetAdjacencyCodimension(dm,adjCodim);CHKERRQ(ierr);
1512dd8e54a2SToby Isaac     }
1513dd8e54a2SToby Isaac   }
1514dd8e54a2SToby Isaac   ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr);
1515dd8e54a2SToby Isaac   ierr = PetscOptionsInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg);CHKERRQ(ierr);
1516dd8e54a2SToby Isaac   if (flg) {
1517dd8e54a2SToby Isaac     ierr = DMForestSetPartitionOverlap(dm,overlap);CHKERRQ(ierr);
1518dd8e54a2SToby Isaac   }
1519a6121fbdSMatthew G. Knepley #if 0
1520a6121fbdSMatthew 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);
1521a6121fbdSMatthew G. Knepley   if (flg) {
1522a6121fbdSMatthew G. Knepley     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1523a6121fbdSMatthew G. Knepley     ierr = DMForestSetInitialRefinement(dm,minRefinement);CHKERRQ(ierr);
1524a6121fbdSMatthew G. Knepley   }
1525a6121fbdSMatthew 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);
1526a6121fbdSMatthew G. Knepley   if (flg) {
1527a6121fbdSMatthew G. Knepley     ierr = DMForestSetMinimumRefinement(dm,0);CHKERRQ(ierr);
1528a6121fbdSMatthew G. Knepley     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
1529a6121fbdSMatthew G. Knepley   }
1530a6121fbdSMatthew G. Knepley #endif
1531dd8e54a2SToby Isaac   ierr = DMForestGetMinimumRefinement(dm,&minRefinement);CHKERRQ(ierr);
1532dd8e54a2SToby Isaac   ierr = PetscOptionsInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg);CHKERRQ(ierr);
1533dd8e54a2SToby Isaac   if (flg) {
1534dd8e54a2SToby Isaac     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1535db4d5e8cSToby Isaac   }
153656ba9f64SToby Isaac   ierr = DMForestGetInitialRefinement(dm,&initRefinement);CHKERRQ(ierr);
153756ba9f64SToby Isaac   ierr = PetscOptionsInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg);CHKERRQ(ierr);
153856ba9f64SToby Isaac   if (flg) {
153956ba9f64SToby Isaac     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
154056ba9f64SToby Isaac   }
1541c7eeac06SToby Isaac   ierr = DMForestGetMaximumRefinement(dm,&maxRefinement);CHKERRQ(ierr);
1542c7eeac06SToby Isaac   ierr = PetscOptionsInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg);CHKERRQ(ierr);
1543c7eeac06SToby Isaac   if (flg) {
1544c7eeac06SToby Isaac     ierr = DMForestSetMaximumRefinement(dm,maxRefinement);CHKERRQ(ierr);
1545c7eeac06SToby Isaac   }
1546c7eeac06SToby Isaac   ierr = DMForestGetAdaptivityStrategy(dm,&adaptStrategy);CHKERRQ(ierr);
1547c7eeac06SToby Isaac   ierr = PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,256,&flg);CHKERRQ(ierr);
1548c7eeac06SToby Isaac   if (flg) {
1549c7eeac06SToby Isaac     ierr = DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);CHKERRQ(ierr);
1550c7eeac06SToby Isaac   }
1551c7eeac06SToby Isaac   ierr = DMForestGetGradeFactor(dm,&grade);CHKERRQ(ierr);
1552c7eeac06SToby Isaac   ierr = PetscOptionsInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg);CHKERRQ(ierr);
1553c7eeac06SToby Isaac   if (flg) {
1554c7eeac06SToby Isaac     ierr = DMForestSetGradeFactor(dm,grade);CHKERRQ(ierr);
1555c7eeac06SToby Isaac   }
1556c7eeac06SToby Isaac   ierr = DMForestGetCellWeightFactor(dm,&weightsFactor);CHKERRQ(ierr);
1557c7eeac06SToby Isaac   ierr = PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);CHKERRQ(ierr);
1558c7eeac06SToby Isaac   if (flg) {
1559c7eeac06SToby Isaac     ierr = DMForestSetCellWeightFactor(dm,weightsFactor);CHKERRQ(ierr);
1560c7eeac06SToby Isaac   }
1561db4d5e8cSToby Isaac   ierr = PetscOptionsTail();CHKERRQ(ierr);
1562db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1563db4d5e8cSToby Isaac }
1564db4d5e8cSToby Isaac 
1565db4d5e8cSToby Isaac #undef __FUNCT__
1566d8984e3bSMatthew G. Knepley #define __FUNCT__ "DMCreateSubDM_Forest"
1567d8984e3bSMatthew G. Knepley PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
1568d8984e3bSMatthew G. Knepley {
1569d8984e3bSMatthew G. Knepley   PetscErrorCode ierr;
1570d8984e3bSMatthew G. Knepley 
1571d8984e3bSMatthew G. Knepley   PetscFunctionBegin;
1572d8984e3bSMatthew G. Knepley   if (subdm) {ierr = DMClone(dm, subdm);CHKERRQ(ierr);}
1573d8984e3bSMatthew G. Knepley   ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr);
1574d8984e3bSMatthew G. Knepley   PetscFunctionReturn(0);
1575d8984e3bSMatthew G. Knepley }
1576d8984e3bSMatthew G. Knepley 
1577d8984e3bSMatthew G. Knepley #undef __FUNCT__
15785421bac9SToby Isaac #define __FUNCT__ "DMRefine_Forest"
15795421bac9SToby Isaac PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
15805421bac9SToby Isaac {
15815421bac9SToby Isaac   DMLabel        refine;
15825421bac9SToby Isaac   DM             fineDM;
15835421bac9SToby Isaac   PetscErrorCode ierr;
15845421bac9SToby Isaac 
15855421bac9SToby Isaac   PetscFunctionBegin;
15865421bac9SToby Isaac   ierr = DMGetFineDM(dm,&fineDM);CHKERRQ(ierr);
15875421bac9SToby Isaac   if (fineDM) {
15885421bac9SToby Isaac     ierr       = PetscObjectReference((PetscObject)fineDM);CHKERRQ(ierr);
15895421bac9SToby Isaac     *dmRefined = fineDM;
15905421bac9SToby Isaac     PetscFunctionReturn(0);
15915421bac9SToby Isaac   }
15925421bac9SToby Isaac   ierr = DMForestTemplate(dm,comm,dmRefined);CHKERRQ(ierr);
15935421bac9SToby Isaac   ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr);
15945421bac9SToby Isaac   if (!refine) {
15955421bac9SToby Isaac     ierr = DMCreateLabel(dm,"refine");CHKERRQ(ierr);
15965421bac9SToby Isaac     ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr);
15975421bac9SToby Isaac     ierr = DMLabelSetDefaultValue(refine,DM_FOREST_REFINE);CHKERRQ(ierr);
15985421bac9SToby Isaac   }
15995421bac9SToby Isaac   ierr = DMForestSetAdaptivityLabel(*dmRefined,"refine");CHKERRQ(ierr);
16005421bac9SToby Isaac   PetscFunctionReturn(0);
16015421bac9SToby Isaac }
16025421bac9SToby Isaac 
16035421bac9SToby Isaac #undef __FUNCT__
16045421bac9SToby Isaac #define __FUNCT__ "DMCoarsen_Forest"
16055421bac9SToby Isaac PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
16065421bac9SToby Isaac {
16075421bac9SToby Isaac   DMLabel        coarsen;
16085421bac9SToby Isaac   DM             coarseDM;
16095421bac9SToby Isaac   PetscErrorCode ierr;
16105421bac9SToby Isaac 
16115421bac9SToby Isaac   PetscFunctionBegin;
16124098eed7SToby Isaac   {
16134098eed7SToby Isaac     PetscMPIInt mpiComparison;
16144098eed7SToby Isaac     MPI_Comm    dmcomm = PetscObjectComm((PetscObject)dm);
16154098eed7SToby Isaac 
16164098eed7SToby Isaac     ierr = MPI_Comm_compare(comm, dmcomm, &mpiComparison);CHKERRQ(ierr);
1617*f885a11aSToby Isaac     if (mpiComparison != MPI_IDENT && mpiComparison != MPI_CONGRUENT) SETERRQ(dmcomm,PETSC_ERR_SUP,"No support for different communicators yet");
16184098eed7SToby Isaac   }
16195421bac9SToby Isaac   ierr = DMGetCoarseDM(dm,&coarseDM);CHKERRQ(ierr);
16205421bac9SToby Isaac   if (coarseDM) {
16215421bac9SToby Isaac     ierr         = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr);
16225421bac9SToby Isaac     *dmCoarsened = coarseDM;
16235421bac9SToby Isaac     PetscFunctionReturn(0);
16245421bac9SToby Isaac   }
16255421bac9SToby Isaac   ierr = DMForestTemplate(dm,comm,dmCoarsened);CHKERRQ(ierr);
16264098eed7SToby Isaac   ierr = DMForestSetAdaptivityPurpose(coarseDM,DM_FOREST_COARSEN);CHKERRQ(ierr);
16275421bac9SToby Isaac   ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr);
16285421bac9SToby Isaac   if (!coarsen) {
16295421bac9SToby Isaac     ierr = DMCreateLabel(dm,"coarsen");CHKERRQ(ierr);
16305421bac9SToby Isaac     ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr);
16315421bac9SToby Isaac     ierr = DMLabelSetDefaultValue(coarsen,DM_FOREST_COARSEN);CHKERRQ(ierr);
16325421bac9SToby Isaac   }
16335421bac9SToby Isaac   ierr = DMForestSetAdaptivityLabel(*dmCoarsened,"coarsen");CHKERRQ(ierr);
16345421bac9SToby Isaac   PetscFunctionReturn(0);
16355421bac9SToby Isaac }
16365421bac9SToby Isaac 
16375421bac9SToby Isaac #undef __FUNCT__
1638d222f98bSToby Isaac #define __FUNCT__ "DMInitialize_Forest"
1639d222f98bSToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm)
1640d222f98bSToby Isaac {
1641d222f98bSToby Isaac   PetscErrorCode ierr;
1642d222f98bSToby Isaac 
1643d222f98bSToby Isaac   PetscFunctionBegin;
1644d222f98bSToby Isaac   ierr = PetscMemzero(dm->ops,sizeof(*(dm->ops)));CHKERRQ(ierr);
1645d222f98bSToby Isaac 
1646d222f98bSToby Isaac   dm->ops->clone          = DMClone_Forest;
1647d222f98bSToby Isaac   dm->ops->setfromoptions = DMSetFromOptions_Forest;
1648d222f98bSToby Isaac   dm->ops->destroy        = DMDestroy_Forest;
1649d8984e3bSMatthew G. Knepley   dm->ops->createsubdm    = DMCreateSubDM_Forest;
16505421bac9SToby Isaac   dm->ops->refine         = DMRefine_Forest;
16515421bac9SToby Isaac   dm->ops->coarsen        = DMCoarsen_Forest;
1652d222f98bSToby Isaac   PetscFunctionReturn(0);
1653d222f98bSToby Isaac }
1654d222f98bSToby Isaac 
16559be51f97SToby Isaac /*MC
16569be51f97SToby Isaac   DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh.  Forests usually have a base DM (see DMForestGetBaseDM()), from which it is refined.  The refinement and partitioning of forests is considered immutable after DMSetUp() is called.  To adapt a mesh, one should call DMForestTemplate() to create a new mesh that will default to being identical to it, specify how that mesh should differ, and then calling DMSetUp() on the new mesh.
16579be51f97SToby Isaac 
16589be51f97SToby Isaac   To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the previous mesh whose values indicate which cells should be refined (DM_FOREST_REFINE) or coarsened (DM_FOREST_COARSEN) and how (subtypes are free to allow additional values for things like anisotropic refinement).  The name of the label should be given to the *new* mesh with DMForestSetAdaptivityLabel().
16599be51f97SToby Isaac 
16609be51f97SToby Isaac   Level: advanced
16619be51f97SToby Isaac 
16629be51f97SToby Isaac .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel()
16639be51f97SToby Isaac M*/
16649be51f97SToby Isaac 
1665d222f98bSToby Isaac #undef __FUNCT__
1666db4d5e8cSToby Isaac #define __FUNCT__ "DMCreate_Forest"
1667db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1668db4d5e8cSToby Isaac {
1669db4d5e8cSToby Isaac   DM_Forest      *forest;
1670db4d5e8cSToby Isaac   PetscErrorCode ierr;
1671db4d5e8cSToby Isaac 
1672db4d5e8cSToby Isaac   PetscFunctionBegin;
1673db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1674db4d5e8cSToby Isaac   ierr                         = PetscNewLog(dm,&forest);CHKERRQ(ierr);
1675db4d5e8cSToby Isaac   dm->dim                      = 0;
1676db4d5e8cSToby Isaac   dm->data                     = forest;
1677db4d5e8cSToby Isaac   forest->refct                = 1;
1678db4d5e8cSToby Isaac   forest->data                 = NULL;
167958762b62SToby Isaac   forest->setfromoptionscalled = PETSC_FALSE;
1680db4d5e8cSToby Isaac   forest->topology             = NULL;
1681db4d5e8cSToby Isaac   forest->base                 = NULL;
1682db4d5e8cSToby Isaac   forest->adjDim               = PETSC_DEFAULT;
1683db4d5e8cSToby Isaac   forest->overlap              = PETSC_DEFAULT;
1684db4d5e8cSToby Isaac   forest->minRefinement        = PETSC_DEFAULT;
1685db4d5e8cSToby Isaac   forest->maxRefinement        = PETSC_DEFAULT;
168656ba9f64SToby Isaac   forest->initRefinement       = PETSC_DEFAULT;
1687c7eeac06SToby Isaac   forest->cStart               = PETSC_DETERMINE;
1688c7eeac06SToby Isaac   forest->cEnd                 = PETSC_DETERMINE;
1689db4d5e8cSToby Isaac   forest->cellSF               = 0;
1690ebdf65a2SToby Isaac   forest->adaptLabel           = NULL;
1691db4d5e8cSToby Isaac   forest->gradeFactor          = 2;
1692db4d5e8cSToby Isaac   forest->cellWeights          = NULL;
1693db4d5e8cSToby Isaac   forest->cellWeightsCopyMode  = PETSC_USE_POINTER;
1694db4d5e8cSToby Isaac   forest->weightsFactor        = 1.;
1695db4d5e8cSToby Isaac   forest->weightCapacity       = 1.;
1696a73e2921SToby Isaac   ierr                         = DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);CHKERRQ(ierr);
1697d222f98bSToby Isaac   ierr                         = DMInitialize_Forest(dm);CHKERRQ(ierr);
1698db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1699db4d5e8cSToby Isaac }
1700db4d5e8cSToby Isaac 
1701