xref: /petsc/src/dm/impls/forest/forest.c (revision 28dfcf7c3f543190bef3ee99789beae6d38a1aec)
19be51f97SToby Isaac #include <petsc/private/dmforestimpl.h> /*I "petscdmforest.h" I*/
29be51f97SToby Isaac #include <petsc/private/dmimpl.h>       /*I "petscdm.h" I*/
3a1b0c543SToby Isaac #include <petsc/private/dmlabelimpl.h>  /*I "petscdmlabel.h" I*/
4ef19d27cSToby Isaac #include <petscsf.h>
5db4d5e8cSToby Isaac 
627d4645fSToby Isaac PetscBool DMForestPackageInitialized = PETSC_FALSE;
727d4645fSToby Isaac 
827d4645fSToby Isaac typedef struct _DMForestTypeLink*DMForestTypeLink;
927d4645fSToby Isaac 
1027d4645fSToby Isaac struct _DMForestTypeLink
1127d4645fSToby Isaac {
1227d4645fSToby Isaac   char             *name;
1327d4645fSToby Isaac   DMForestTypeLink next;
1427d4645fSToby Isaac };
1527d4645fSToby Isaac 
1627d4645fSToby Isaac DMForestTypeLink DMForestTypeList;
1727d4645fSToby Isaac 
1827d4645fSToby Isaac static PetscErrorCode DMForestPackageFinalize(void)
1927d4645fSToby Isaac {
2027d4645fSToby Isaac   DMForestTypeLink oldLink, link = DMForestTypeList;
2127d4645fSToby Isaac   PetscErrorCode   ierr;
2227d4645fSToby Isaac 
2327d4645fSToby Isaac   PetscFunctionBegin;
2427d4645fSToby Isaac   while (link) {
2527d4645fSToby Isaac     oldLink = link;
26f416af30SBarry Smith     ierr    = PetscFree(oldLink->name);CHKERRQ(ierr);
2727d4645fSToby Isaac     link    = oldLink->next;
2827d4645fSToby Isaac     ierr    = PetscFree(oldLink);CHKERRQ(ierr);
2927d4645fSToby Isaac   }
3027d4645fSToby Isaac   PetscFunctionReturn(0);
3127d4645fSToby Isaac }
3227d4645fSToby Isaac 
3327d4645fSToby Isaac static PetscErrorCode DMForestPackageInitialize(void)
3427d4645fSToby Isaac {
3527d4645fSToby Isaac   PetscErrorCode ierr;
3627d4645fSToby Isaac 
3727d4645fSToby Isaac   PetscFunctionBegin;
3827d4645fSToby Isaac   if (DMForestPackageInitialized) PetscFunctionReturn(0);
3927d4645fSToby Isaac   DMForestPackageInitialized = PETSC_TRUE;
40f885a11aSToby Isaac 
4127d4645fSToby Isaac   ierr = DMForestRegisterType(DMFOREST);CHKERRQ(ierr);
4227d4645fSToby Isaac   ierr = PetscRegisterFinalize(DMForestPackageFinalize);CHKERRQ(ierr);
4327d4645fSToby Isaac   PetscFunctionReturn(0);
4427d4645fSToby Isaac }
4527d4645fSToby Isaac 
469be51f97SToby Isaac /*@C
479be51f97SToby Isaac   DMForestRegisterType - Registers a DMType as a subtype of DMFOREST (so that DMIsForest() will be correct)
489be51f97SToby Isaac 
499be51f97SToby Isaac   Not Collective
509be51f97SToby Isaac 
519be51f97SToby Isaac   Input parameter:
529be51f97SToby Isaac . name - the name of the type
539be51f97SToby Isaac 
549be51f97SToby Isaac   Level: advanced
559be51f97SToby Isaac 
569be51f97SToby Isaac .seealso: DMFOREST, DMIsForest()
579be51f97SToby Isaac @*/
5827d4645fSToby Isaac PetscErrorCode DMForestRegisterType(DMType name)
5927d4645fSToby Isaac {
6027d4645fSToby Isaac   DMForestTypeLink link;
6127d4645fSToby Isaac   PetscErrorCode   ierr;
6227d4645fSToby Isaac 
6327d4645fSToby Isaac   PetscFunctionBegin;
6427d4645fSToby Isaac   ierr             = DMForestPackageInitialize();CHKERRQ(ierr);
6527d4645fSToby Isaac   ierr             = PetscNew(&link);CHKERRQ(ierr);
6627d4645fSToby Isaac   ierr             = PetscStrallocpy(name,&link->name);CHKERRQ(ierr);
6727d4645fSToby Isaac   link->next       = DMForestTypeList;
6827d4645fSToby Isaac   DMForestTypeList = link;
6927d4645fSToby Isaac   PetscFunctionReturn(0);
7027d4645fSToby Isaac }
7127d4645fSToby Isaac 
729be51f97SToby Isaac /*@
739be51f97SToby Isaac   DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes
749be51f97SToby Isaac 
759be51f97SToby Isaac   Not Collective
769be51f97SToby Isaac 
779be51f97SToby Isaac   Input parameter:
789be51f97SToby Isaac . dm - the DM object
799be51f97SToby Isaac 
809be51f97SToby Isaac   Output parameter:
819be51f97SToby Isaac . isForest - whether dm is a subtype of DMFOREST
829be51f97SToby Isaac 
839be51f97SToby Isaac   Level: intermediate
849be51f97SToby Isaac 
859be51f97SToby Isaac .seealso: DMFOREST, DMForestRegisterType()
869be51f97SToby Isaac @*/
8727d4645fSToby Isaac PetscErrorCode DMIsForest(DM dm, PetscBool *isForest)
8827d4645fSToby Isaac {
8927d4645fSToby Isaac   DMForestTypeLink link = DMForestTypeList;
9027d4645fSToby Isaac   PetscErrorCode   ierr;
9127d4645fSToby Isaac 
9227d4645fSToby Isaac   PetscFunctionBegin;
9327d4645fSToby Isaac   while (link) {
9427d4645fSToby Isaac     PetscBool sameType;
9527d4645fSToby Isaac     ierr = PetscObjectTypeCompare((PetscObject)dm,link->name,&sameType);CHKERRQ(ierr);
9627d4645fSToby Isaac     if (sameType) {
9727d4645fSToby Isaac       *isForest = PETSC_TRUE;
9827d4645fSToby Isaac       PetscFunctionReturn(0);
9927d4645fSToby Isaac     }
10027d4645fSToby Isaac     link = link->next;
10127d4645fSToby Isaac   }
10227d4645fSToby Isaac   *isForest = PETSC_FALSE;
10327d4645fSToby Isaac   PetscFunctionReturn(0);
10427d4645fSToby Isaac }
10527d4645fSToby Isaac 
1069be51f97SToby Isaac /*@
1079be51f97SToby Isaac   DMForestTemplate - Create a new DM that will be adapted from a source DM.  The new DM reproduces the configuration
1089be51f97SToby Isaac   of the source, but is not yet setup, so that the user can then define only the ways that the new DM should differ
1099be51f97SToby Isaac   (by, e.g., refinement or repartitioning).  The source DM is also set as the adaptivity source DM of the new DM (see
1109be51f97SToby Isaac   DMForestSetAdaptivityForest()).
1119be51f97SToby Isaac 
1129be51f97SToby Isaac   Collective on dm
1139be51f97SToby Isaac 
1149be51f97SToby Isaac   Input Parameters:
1159be51f97SToby Isaac + dm - the source DM object
1169be51f97SToby Isaac - comm - the communicator for the new DM (this communicator is currently ignored, but is present so that DMForestTemplate() can be used within DMCoarsen())
1179be51f97SToby Isaac 
1189be51f97SToby Isaac   Output Parameter:
1199be51f97SToby Isaac . tdm - the new DM object
1209be51f97SToby Isaac 
1219be51f97SToby Isaac   Level: intermediate
1229be51f97SToby Isaac 
1239be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest()
1249be51f97SToby Isaac @*/
1259be51f97SToby Isaac PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm)
126a0452a8eSToby Isaac {
127a0452a8eSToby Isaac   DM_Forest                  *forest = (DM_Forest*) dm->data;
12820e8089bSToby Isaac   DMType                     type;
129a0452a8eSToby Isaac   DM                         base;
130a0452a8eSToby Isaac   DMForestTopology           topology;
13105e99e11SStefano Zampini   MatType                    mtype;
132a0452a8eSToby Isaac   PetscInt                   dim, overlap, ref, factor;
133a0452a8eSToby Isaac   DMForestAdaptivityStrategy strat;
134795844e7SToby Isaac   PetscDS                    ds;
135795844e7SToby Isaac   void                       *ctx;
13649fc9a2fSToby Isaac   PetscErrorCode             (*map)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void*);
1373e58adeeSToby Isaac   void                       *mapCtx;
138a0452a8eSToby Isaac   PetscErrorCode             ierr;
139a0452a8eSToby Isaac 
140a0452a8eSToby Isaac   PetscFunctionBegin;
141a0452a8eSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14220e8089bSToby Isaac   ierr = DMCreate(PetscObjectComm((PetscObject)dm),tdm);CHKERRQ(ierr);
14320e8089bSToby Isaac   ierr = DMGetType(dm,&type);CHKERRQ(ierr);
14420e8089bSToby Isaac   ierr = DMSetType(*tdm,type);CHKERRQ(ierr);
145a0452a8eSToby Isaac   ierr = DMForestGetBaseDM(dm,&base);CHKERRQ(ierr);
14620e8089bSToby Isaac   ierr = DMForestSetBaseDM(*tdm,base);CHKERRQ(ierr);
147a0452a8eSToby Isaac   ierr = DMForestGetTopology(dm,&topology);CHKERRQ(ierr);
14820e8089bSToby Isaac   ierr = DMForestSetTopology(*tdm,topology);CHKERRQ(ierr);
149a0452a8eSToby Isaac   ierr = DMForestGetAdjacencyDimension(dm,&dim);CHKERRQ(ierr);
15020e8089bSToby Isaac   ierr = DMForestSetAdjacencyDimension(*tdm,dim);CHKERRQ(ierr);
151a0452a8eSToby Isaac   ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr);
15220e8089bSToby Isaac   ierr = DMForestSetPartitionOverlap(*tdm,overlap);CHKERRQ(ierr);
153a0452a8eSToby Isaac   ierr = DMForestGetMinimumRefinement(dm,&ref);CHKERRQ(ierr);
15420e8089bSToby Isaac   ierr = DMForestSetMinimumRefinement(*tdm,ref);CHKERRQ(ierr);
155a0452a8eSToby Isaac   ierr = DMForestGetMaximumRefinement(dm,&ref);CHKERRQ(ierr);
15620e8089bSToby Isaac   ierr = DMForestSetMaximumRefinement(*tdm,ref);CHKERRQ(ierr);
157a0452a8eSToby Isaac   ierr = DMForestGetAdaptivityStrategy(dm,&strat);CHKERRQ(ierr);
15820e8089bSToby Isaac   ierr = DMForestSetAdaptivityStrategy(*tdm,strat);CHKERRQ(ierr);
159a0452a8eSToby Isaac   ierr = DMForestGetGradeFactor(dm,&factor);CHKERRQ(ierr);
16020e8089bSToby Isaac   ierr = DMForestSetGradeFactor(*tdm,factor);CHKERRQ(ierr);
1613e58adeeSToby Isaac   ierr = DMForestGetBaseCoordinateMapping(dm,&map,&mapCtx);CHKERRQ(ierr);
1625e8c540aSToby Isaac   ierr = DMForestSetBaseCoordinateMapping(*tdm,map,mapCtx);CHKERRQ(ierr);
163a0452a8eSToby Isaac   if (forest->ftemplate) {
16420e8089bSToby Isaac     ierr = (forest->ftemplate)(dm, *tdm);CHKERRQ(ierr);
165a0452a8eSToby Isaac   }
16620e8089bSToby Isaac   ierr = DMForestSetAdaptivityForest(*tdm,dm);CHKERRQ(ierr);
167795844e7SToby Isaac   ierr = DMGetDS(dm,&ds);CHKERRQ(ierr);
168795844e7SToby Isaac   ierr = DMSetDS(*tdm,ds);CHKERRQ(ierr);
169795844e7SToby Isaac   ierr = DMGetApplicationContext(dm,&ctx);CHKERRQ(ierr);
170795844e7SToby Isaac   ierr = DMSetApplicationContext(*tdm,&ctx);CHKERRQ(ierr);
17190b157c4SStefano Zampini   {
17290b157c4SStefano Zampini     PetscBool            isper;
173795844e7SToby Isaac     const PetscReal      *maxCell, *L;
174795844e7SToby Isaac     const DMBoundaryType *bd;
175795844e7SToby Isaac 
17690b157c4SStefano Zampini     ierr = DMGetPeriodicity(dm,&isper,&maxCell,&L,&bd);CHKERRQ(ierr);
17790b157c4SStefano Zampini     ierr = DMSetPeriodicity(*tdm,isper,maxCell,L,bd);CHKERRQ(ierr);
178795844e7SToby Isaac   }
179bff67a9bSToby Isaac   ierr = DMCopyBoundary(dm,*tdm);CHKERRQ(ierr);
18005e99e11SStefano Zampini   ierr = DMGetMatType(dm,&mtype);CHKERRQ(ierr);
18105e99e11SStefano Zampini   ierr = DMSetMatType(*tdm,mtype);CHKERRQ(ierr);
182a0452a8eSToby Isaac   PetscFunctionReturn(0);
183a0452a8eSToby Isaac }
184a0452a8eSToby Isaac 
18501d9d024SToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm);
18601d9d024SToby Isaac 
187db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm)
188db4d5e8cSToby Isaac {
189db4d5e8cSToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
190db4d5e8cSToby Isaac   const char     *type;
191db4d5e8cSToby Isaac   PetscErrorCode ierr;
192db4d5e8cSToby Isaac 
193db4d5e8cSToby Isaac   PetscFunctionBegin;
194db4d5e8cSToby Isaac   forest->refct++;
195db4d5e8cSToby Isaac   (*newdm)->data = forest;
196db4d5e8cSToby Isaac   ierr           = PetscObjectGetType((PetscObject) dm, &type);CHKERRQ(ierr);
197db4d5e8cSToby Isaac   ierr           = PetscObjectChangeTypeName((PetscObject) *newdm, type);CHKERRQ(ierr);
19801d9d024SToby Isaac   ierr           = DMInitialize_Forest(*newdm);CHKERRQ(ierr);
199db4d5e8cSToby Isaac   PetscFunctionReturn(0);
200db4d5e8cSToby Isaac }
201db4d5e8cSToby Isaac 
202d222f98bSToby Isaac static PetscErrorCode DMDestroy_Forest(DM dm)
203db4d5e8cSToby Isaac {
204db4d5e8cSToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
205db4d5e8cSToby Isaac   PetscErrorCode ierr;
206db4d5e8cSToby Isaac 
207db4d5e8cSToby Isaac   PetscFunctionBegin;
208db4d5e8cSToby Isaac   if (--forest->refct > 0) PetscFunctionReturn(0);
209d222f98bSToby Isaac   if (forest->destroy) {ierr = forest->destroy(dm);CHKERRQ(ierr);}
210db4d5e8cSToby Isaac   ierr = PetscSFDestroy(&forest->cellSF);CHKERRQ(ierr);
2110f17b9e3SToby Isaac   ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr);
2120f17b9e3SToby Isaac   ierr = PetscSFDestroy(&forest->coarseToPreFine);CHKERRQ(ierr);
213a1b0c543SToby Isaac   ierr = DMLabelDestroy(&forest->adaptLabel);CHKERRQ(ierr);
2149a81d013SToby Isaac   ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr);
21556ba9f64SToby Isaac   ierr = DMDestroy(&forest->base);CHKERRQ(ierr);
216ba936b91SToby Isaac   ierr = DMDestroy(&forest->adapt);CHKERRQ(ierr);
21730f902e7SToby Isaac   ierr = PetscFree(forest->topology);CHKERRQ(ierr);
21830f902e7SToby Isaac   ierr = PetscFree(forest);CHKERRQ(ierr);
219db4d5e8cSToby Isaac   PetscFunctionReturn(0);
220db4d5e8cSToby Isaac }
221db4d5e8cSToby Isaac 
2229be51f97SToby Isaac /*@C
2239be51f97SToby Isaac   DMForestSetTopology - Set the topology of a DMForest during the pre-setup phase.  The topology is a string (e.g.
22468d54884SBarry Smith   "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest during
2259be51f97SToby Isaac   DMSetUp().
2269be51f97SToby Isaac 
2279be51f97SToby Isaac   Logically collective on dm
2289be51f97SToby Isaac 
2299be51f97SToby Isaac   Input parameters:
2309be51f97SToby Isaac + dm - the forest
2319be51f97SToby Isaac - topology - the topology of the forest
2329be51f97SToby Isaac 
2339be51f97SToby Isaac   Level: intermediate
2349be51f97SToby Isaac 
2359be51f97SToby Isaac .seealso(): DMForestGetTopology(), DMForestSetBaseDM()
2369be51f97SToby Isaac @*/
237dd8e54a2SToby Isaac PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology)
238db4d5e8cSToby Isaac {
239db4d5e8cSToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
240db4d5e8cSToby Isaac   PetscErrorCode ierr;
241db4d5e8cSToby Isaac 
242db4d5e8cSToby Isaac   PetscFunctionBegin;
243db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
244ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the topology after setup");
245dd8e54a2SToby Isaac   ierr = PetscFree(forest->topology);CHKERRQ(ierr);
246dd8e54a2SToby Isaac   ierr = PetscStrallocpy((const char*)topology,(char**) &forest->topology);CHKERRQ(ierr);
247db4d5e8cSToby Isaac   PetscFunctionReturn(0);
248db4d5e8cSToby Isaac }
249db4d5e8cSToby Isaac 
2509be51f97SToby Isaac /*@C
2519be51f97SToby Isaac   DMForestGetTopology - Get a string describing the topology of a DMForest.
2529be51f97SToby Isaac 
2539be51f97SToby Isaac   Not collective
2549be51f97SToby Isaac 
2559be51f97SToby Isaac   Input parameter:
2569be51f97SToby Isaac . dm - the forest
2579be51f97SToby Isaac 
2589be51f97SToby Isaac   Output parameter:
2599be51f97SToby Isaac . topology - the topology of the forest (e.g., 'cube', 'shell')
2609be51f97SToby Isaac 
2619be51f97SToby Isaac   Level: intermediate
2629be51f97SToby Isaac 
2639be51f97SToby Isaac .seealso: DMForestSetTopology()
2649be51f97SToby Isaac @*/
265dd8e54a2SToby Isaac PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology)
266dd8e54a2SToby Isaac {
267dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
268dd8e54a2SToby Isaac 
269dd8e54a2SToby Isaac   PetscFunctionBegin;
270dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
271dd8e54a2SToby Isaac   PetscValidPointer(topology,2);
272dd8e54a2SToby Isaac   *topology = forest->topology;
273dd8e54a2SToby Isaac   PetscFunctionReturn(0);
274dd8e54a2SToby Isaac }
275dd8e54a2SToby Isaac 
2769be51f97SToby Isaac /*@
2779be51f97SToby Isaac   DMForestSetBaseDM - During the pre-setup phase, set the DM that defines the base mesh of a DMForest forest.  The
2789be51f97SToby Isaac   forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its
279765b024eSBarry Smith   base.  In general, two forest must share a base to be comparable, to do things like construct interpolators.
2809be51f97SToby Isaac 
2819be51f97SToby Isaac   Logically collective on dm
2829be51f97SToby Isaac 
2839be51f97SToby Isaac   Input Parameters:
2849be51f97SToby Isaac + dm - the forest
2859be51f97SToby Isaac - base - the base DM of the forest
2869be51f97SToby Isaac 
28795452b02SPatrick Sanan   Notes:
28895452b02SPatrick Sanan     Currently the base DM must be a DMPLEX
289765b024eSBarry Smith 
2909be51f97SToby Isaac   Level: intermediate
2919be51f97SToby Isaac 
2929be51f97SToby Isaac .seealso(): DMForestGetBaseDM()
2939be51f97SToby Isaac @*/
294dd8e54a2SToby Isaac PetscErrorCode DMForestSetBaseDM(DM dm, DM base)
295dd8e54a2SToby Isaac {
296dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
297dd8e54a2SToby Isaac   PetscInt       dim, dimEmbed;
298dd8e54a2SToby Isaac   PetscErrorCode ierr;
299dd8e54a2SToby Isaac 
300dd8e54a2SToby Isaac   PetscFunctionBegin;
301dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
302ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the base after setup");
303dd8e54a2SToby Isaac   ierr         = PetscObjectReference((PetscObject)base);CHKERRQ(ierr);
304dd8e54a2SToby Isaac   ierr         = DMDestroy(&forest->base);CHKERRQ(ierr);
305dd8e54a2SToby Isaac   forest->base = base;
306a0452a8eSToby Isaac   if (base) {
307*28dfcf7cSStefano Zampini     PetscBool        isper;
308*28dfcf7cSStefano Zampini     const PetscReal *maxCell, *L;
309*28dfcf7cSStefano Zampini     const DMBoundaryType *bd;
310*28dfcf7cSStefano Zampini 
311a0452a8eSToby Isaac     PetscValidHeaderSpecific(base, DM_CLASSID, 2);
312dd8e54a2SToby Isaac     ierr = DMGetDimension(base,&dim);CHKERRQ(ierr);
313dd8e54a2SToby Isaac     ierr = DMSetDimension(dm,dim);CHKERRQ(ierr);
314dd8e54a2SToby Isaac     ierr = DMGetCoordinateDim(base,&dimEmbed);CHKERRQ(ierr);
315dd8e54a2SToby Isaac     ierr = DMSetCoordinateDim(dm,dimEmbed);CHKERRQ(ierr);
316*28dfcf7cSStefano Zampini     ierr = DMGetPeriodicity(base,&isper,&maxCell,&L,&bd);CHKERRQ(ierr);
317*28dfcf7cSStefano Zampini     ierr = DMSetPeriodicity(dm,isper,maxCell,L,bd);CHKERRQ(ierr);
318*28dfcf7cSStefano Zampini   } else {
319*28dfcf7cSStefano Zampini     ierr = DMSetPeriodicity(dm,PETSC_FALSE,NULL,NULL,NULL);CHKERRQ(ierr);
320a0452a8eSToby Isaac   }
321dd8e54a2SToby Isaac   PetscFunctionReturn(0);
322dd8e54a2SToby Isaac }
323dd8e54a2SToby Isaac 
3249be51f97SToby Isaac /*@
3259be51f97SToby Isaac   DMForestGetBaseDM - Get the base DM of a DMForest forest.  The forest will be hierarchically refined from the base,
32668d54884SBarry Smith   and all refinements/coarsenings of the forest will share its base.  In general, two forest must share a base to be
3279be51f97SToby Isaac   comparable, to do things like construct interpolators.
3289be51f97SToby Isaac 
3299be51f97SToby Isaac   Not collective
3309be51f97SToby Isaac 
3319be51f97SToby Isaac   Input Parameter:
3329be51f97SToby Isaac . dm - the forest
3339be51f97SToby Isaac 
3349be51f97SToby Isaac   Output Parameter:
3359be51f97SToby Isaac . base - the base DM of the forest
3369be51f97SToby Isaac 
337367003a6SStefano Zampini   Notes:
338367003a6SStefano Zampini     After DMSetUp(), the base DM will be redundantly distributed across MPI processes
339367003a6SStefano Zampini 
3409be51f97SToby Isaac   Level: intermediate
3419be51f97SToby Isaac 
3429be51f97SToby Isaac .seealso(); DMForestSetBaseDM()
3439be51f97SToby Isaac @*/
344dd8e54a2SToby Isaac PetscErrorCode DMForestGetBaseDM(DM dm, DM *base)
345dd8e54a2SToby Isaac {
346dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
347dd8e54a2SToby Isaac 
348dd8e54a2SToby Isaac   PetscFunctionBegin;
349dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
350dd8e54a2SToby Isaac   PetscValidPointer(base, 2);
351dd8e54a2SToby Isaac   *base = forest->base;
352dd8e54a2SToby Isaac   PetscFunctionReturn(0);
353dd8e54a2SToby Isaac }
354dd8e54a2SToby Isaac 
35599478f86SToby Isaac PetscErrorCode DMForestSetBaseCoordinateMapping(DM dm, PetscErrorCode (*func)(DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx)
356cf38a08cSToby Isaac {
357cf38a08cSToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
358cf38a08cSToby Isaac 
359cf38a08cSToby Isaac   PetscFunctionBegin;
360cf38a08cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
361cf38a08cSToby Isaac   forest->mapcoordinates    = func;
362cf38a08cSToby Isaac   forest->mapcoordinatesctx = ctx;
363cf38a08cSToby Isaac   PetscFunctionReturn(0);
364cf38a08cSToby Isaac }
365cf38a08cSToby Isaac 
36699478f86SToby Isaac PetscErrorCode DMForestGetBaseCoordinateMapping(DM dm, PetscErrorCode (**func) (DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx)
367cf38a08cSToby Isaac {
368cf38a08cSToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
369cf38a08cSToby Isaac 
370cf38a08cSToby Isaac   PetscFunctionBegin;
371cf38a08cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
372cf38a08cSToby Isaac   if (func) *func = forest->mapcoordinates;
373cf38a08cSToby Isaac   if (ctx) *((void**) ctx) = forest->mapcoordinatesctx;
374cf38a08cSToby Isaac   PetscFunctionReturn(0);
375cf38a08cSToby Isaac }
376cf38a08cSToby Isaac 
3779be51f97SToby Isaac /*@
3789be51f97SToby Isaac   DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be
3799be51f97SToby Isaac   adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) im DMSetUp().  Usually not needed
3809be51f97SToby Isaac   by users directly: DMForestTemplate() constructs a new forest to be adapted from an old forest and calls this
3819be51f97SToby Isaac   routine.
3829be51f97SToby Isaac 
383dffe73a3SToby Isaac   Note that this can be called after setup with adapt = NULL, which will clear all internal data related to the
384dffe73a3SToby Isaac   adaptivity forest from dm.  This way, repeatedly adapting does not leave stale DM objects in memory.
385dffe73a3SToby Isaac 
3869be51f97SToby Isaac   Logically collective on dm
3879be51f97SToby Isaac 
3889be51f97SToby Isaac   Input Parameter:
3899be51f97SToby Isaac + dm - the new forest, which will be constructed from adapt
3909be51f97SToby Isaac - adapt - the old forest
3919be51f97SToby Isaac 
3929be51f97SToby Isaac   Level: intermediate
3939be51f97SToby Isaac 
3949be51f97SToby Isaac .seealso: DMForestGetAdaptivityForest(), DMForestSetAdaptivityPurpose()
3959be51f97SToby Isaac @*/
396ba936b91SToby Isaac PetscErrorCode DMForestSetAdaptivityForest(DM dm,DM adapt)
397dd8e54a2SToby Isaac {
398dffe73a3SToby Isaac   DM_Forest      *forest, *adaptForest, *oldAdaptForest;
399dffe73a3SToby Isaac   DM             oldAdapt;
400456cc5b7SMatthew G. Knepley   PetscBool      isForest;
401dd8e54a2SToby Isaac   PetscErrorCode ierr;
402dd8e54a2SToby Isaac 
403dd8e54a2SToby Isaac   PetscFunctionBegin;
404dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
405ba936b91SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
406456cc5b7SMatthew G. Knepley   ierr = DMIsForest(dm, &isForest);CHKERRQ(ierr);
407456cc5b7SMatthew G. Knepley   if (!isForest) PetscFunctionReturn(0);
408ba936b91SToby Isaac   forest   = (DM_Forest*) dm->data;
409dffe73a3SToby Isaac   ierr     = DMForestGetAdaptivityForest(dm,&oldAdapt);CHKERRQ(ierr);
410dffe73a3SToby Isaac   if (adapt != NULL && dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
411193eb951SToby Isaac   adaptForest    = (DM_Forest*) (adapt ? adapt->data : NULL);
412193eb951SToby Isaac   oldAdaptForest = (DM_Forest*) (oldAdapt ? oldAdapt->data : NULL);
413dffe73a3SToby Isaac   if (adaptForest != oldAdaptForest) {
414dffe73a3SToby Isaac     ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr);
415dffe73a3SToby Isaac     ierr = PetscSFDestroy(&forest->coarseToPreFine);CHKERRQ(ierr);
416dffe73a3SToby Isaac     if (forest->clearadaptivityforest) {ierr = (forest->clearadaptivityforest)(dm);CHKERRQ(ierr);}
417dffe73a3SToby Isaac   }
41826d9498aSToby Isaac   switch (forest->adaptPurpose) {
419cd3c525cSToby Isaac   case DM_ADAPT_DETERMINE:
420ba936b91SToby Isaac     ierr          = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr);
421ba936b91SToby Isaac     ierr          = DMDestroy(&(forest->adapt));CHKERRQ(ierr);
422ba936b91SToby Isaac     forest->adapt = adapt;
42326d9498aSToby Isaac     break;
424a1b0c543SToby Isaac   case DM_ADAPT_REFINE:
42526d9498aSToby Isaac     ierr = DMSetCoarseDM(dm,adapt);CHKERRQ(ierr);
42626d9498aSToby Isaac     break;
427a1b0c543SToby Isaac   case DM_ADAPT_COARSEN:
42826d9498aSToby Isaac     ierr = DMSetFineDM(dm,adapt);CHKERRQ(ierr);
42926d9498aSToby Isaac     break;
43026d9498aSToby Isaac   default:
43126d9498aSToby Isaac     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
43226d9498aSToby Isaac   }
433dd8e54a2SToby Isaac   PetscFunctionReturn(0);
434dd8e54a2SToby Isaac }
435dd8e54a2SToby Isaac 
4369be51f97SToby Isaac /*@
4379be51f97SToby Isaac   DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted.
4389be51f97SToby Isaac 
4399be51f97SToby Isaac   Not collective
4409be51f97SToby Isaac 
4419be51f97SToby Isaac   Input Parameter:
4429be51f97SToby Isaac . dm - the forest
4439be51f97SToby Isaac 
4449be51f97SToby Isaac   Output Parameter:
4459be51f97SToby Isaac . adapt - the forest from which dm is/was adapted
4469be51f97SToby Isaac 
4479be51f97SToby Isaac   Level: intermediate
4489be51f97SToby Isaac 
4499be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose()
4509be51f97SToby Isaac @*/
451ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt)
452dd8e54a2SToby Isaac {
453ba936b91SToby Isaac   DM_Forest      *forest;
45426d9498aSToby Isaac   PetscErrorCode ierr;
455dd8e54a2SToby Isaac 
456dd8e54a2SToby Isaac   PetscFunctionBegin;
457dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
458ba936b91SToby Isaac   forest = (DM_Forest*) dm->data;
45926d9498aSToby Isaac   switch (forest->adaptPurpose) {
460cd3c525cSToby Isaac   case DM_ADAPT_DETERMINE:
461ba936b91SToby Isaac     *adapt = forest->adapt;
46226d9498aSToby Isaac     break;
463a1b0c543SToby Isaac   case DM_ADAPT_REFINE:
46426d9498aSToby Isaac     ierr = DMGetCoarseDM(dm,adapt);CHKERRQ(ierr);
46526d9498aSToby Isaac     break;
466a1b0c543SToby Isaac   case DM_ADAPT_COARSEN:
46726d9498aSToby Isaac     ierr = DMGetFineDM(dm,adapt);CHKERRQ(ierr);
46826d9498aSToby Isaac     break;
46926d9498aSToby Isaac   default:
47026d9498aSToby Isaac     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
47126d9498aSToby Isaac   }
47226d9498aSToby Isaac   PetscFunctionReturn(0);
47326d9498aSToby Isaac }
47426d9498aSToby Isaac 
4759be51f97SToby Isaac /*@
4769be51f97SToby Isaac   DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its
477a1b0c543SToby Isaac   source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening
478cd3c525cSToby Isaac   (DM_ADAPT_COARSEN), or undefined (DM_ADAPT_DETERMINE).  This only matters for the purposes of reference counting:
4799be51f97SToby Isaac   during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse
4809be51f97SToby Isaac   relationship (see DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does
4819be51f97SToby Isaac   not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can
4829be51f97SToby Isaac   cause a memory leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
4839be51f97SToby Isaac 
4849be51f97SToby Isaac   Logically collective on dm
4859be51f97SToby Isaac 
4869be51f97SToby Isaac   Input Parameters:
4879be51f97SToby Isaac + dm - the forest
488cd3c525cSToby Isaac - purpose - the adaptivity purpose (DM_ADAPT_DETERMINE/DM_ADAPT_REFINE/DM_ADAPT_COARSEN)
4899be51f97SToby Isaac 
4909be51f97SToby Isaac   Level: advanced
4919be51f97SToby Isaac 
4929be51f97SToby Isaac .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest()
4939be51f97SToby Isaac @*/
494a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMAdaptFlag purpose)
49526d9498aSToby Isaac {
49626d9498aSToby Isaac   DM_Forest      *forest;
49726d9498aSToby Isaac   PetscErrorCode ierr;
49826d9498aSToby Isaac 
49926d9498aSToby Isaac   PetscFunctionBegin;
50026d9498aSToby Isaac   forest = (DM_Forest*) dm->data;
5019be51f97SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
50226d9498aSToby Isaac   if (purpose != forest->adaptPurpose) {
50326d9498aSToby Isaac     DM adapt;
50426d9498aSToby Isaac 
50526d9498aSToby Isaac     ierr = DMForestGetAdaptivityForest(dm,&adapt);CHKERRQ(ierr);
50626d9498aSToby Isaac     ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr);
50726d9498aSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
508f885a11aSToby Isaac 
50926d9498aSToby Isaac     forest->adaptPurpose = purpose;
510f885a11aSToby Isaac 
51126d9498aSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,adapt);CHKERRQ(ierr);
51226d9498aSToby Isaac     ierr = DMDestroy(&adapt);CHKERRQ(ierr);
51326d9498aSToby Isaac   }
514dd8e54a2SToby Isaac   PetscFunctionReturn(0);
515dd8e54a2SToby Isaac }
516dd8e54a2SToby Isaac 
51756c0450aSToby Isaac /*@
51856c0450aSToby Isaac   DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with
519a1b0c543SToby Isaac   DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening (DM_ADAPT_COARSEN), or
520cd3c525cSToby Isaac   undefined (DM_ADAPT_DETERMINE).  This only matters for the purposes of reference counting: during DMDestroy(), cyclic
52156c0450aSToby Isaac   references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see
52256c0450aSToby Isaac   DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does not maintain a
52356c0450aSToby Isaac   reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory
52456c0450aSToby Isaac   leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
52556c0450aSToby Isaac 
52656c0450aSToby Isaac   Not collective
52756c0450aSToby Isaac 
52856c0450aSToby Isaac   Input Parameter:
52956c0450aSToby Isaac . dm - the forest
53056c0450aSToby Isaac 
53156c0450aSToby Isaac   Output Parameter:
532cd3c525cSToby Isaac . purpose - the adaptivity purpose (DM_ADAPT_DETERMINE/DM_ADAPT_REFINE/DM_ADAPT_COARSEN)
53356c0450aSToby Isaac 
53456c0450aSToby Isaac   Level: advanced
53556c0450aSToby Isaac 
53656c0450aSToby Isaac .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest()
53756c0450aSToby Isaac @*/
538a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMAdaptFlag *purpose)
53956c0450aSToby Isaac {
54056c0450aSToby Isaac   DM_Forest *forest;
54156c0450aSToby Isaac 
54256c0450aSToby Isaac   PetscFunctionBegin;
54356c0450aSToby Isaac   forest   = (DM_Forest*) dm->data;
54456c0450aSToby Isaac   *purpose = forest->adaptPurpose;
54556c0450aSToby Isaac   PetscFunctionReturn(0);
54656c0450aSToby Isaac }
54756c0450aSToby Isaac 
5489be51f97SToby Isaac /*@
5499be51f97SToby Isaac   DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine
5509be51f97SToby Isaac   cell adjacency (for the purposes of partitioning and overlap).
5519be51f97SToby Isaac 
5529be51f97SToby Isaac   Logically collective on dm
5539be51f97SToby Isaac 
5549be51f97SToby Isaac   Input Parameters:
5559be51f97SToby Isaac + dm - the forest
5569be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency)
5579be51f97SToby Isaac 
5589be51f97SToby Isaac   Level: intermediate
5599be51f97SToby Isaac 
5609be51f97SToby Isaac .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap()
5619be51f97SToby Isaac @*/
562dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim)
563dd8e54a2SToby Isaac {
564dd8e54a2SToby Isaac   PetscInt       dim;
565dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
566dd8e54a2SToby Isaac   PetscErrorCode ierr;
567dd8e54a2SToby Isaac 
568dd8e54a2SToby Isaac   PetscFunctionBegin;
569dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
570ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup");
571dd8e54a2SToby Isaac   if (adjDim < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %d", adjDim);
572dd8e54a2SToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
573dd8e54a2SToby Isaac   if (adjDim > dim) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %d: %d", dim, adjDim);
574dd8e54a2SToby Isaac   forest->adjDim = adjDim;
575dd8e54a2SToby Isaac   PetscFunctionReturn(0);
576dd8e54a2SToby Isaac }
577dd8e54a2SToby Isaac 
5789be51f97SToby Isaac /*@
5799be51f97SToby Isaac   DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that,
5809be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
5819be51f97SToby Isaac 
5829be51f97SToby Isaac   Logically collective on dm
5839be51f97SToby Isaac 
5849be51f97SToby Isaac   Input Parameters:
5859be51f97SToby Isaac + dm - the forest
5869be51f97SToby Isaac - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
5879be51f97SToby Isaac 
5889be51f97SToby Isaac   Level: intermediate
5899be51f97SToby Isaac 
5909be51f97SToby Isaac .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension()
5919be51f97SToby Isaac @*/
592dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim)
593dd8e54a2SToby Isaac {
594dd8e54a2SToby Isaac   PetscInt       dim;
595dd8e54a2SToby Isaac   PetscErrorCode ierr;
596dd8e54a2SToby Isaac 
597dd8e54a2SToby Isaac   PetscFunctionBegin;
598dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
599dd8e54a2SToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
600dd8e54a2SToby Isaac   ierr = DMForestSetAdjacencyDimension(dm,dim-adjCodim);CHKERRQ(ierr);
601dd8e54a2SToby Isaac   PetscFunctionReturn(0);
602dd8e54a2SToby Isaac }
603dd8e54a2SToby Isaac 
6049be51f97SToby Isaac /*@
6059be51f97SToby Isaac   DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the
6069be51f97SToby Isaac   purposes of partitioning and overlap).
6079be51f97SToby Isaac 
6089be51f97SToby Isaac   Not collective
6099be51f97SToby Isaac 
6109be51f97SToby Isaac   Input Parameter:
6119be51f97SToby Isaac . dm - the forest
6129be51f97SToby Isaac 
6139be51f97SToby Isaac   Output Parameter:
6149be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency)
6159be51f97SToby Isaac 
6169be51f97SToby Isaac   Level: intermediate
6179be51f97SToby Isaac 
6189be51f97SToby Isaac .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap()
6199be51f97SToby Isaac @*/
620dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim)
621dd8e54a2SToby Isaac {
622dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
623dd8e54a2SToby Isaac 
624dd8e54a2SToby Isaac   PetscFunctionBegin;
625dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
626dd8e54a2SToby Isaac   PetscValidIntPointer(adjDim,2);
627dd8e54a2SToby Isaac   *adjDim = forest->adjDim;
628dd8e54a2SToby Isaac   PetscFunctionReturn(0);
629dd8e54a2SToby Isaac }
630dd8e54a2SToby Isaac 
6319be51f97SToby Isaac /*@
6329be51f97SToby Isaac   DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that,
6339be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
6349be51f97SToby Isaac 
6359be51f97SToby Isaac   Not collective
6369be51f97SToby Isaac 
6379be51f97SToby Isaac   Input Parameter:
6389be51f97SToby Isaac . dm - the forest
6399be51f97SToby Isaac 
6409be51f97SToby Isaac   Output Parameter:
6419be51f97SToby Isaac . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
6429be51f97SToby Isaac 
6439be51f97SToby Isaac   Level: intermediate
6449be51f97SToby Isaac 
6459be51f97SToby Isaac .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension()
6469be51f97SToby Isaac @*/
647dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim)
648dd8e54a2SToby Isaac {
649dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
650dd8e54a2SToby Isaac   PetscInt       dim;
651dd8e54a2SToby Isaac   PetscErrorCode ierr;
652dd8e54a2SToby Isaac 
653dd8e54a2SToby Isaac   PetscFunctionBegin;
654dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
655dd8e54a2SToby Isaac   PetscValidIntPointer(adjCodim,2);
656dd8e54a2SToby Isaac   ierr      = DMGetDimension(dm,&dim);CHKERRQ(ierr);
657dd8e54a2SToby Isaac   *adjCodim = dim - forest->adjDim;
658dd8e54a2SToby Isaac   PetscFunctionReturn(0);
659dd8e54a2SToby Isaac }
660dd8e54a2SToby Isaac 
6619be51f97SToby Isaac /*@
6629be51f97SToby Isaac   DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel
6639be51f97SToby Isaac   partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding
6649be51f97SToby Isaac   adjacent cells
6659be51f97SToby Isaac 
6669be51f97SToby Isaac   Logically collective on dm
6679be51f97SToby Isaac 
6689be51f97SToby Isaac   Input Parameters:
6699be51f97SToby Isaac + dm - the forest
6709be51f97SToby Isaac - overlap - default 0
6719be51f97SToby Isaac 
6729be51f97SToby Isaac   Level: intermediate
6739be51f97SToby Isaac 
6749be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
6759be51f97SToby Isaac @*/
676dd8e54a2SToby Isaac PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap)
677dd8e54a2SToby Isaac {
678dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
679dd8e54a2SToby Isaac 
680dd8e54a2SToby Isaac   PetscFunctionBegin;
681dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
682ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup");
683dd8e54a2SToby Isaac   if (overlap < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %d", overlap);
684dd8e54a2SToby Isaac   forest->overlap = overlap;
685dd8e54a2SToby Isaac   PetscFunctionReturn(0);
686dd8e54a2SToby Isaac }
687dd8e54a2SToby Isaac 
6889be51f97SToby Isaac /*@
6899be51f97SToby Isaac   DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values
6909be51f97SToby Isaac   > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells
6919be51f97SToby Isaac 
6929be51f97SToby Isaac   Not collective
6939be51f97SToby Isaac 
6949be51f97SToby Isaac   Input Parameter:
6959be51f97SToby Isaac . dm - the forest
6969be51f97SToby Isaac 
6979be51f97SToby Isaac   Output Parameter:
6989be51f97SToby Isaac . overlap - default 0
6999be51f97SToby Isaac 
7009be51f97SToby Isaac   Level: intermediate
7019be51f97SToby Isaac 
7029be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
7039be51f97SToby Isaac @*/
704dd8e54a2SToby Isaac PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap)
705dd8e54a2SToby Isaac {
706dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
707dd8e54a2SToby Isaac 
708dd8e54a2SToby Isaac   PetscFunctionBegin;
709dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
710dd8e54a2SToby Isaac   PetscValidIntPointer(overlap,2);
711dd8e54a2SToby Isaac   *overlap = forest->overlap;
712dd8e54a2SToby Isaac   PetscFunctionReturn(0);
713dd8e54a2SToby Isaac }
714dd8e54a2SToby Isaac 
7159be51f97SToby Isaac /*@
7169be51f97SToby Isaac   DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base
7179be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest
7189be51f97SToby Isaac   (see DMForestGetAdaptivityForest()) this limits the amount of coarsening.
7199be51f97SToby Isaac 
7209be51f97SToby Isaac   Logically collective on dm
7219be51f97SToby Isaac 
7229be51f97SToby Isaac   Input Parameters:
7239be51f97SToby Isaac + dm - the forest
7249be51f97SToby Isaac - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7259be51f97SToby Isaac 
7269be51f97SToby Isaac   Level: intermediate
7279be51f97SToby Isaac 
7289be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
7299be51f97SToby Isaac @*/
730dd8e54a2SToby Isaac PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement)
731dd8e54a2SToby Isaac {
732dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
733dd8e54a2SToby Isaac 
734dd8e54a2SToby Isaac   PetscFunctionBegin;
735dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
736ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup");
737dd8e54a2SToby Isaac   forest->minRefinement = minRefinement;
738dd8e54a2SToby Isaac   PetscFunctionReturn(0);
739dd8e54a2SToby Isaac }
740dd8e54a2SToby Isaac 
7419be51f97SToby Isaac /*@
7429be51f97SToby Isaac   DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see
7439be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest (see
7449be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of coarsening.
7459be51f97SToby Isaac 
7469be51f97SToby Isaac   Not collective
7479be51f97SToby Isaac 
7489be51f97SToby Isaac   Input Parameter:
7499be51f97SToby Isaac . dm - the forest
7509be51f97SToby Isaac 
7519be51f97SToby Isaac   Output Parameter:
7529be51f97SToby Isaac . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7539be51f97SToby Isaac 
7549be51f97SToby Isaac   Level: intermediate
7559be51f97SToby Isaac 
7569be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestGetMaximumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
7579be51f97SToby Isaac @*/
758dd8e54a2SToby Isaac PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement)
759dd8e54a2SToby Isaac {
760dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
761dd8e54a2SToby Isaac 
762dd8e54a2SToby Isaac   PetscFunctionBegin;
763dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
764dd8e54a2SToby Isaac   PetscValidIntPointer(minRefinement,2);
765dd8e54a2SToby Isaac   *minRefinement = forest->minRefinement;
766dd8e54a2SToby Isaac   PetscFunctionReturn(0);
767dd8e54a2SToby Isaac }
768dd8e54a2SToby Isaac 
7699be51f97SToby Isaac /*@
7709be51f97SToby Isaac   DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base
7719be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.
7729be51f97SToby Isaac 
7739be51f97SToby Isaac   Logically collective on dm
7749be51f97SToby Isaac 
7759be51f97SToby Isaac   Input Parameters:
7769be51f97SToby Isaac + dm - the forest
7779be51f97SToby Isaac - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7789be51f97SToby Isaac 
7799be51f97SToby Isaac   Level: intermediate
7809be51f97SToby Isaac 
7819be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
7829be51f97SToby Isaac @*/
78356ba9f64SToby Isaac PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement)
78456ba9f64SToby Isaac {
78556ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
78656ba9f64SToby Isaac 
78756ba9f64SToby Isaac   PetscFunctionBegin;
78856ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
78956ba9f64SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup");
79056ba9f64SToby Isaac   forest->initRefinement = initRefinement;
79156ba9f64SToby Isaac   PetscFunctionReturn(0);
79256ba9f64SToby Isaac }
79356ba9f64SToby Isaac 
7949be51f97SToby Isaac /*@
7959be51f97SToby Isaac   DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see
7969be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.
7979be51f97SToby Isaac 
7989be51f97SToby Isaac   Not collective
7999be51f97SToby Isaac 
8009be51f97SToby Isaac   Input Parameter:
8019be51f97SToby Isaac . dm - the forest
8029be51f97SToby Isaac 
8039be51f97SToby Isaac   Output Paramater:
8049be51f97SToby Isaac . initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8059be51f97SToby Isaac 
8069be51f97SToby Isaac   Level: intermediate
8079be51f97SToby Isaac 
8089be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
8099be51f97SToby Isaac @*/
81056ba9f64SToby Isaac PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement)
81156ba9f64SToby Isaac {
81256ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
81356ba9f64SToby Isaac 
81456ba9f64SToby Isaac   PetscFunctionBegin;
81556ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
81656ba9f64SToby Isaac   PetscValidIntPointer(initRefinement,2);
81756ba9f64SToby Isaac   *initRefinement = forest->initRefinement;
81856ba9f64SToby Isaac   PetscFunctionReturn(0);
81956ba9f64SToby Isaac }
82056ba9f64SToby Isaac 
8219be51f97SToby Isaac /*@
8229be51f97SToby Isaac   DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base
8239be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest
8249be51f97SToby Isaac   (see DMForestGetAdaptivityForest()), this limits the amount of refinement.
8259be51f97SToby Isaac 
8269be51f97SToby Isaac   Logically collective on dm
8279be51f97SToby Isaac 
8289be51f97SToby Isaac   Input Parameters:
8299be51f97SToby Isaac + dm - the forest
8309be51f97SToby Isaac - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8319be51f97SToby Isaac 
8329be51f97SToby Isaac   Level: intermediate
8339be51f97SToby Isaac 
8349be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityDM()
8359be51f97SToby Isaac @*/
836c7eeac06SToby Isaac PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement)
837dd8e54a2SToby Isaac {
838dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
839dd8e54a2SToby Isaac 
840dd8e54a2SToby Isaac   PetscFunctionBegin;
841dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
842ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup");
843c7eeac06SToby Isaac   forest->maxRefinement = maxRefinement;
844dd8e54a2SToby Isaac   PetscFunctionReturn(0);
845dd8e54a2SToby Isaac }
846dd8e54a2SToby Isaac 
8479be51f97SToby Isaac /*@
8489be51f97SToby Isaac   DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see
8499be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest (see
8509be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of refinement.
8519be51f97SToby Isaac 
8529be51f97SToby Isaac   Not collective
8539be51f97SToby Isaac 
8549be51f97SToby Isaac   Input Parameter:
8559be51f97SToby Isaac . dm - the forest
8569be51f97SToby Isaac 
8579be51f97SToby Isaac   Output Parameter:
8589be51f97SToby Isaac . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8599be51f97SToby Isaac 
8609be51f97SToby Isaac   Level: intermediate
8619be51f97SToby Isaac 
8629be51f97SToby Isaac .seealso: DMForestSetMaximumRefinement(), DMForestGetMinimumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
8639be51f97SToby Isaac @*/
864c7eeac06SToby Isaac PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement)
865dd8e54a2SToby Isaac {
866dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
867dd8e54a2SToby Isaac 
868dd8e54a2SToby Isaac   PetscFunctionBegin;
869dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
870c7eeac06SToby Isaac   PetscValidIntPointer(maxRefinement,2);
871c7eeac06SToby Isaac   *maxRefinement = forest->maxRefinement;
872dd8e54a2SToby Isaac   PetscFunctionReturn(0);
873dd8e54a2SToby Isaac }
874c7eeac06SToby Isaac 
8759be51f97SToby Isaac /*@C
8769be51f97SToby Isaac   DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes.
8779be51f97SToby Isaac   Subtypes of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree
8789be51f97SToby Isaac   for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to
8799be51f97SToby Isaac   specify refinement/coarsening.
8809be51f97SToby Isaac 
8819be51f97SToby Isaac   Logically collective on dm
8829be51f97SToby Isaac 
8839be51f97SToby Isaac   Input Parameters:
8849be51f97SToby Isaac + dm - the forest
8859be51f97SToby Isaac - adaptStrategy - default DMFORESTADAPTALL
8869be51f97SToby Isaac 
8879be51f97SToby Isaac   Level: advanced
8889be51f97SToby Isaac 
8899be51f97SToby Isaac .seealso: DMForestGetAdaptivityStrategy()
8909be51f97SToby Isaac @*/
891c7eeac06SToby Isaac PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy)
892c7eeac06SToby Isaac {
893c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
894c7eeac06SToby Isaac   PetscErrorCode ierr;
895c7eeac06SToby Isaac 
896c7eeac06SToby Isaac   PetscFunctionBegin;
897c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
898c7eeac06SToby Isaac   ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr);
899a73e2921SToby Isaac   ierr = PetscStrallocpy((const char*) adaptStrategy,(char**)&forest->adaptStrategy);CHKERRQ(ierr);
900c7eeac06SToby Isaac   PetscFunctionReturn(0);
901c7eeac06SToby Isaac }
902c7eeac06SToby Isaac 
9039be51f97SToby Isaac /*@C
9049be51f97SToby Isaac   DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes.  Subtypes
9059be51f97SToby Isaac   of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all
9069be51f97SToby Isaac   processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only
9079be51f97SToby Isaac   one process needs to specify refinement/coarsening.
9089be51f97SToby Isaac 
9099be51f97SToby Isaac   Not collective
9109be51f97SToby Isaac 
9119be51f97SToby Isaac   Input Parameter:
9129be51f97SToby Isaac . dm - the forest
9139be51f97SToby Isaac 
9149be51f97SToby Isaac   Output Parameter:
9159be51f97SToby Isaac . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL)
9169be51f97SToby Isaac 
9179be51f97SToby Isaac   Level: advanced
9189be51f97SToby Isaac 
9199be51f97SToby Isaac .seealso: DMForestSetAdaptivityStrategy()
9209be51f97SToby Isaac @*/
921c7eeac06SToby Isaac PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy)
922c7eeac06SToby Isaac {
923c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
924c7eeac06SToby Isaac 
925c7eeac06SToby Isaac   PetscFunctionBegin;
926c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
927c7eeac06SToby Isaac   PetscValidPointer(adaptStrategy,2);
928c7eeac06SToby Isaac   *adaptStrategy = forest->adaptStrategy;
929c7eeac06SToby Isaac   PetscFunctionReturn(0);
930c7eeac06SToby Isaac }
931c7eeac06SToby Isaac 
9322a133e43SToby Isaac /*@
9332a133e43SToby Isaac   DMForestGetAdaptivitySuccess - Return whether the requested adaptation (refinement, coarsening, repartitioning,
9342a133e43SToby Isaac   etc.) was successful.  PETSC_FALSE indicates that the post-adaptation forest is the same as the pre-adpatation
9352a133e43SToby Isaac   forest.  A requested adaptation may have been unsuccessful if, for example, the requested refinement would have
9362a133e43SToby Isaac   exceeded the maximum refinement level.
9372a133e43SToby Isaac 
9382a133e43SToby Isaac   Collective on dm
9392a133e43SToby Isaac 
9402a133e43SToby Isaac   Input Parameter:
9412a133e43SToby Isaac 
9422a133e43SToby Isaac . dm - the post-adaptation forest
9432a133e43SToby Isaac 
9442a133e43SToby Isaac   Output Parameter:
9452a133e43SToby Isaac 
9462a133e43SToby Isaac . success - PETSC_TRUE if the post-adaptation forest is different from the pre-adaptation forest.
9472a133e43SToby Isaac 
9482a133e43SToby Isaac   Level: intermediate
9492a133e43SToby Isaac 
9502a133e43SToby Isaac .see
9512a133e43SToby Isaac @*/
9522a133e43SToby Isaac PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success)
9532a133e43SToby Isaac {
9542a133e43SToby Isaac   DM_Forest      *forest;
9552a133e43SToby Isaac   PetscErrorCode ierr;
9562a133e43SToby Isaac 
9572a133e43SToby Isaac   PetscFunctionBegin;
9582a133e43SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9592a133e43SToby Isaac   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() has not been called yet.");
9602a133e43SToby Isaac   forest = (DM_Forest *) dm->data;
9612a133e43SToby Isaac   ierr = (forest->getadaptivitysuccess)(dm,success);CHKERRQ(ierr);
9622a133e43SToby Isaac   PetscFunctionReturn(0);
9632a133e43SToby Isaac }
9642a133e43SToby Isaac 
965bf9b5d84SToby Isaac /*@
966bf9b5d84SToby Isaac   DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed
967bf9b5d84SToby 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().
968bf9b5d84SToby Isaac 
969bf9b5d84SToby Isaac   Logically collective on dm
970bf9b5d84SToby Isaac 
971bf9b5d84SToby Isaac   Input Parameters:
972bf9b5d84SToby Isaac + dm - the post-adaptation forest
973bf9b5d84SToby Isaac - computeSF - default PETSC_TRUE
974bf9b5d84SToby Isaac 
975bf9b5d84SToby Isaac   Level: advanced
976bf9b5d84SToby Isaac 
977bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
978bf9b5d84SToby Isaac @*/
979bf9b5d84SToby Isaac PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF)
980bf9b5d84SToby Isaac {
981bf9b5d84SToby Isaac   DM_Forest *forest;
982bf9b5d84SToby Isaac 
983bf9b5d84SToby Isaac   PetscFunctionBegin;
984bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
985bf9b5d84SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called");
986bf9b5d84SToby Isaac   forest                 = (DM_Forest*) dm->data;
987bf9b5d84SToby Isaac   forest->computeAdaptSF = computeSF;
988bf9b5d84SToby Isaac   PetscFunctionReturn(0);
989bf9b5d84SToby Isaac }
990bf9b5d84SToby Isaac 
9910eb7e1eaSToby Isaac PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time)
99280b27e07SToby Isaac {
99380b27e07SToby Isaac   DM_Forest      *forest;
99480b27e07SToby Isaac   PetscErrorCode ierr;
99580b27e07SToby Isaac 
99680b27e07SToby Isaac   PetscFunctionBegin;
99780b27e07SToby Isaac   PetscValidHeaderSpecific(dmIn   ,DM_CLASSID  ,1);
99880b27e07SToby Isaac   PetscValidHeaderSpecific(vecIn  ,VEC_CLASSID ,2);
99980b27e07SToby Isaac   PetscValidHeaderSpecific(dmOut  ,DM_CLASSID  ,3);
100080b27e07SToby Isaac   PetscValidHeaderSpecific(vecOut ,VEC_CLASSID ,4);
100180b27e07SToby Isaac   forest = (DM_Forest *) dmIn->data;
100280b27e07SToby Isaac   if (!forest->transfervec) SETERRQ(PetscObjectComm((PetscObject)dmIn),PETSC_ERR_SUP,"DMForestTransferVec() not implemented");
10030eb7e1eaSToby Isaac   ierr = (forest->transfervec)(dmIn,vecIn,dmOut,vecOut,useBCs,time);CHKERRQ(ierr);
100480b27e07SToby Isaac   PetscFunctionReturn(0);
100580b27e07SToby Isaac }
100680b27e07SToby Isaac 
1007bf9b5d84SToby Isaac /*@
1008bf9b5d84SToby Isaac   DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the
1009bf9b5d84SToby Isaac   pre-adaptation forest to the post-adaptiation forest.  After DMSetUp() is called, these transfer PetscSFs can be
1010bf9b5d84SToby Isaac   accessed with DMForestGetAdaptivitySF().
1011bf9b5d84SToby Isaac 
1012bf9b5d84SToby Isaac   Not collective
1013bf9b5d84SToby Isaac 
1014bf9b5d84SToby Isaac   Input Parameter:
1015bf9b5d84SToby Isaac . dm - the post-adaptation forest
1016bf9b5d84SToby Isaac 
1017bf9b5d84SToby Isaac   Output Parameter:
1018bf9b5d84SToby Isaac . computeSF - default PETSC_TRUE
1019bf9b5d84SToby Isaac 
1020bf9b5d84SToby Isaac   Level: advanced
1021bf9b5d84SToby Isaac 
1022bf9b5d84SToby Isaac .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
1023bf9b5d84SToby Isaac @*/
1024bf9b5d84SToby Isaac PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF)
1025bf9b5d84SToby Isaac {
1026bf9b5d84SToby Isaac   DM_Forest *forest;
1027bf9b5d84SToby Isaac 
1028bf9b5d84SToby Isaac   PetscFunctionBegin;
1029bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1030bf9b5d84SToby Isaac   forest     = (DM_Forest*) dm->data;
1031bf9b5d84SToby Isaac   *computeSF = forest->computeAdaptSF;
1032bf9b5d84SToby Isaac   PetscFunctionReturn(0);
1033bf9b5d84SToby Isaac }
1034bf9b5d84SToby Isaac 
1035bf9b5d84SToby Isaac /*@
1036bf9b5d84SToby Isaac   DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest.
1037bf9b5d84SToby Isaac   Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be
1038bf9b5d84SToby Isaac   some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa.  Therefore there are two
1039bf9b5d84SToby Isaac   PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates
1040bf9b5d84SToby Isaac   pre-adaptation fine cells to post-adaptation coarse cells.
1041bf9b5d84SToby Isaac 
1042bf9b5d84SToby Isaac   Not collective
1043bf9b5d84SToby Isaac 
1044bf9b5d84SToby Isaac   Input Parameter:
1045bf9b5d84SToby Isaac   dm - the post-adaptation forest
1046bf9b5d84SToby Isaac 
1047bf9b5d84SToby Isaac   Output Parameter:
10480f17b9e3SToby Isaac   preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post-
10490f17b9e3SToby Isaac   coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre-
1050bf9b5d84SToby Isaac 
1051bf9b5d84SToby Isaac   Level: advanced
1052bf9b5d84SToby Isaac 
1053bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF()
1054bf9b5d84SToby Isaac @*/
10550f17b9e3SToby Isaac PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine)
1056bf9b5d84SToby Isaac {
1057bf9b5d84SToby Isaac   DM_Forest      *forest;
1058bf9b5d84SToby Isaac   PetscErrorCode ierr;
1059bf9b5d84SToby Isaac 
1060bf9b5d84SToby Isaac   PetscFunctionBegin;
1061bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1062bf9b5d84SToby Isaac   ierr   = DMSetUp(dm);CHKERRQ(ierr);
1063bf9b5d84SToby Isaac   forest = (DM_Forest*) dm->data;
1064f885a11aSToby Isaac   if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine;
1065f885a11aSToby Isaac   if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine;
1066bf9b5d84SToby Isaac   PetscFunctionReturn(0);
1067bf9b5d84SToby Isaac }
1068bf9b5d84SToby Isaac 
10699be51f97SToby Isaac /*@
10709be51f97SToby Isaac   DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to
10719be51f97SToby Isaac   indicate that the diameter of neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may
10729be51f97SToby Isaac   only support one particular choice of grading factor.
10739be51f97SToby Isaac 
10749be51f97SToby Isaac   Logically collective on dm
10759be51f97SToby Isaac 
10769be51f97SToby Isaac   Input Parameters:
10779be51f97SToby Isaac + dm - the forest
10789be51f97SToby Isaac - grade - the grading factor
10799be51f97SToby Isaac 
10809be51f97SToby Isaac   Level: advanced
10819be51f97SToby Isaac 
10829be51f97SToby Isaac .seealso: DMForestGetGradeFactor()
10839be51f97SToby Isaac @*/
1084c7eeac06SToby Isaac PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade)
1085c7eeac06SToby Isaac {
1086c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1087c7eeac06SToby Isaac 
1088c7eeac06SToby Isaac   PetscFunctionBegin;
1089c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1090ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup");
1091c7eeac06SToby Isaac   forest->gradeFactor = grade;
1092c7eeac06SToby Isaac   PetscFunctionReturn(0);
1093c7eeac06SToby Isaac }
1094c7eeac06SToby Isaac 
10959be51f97SToby Isaac /*@
10969be51f97SToby Isaac   DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of
10979be51f97SToby Isaac   neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may only support one particular
10989be51f97SToby Isaac   choice of grading factor.
10999be51f97SToby Isaac 
11009be51f97SToby Isaac   Not collective
11019be51f97SToby Isaac 
11029be51f97SToby Isaac   Input Parameter:
11039be51f97SToby Isaac . dm - the forest
11049be51f97SToby Isaac 
11059be51f97SToby Isaac   Output Parameter:
11069be51f97SToby Isaac . grade - the grading factor
11079be51f97SToby Isaac 
11089be51f97SToby Isaac   Level: advanced
11099be51f97SToby Isaac 
11109be51f97SToby Isaac .seealso: DMForestSetGradeFactor()
11119be51f97SToby Isaac @*/
1112c7eeac06SToby Isaac PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade)
1113c7eeac06SToby Isaac {
1114c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1115c7eeac06SToby Isaac 
1116c7eeac06SToby Isaac   PetscFunctionBegin;
1117c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1118c7eeac06SToby Isaac   PetscValidIntPointer(grade,2);
1119c7eeac06SToby Isaac   *grade = forest->gradeFactor;
1120c7eeac06SToby Isaac   PetscFunctionReturn(0);
1121c7eeac06SToby Isaac }
1122c7eeac06SToby Isaac 
11239be51f97SToby Isaac /*@
11249be51f97SToby Isaac   DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes
11259be51f97SToby Isaac   the cell weight (see DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be
11269be51f97SToby Isaac   (cellWeight) * (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on
11279be51f97SToby Isaac   its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the
11289be51f97SToby Isaac   computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
11299be51f97SToby Isaac 
11309be51f97SToby Isaac   Logically collective on dm
11319be51f97SToby Isaac 
11329be51f97SToby Isaac   Input Parameters:
11339be51f97SToby Isaac + dm - the forest
11349be51f97SToby Isaac - weightsFactors - default 1.
11359be51f97SToby Isaac 
11369be51f97SToby Isaac   Level: advanced
11379be51f97SToby Isaac 
11389be51f97SToby Isaac .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights()
11399be51f97SToby Isaac @*/
1140ef51cf95SToby Isaac PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor)
1141c7eeac06SToby Isaac {
1142c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1143c7eeac06SToby Isaac 
1144c7eeac06SToby Isaac   PetscFunctionBegin;
1145c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1146ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup");
1147c7eeac06SToby Isaac   forest->weightsFactor = weightsFactor;
1148c7eeac06SToby Isaac   PetscFunctionReturn(0);
1149c7eeac06SToby Isaac }
1150c7eeac06SToby Isaac 
11519be51f97SToby Isaac /*@
11529be51f97SToby Isaac   DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see
11539be51f97SToby Isaac   DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be (cellWeight) *
11549be51f97SToby Isaac   (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on its level; a
11559be51f97SToby Isaac   factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation
11569be51f97SToby Isaac   associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
11579be51f97SToby Isaac 
11589be51f97SToby Isaac   Not collective
11599be51f97SToby Isaac 
11609be51f97SToby Isaac   Input Parameter:
11619be51f97SToby Isaac . dm - the forest
11629be51f97SToby Isaac 
11639be51f97SToby Isaac   Output Parameter:
11649be51f97SToby Isaac . weightsFactors - default 1.
11659be51f97SToby Isaac 
11669be51f97SToby Isaac   Level: advanced
11679be51f97SToby Isaac 
11689be51f97SToby Isaac .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights()
11699be51f97SToby Isaac @*/
1170ef51cf95SToby Isaac PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor)
1171c7eeac06SToby Isaac {
1172c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1173c7eeac06SToby Isaac 
1174c7eeac06SToby Isaac   PetscFunctionBegin;
1175c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1176c7eeac06SToby Isaac   PetscValidRealPointer(weightsFactor,2);
1177c7eeac06SToby Isaac   *weightsFactor = forest->weightsFactor;
1178c7eeac06SToby Isaac   PetscFunctionReturn(0);
1179c7eeac06SToby Isaac }
1180c7eeac06SToby Isaac 
11819be51f97SToby Isaac /*@
11829be51f97SToby Isaac   DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process
11839be51f97SToby Isaac 
11849be51f97SToby Isaac   Not collective
11859be51f97SToby Isaac 
11869be51f97SToby Isaac   Input Parameter:
11879be51f97SToby Isaac . dm - the forest
11889be51f97SToby Isaac 
11899be51f97SToby Isaac   Output Parameters:
11909be51f97SToby Isaac + cStart - the first cell on this process
11919be51f97SToby Isaac - cEnd - one after the final cell on this process
11929be51f97SToby Isaac 
11931a244344SSatish Balay   Level: intermediate
11949be51f97SToby Isaac 
11959be51f97SToby Isaac .seealso: DMForestGetCellSF()
11969be51f97SToby Isaac @*/
1197c7eeac06SToby Isaac PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd)
1198c7eeac06SToby Isaac {
1199c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1200c7eeac06SToby Isaac   PetscErrorCode ierr;
1201c7eeac06SToby Isaac 
1202c7eeac06SToby Isaac   PetscFunctionBegin;
1203c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1204c7eeac06SToby Isaac   PetscValidIntPointer(cStart,2);
1205c7eeac06SToby Isaac   PetscValidIntPointer(cEnd,2);
1206c7eeac06SToby Isaac   if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) {
1207c7eeac06SToby Isaac     ierr = forest->createcellchart(dm,&forest->cStart,&forest->cEnd);CHKERRQ(ierr);
1208c7eeac06SToby Isaac   }
1209c7eeac06SToby Isaac   *cStart =  forest->cStart;
1210c7eeac06SToby Isaac   *cEnd   =  forest->cEnd;
1211c7eeac06SToby Isaac   PetscFunctionReturn(0);
1212c7eeac06SToby Isaac }
1213c7eeac06SToby Isaac 
12149be51f97SToby Isaac /*@
12159be51f97SToby Isaac   DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes
12169be51f97SToby Isaac 
12179be51f97SToby Isaac   Not collective
12189be51f97SToby Isaac 
12199be51f97SToby Isaac   Input Parameter:
12209be51f97SToby Isaac . dm - the forest
12219be51f97SToby Isaac 
12229be51f97SToby Isaac   Output Parameter:
12239be51f97SToby Isaac . cellSF - the PetscSF
12249be51f97SToby Isaac 
12251a244344SSatish Balay   Level: intermediate
12269be51f97SToby Isaac 
12279be51f97SToby Isaac .seealso: DMForestGetCellChart()
12289be51f97SToby Isaac @*/
1229c7eeac06SToby Isaac PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF)
1230c7eeac06SToby Isaac {
1231c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1232c7eeac06SToby Isaac   PetscErrorCode ierr;
1233c7eeac06SToby Isaac 
1234c7eeac06SToby Isaac   PetscFunctionBegin;
1235c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1236c7eeac06SToby Isaac   PetscValidPointer(cellSF,2);
1237c7eeac06SToby Isaac   if ((!forest->cellSF) && forest->createcellsf) {
1238c7eeac06SToby Isaac     ierr = forest->createcellsf(dm,&forest->cellSF);CHKERRQ(ierr);
1239c7eeac06SToby Isaac   }
1240c7eeac06SToby Isaac   *cellSF = forest->cellSF;
1241c7eeac06SToby Isaac   PetscFunctionReturn(0);
1242c7eeac06SToby Isaac }
1243c7eeac06SToby Isaac 
12449be51f97SToby Isaac /*@C
12459be51f97SToby Isaac   DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see
12469be51f97SToby Isaac   DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination).  The
1247cd3c525cSToby Isaac   interpretation of the label values is up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP,
1248cd3c525cSToby Isaac   DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have been reserved as choices that should be accepted by all subtypes.
12499be51f97SToby Isaac 
12509be51f97SToby Isaac   Logically collective on dm
12519be51f97SToby Isaac 
12529be51f97SToby Isaac   Input Parameters:
12539be51f97SToby Isaac - dm - the forest
1254a1b0c543SToby Isaac + adaptLabel - the label in the pre-adaptation forest
12559be51f97SToby Isaac 
12569be51f97SToby Isaac   Level: intermediate
12579be51f97SToby Isaac 
12589be51f97SToby Isaac .seealso DMForestGetAdaptivityLabel()
12599be51f97SToby Isaac @*/
1260a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityLabel(DM dm, DMLabel adaptLabel)
1261c7eeac06SToby Isaac {
1262c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1263c7eeac06SToby Isaac   PetscErrorCode ierr;
1264c7eeac06SToby Isaac 
1265c7eeac06SToby Isaac   PetscFunctionBegin;
1266c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1267a1b0c543SToby Isaac   adaptLabel->refct++;
1268a1b0c543SToby Isaac   if (forest->adaptLabel) {ierr = DMLabelDestroy(&forest->adaptLabel);CHKERRQ(ierr);}
1269a1b0c543SToby Isaac   forest->adaptLabel = adaptLabel;
1270c7eeac06SToby Isaac   PetscFunctionReturn(0);
1271c7eeac06SToby Isaac }
1272c7eeac06SToby Isaac 
12739be51f97SToby Isaac /*@C
12749be51f97SToby Isaac   DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that
12759be51f97SToby Isaac   holds the adaptation flags (refinement, coarsening, or some combination).  The interpretation of the label values is
1276cd3c525cSToby Isaac   up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have
1277cd3c525cSToby Isaac   been reserved as choices that should be accepted by all subtypes.
12789be51f97SToby Isaac 
12799be51f97SToby Isaac   Not collective
12809be51f97SToby Isaac 
12819be51f97SToby Isaac   Input Parameter:
12829be51f97SToby Isaac . dm - the forest
12839be51f97SToby Isaac 
12849be51f97SToby Isaac   Output Parameter:
12859be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest
12869be51f97SToby Isaac 
12879be51f97SToby Isaac   Level: intermediate
12889be51f97SToby Isaac 
12899be51f97SToby Isaac .seealso DMForestSetAdaptivityLabel()
12909be51f97SToby Isaac @*/
1291a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityLabel(DM dm, DMLabel *adaptLabel)
1292c7eeac06SToby Isaac {
1293c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1294c7eeac06SToby Isaac 
1295c7eeac06SToby Isaac   PetscFunctionBegin;
1296c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1297ba936b91SToby Isaac   *adaptLabel = forest->adaptLabel;
1298c7eeac06SToby Isaac   PetscFunctionReturn(0);
1299c7eeac06SToby Isaac }
1300c7eeac06SToby Isaac 
13019be51f97SToby Isaac /*@
13029be51f97SToby Isaac   DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
13039be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
13049be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
13059be51f97SToby Isaac   of 1.
13069be51f97SToby Isaac 
13079be51f97SToby Isaac   Logically collective on dm
13089be51f97SToby Isaac 
13099be51f97SToby Isaac   Input Parameters:
13109be51f97SToby Isaac + dm - the forest
13119be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
13129be51f97SToby Isaac - copyMode - how weights should reference weights
13139be51f97SToby Isaac 
13149be51f97SToby Isaac   Level: advanced
13159be51f97SToby Isaac 
13169be51f97SToby Isaac .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity()
13179be51f97SToby Isaac @*/
1318c7eeac06SToby Isaac PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode)
1319c7eeac06SToby Isaac {
1320c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1321c7eeac06SToby Isaac   PetscInt       cStart, cEnd;
1322c7eeac06SToby Isaac   PetscErrorCode ierr;
1323c7eeac06SToby Isaac 
1324c7eeac06SToby Isaac   PetscFunctionBegin;
1325c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1326c7eeac06SToby Isaac   ierr = DMForestGetCellChart(dm,&cStart,&cEnd);CHKERRQ(ierr);
1327c7eeac06SToby Isaac   if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd);
1328c7eeac06SToby Isaac   if (copyMode == PETSC_COPY_VALUES) {
1329c7eeac06SToby Isaac     if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) {
1330c7eeac06SToby Isaac       ierr = PetscMalloc1(cEnd-cStart,&forest->cellWeights);CHKERRQ(ierr);
1331c7eeac06SToby Isaac     }
1332c7eeac06SToby Isaac     ierr                        = PetscMemcpy(forest->cellWeights,weights,(cEnd-cStart)*sizeof(*weights));CHKERRQ(ierr);
1333c7eeac06SToby Isaac     forest->cellWeightsCopyMode = PETSC_OWN_POINTER;
1334c7eeac06SToby Isaac     PetscFunctionReturn(0);
1335c7eeac06SToby Isaac   }
1336c7eeac06SToby Isaac   if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) {
1337c7eeac06SToby Isaac     ierr = PetscFree(forest->cellWeights);CHKERRQ(ierr);
1338c7eeac06SToby Isaac   }
1339c7eeac06SToby Isaac   forest->cellWeights         = weights;
1340c7eeac06SToby Isaac   forest->cellWeightsCopyMode = copyMode;
1341c7eeac06SToby Isaac   PetscFunctionReturn(0);
1342c7eeac06SToby Isaac }
1343c7eeac06SToby Isaac 
13449be51f97SToby Isaac /*@
13459be51f97SToby Isaac   DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
13469be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
13479be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
13489be51f97SToby Isaac   of 1.
13499be51f97SToby Isaac 
13509be51f97SToby Isaac   Not collective
13519be51f97SToby Isaac 
13529be51f97SToby Isaac   Input Parameter:
13539be51f97SToby Isaac . dm - the forest
13549be51f97SToby Isaac 
13559be51f97SToby Isaac   Output Parameter:
13569be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
13579be51f97SToby Isaac 
13589be51f97SToby Isaac   Level: advanced
13599be51f97SToby Isaac 
13609be51f97SToby Isaac .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity()
13619be51f97SToby Isaac @*/
1362c7eeac06SToby Isaac PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights)
1363c7eeac06SToby Isaac {
1364c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1365c7eeac06SToby Isaac 
1366c7eeac06SToby Isaac   PetscFunctionBegin;
1367c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1368c7eeac06SToby Isaac   PetscValidPointer(weights,2);
1369c7eeac06SToby Isaac   *weights = forest->cellWeights;
1370c7eeac06SToby Isaac   PetscFunctionReturn(0);
1371c7eeac06SToby Isaac }
1372c7eeac06SToby Isaac 
13739be51f97SToby Isaac /*@
13749be51f97SToby Isaac   DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning
13759be51f97SToby Isaac   a pre-adaptation forest (see DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each
13769be51f97SToby Isaac   process's cells to the process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that
13779be51f97SToby Isaac   the current process should not have any cells after repartitioning.
13789be51f97SToby Isaac 
13799be51f97SToby Isaac   Logically Collective on dm
13809be51f97SToby Isaac 
13819be51f97SToby Isaac   Input parameters:
13829be51f97SToby Isaac + dm - the forest
13839be51f97SToby Isaac - capacity - this process's capacity
13849be51f97SToby Isaac 
13859be51f97SToby Isaac   Level: advanced
13869be51f97SToby Isaac 
13879be51f97SToby Isaac .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
13889be51f97SToby Isaac @*/
1389c7eeac06SToby Isaac PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity)
1390c7eeac06SToby Isaac {
1391c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1392c7eeac06SToby Isaac 
1393c7eeac06SToby Isaac   PetscFunctionBegin;
1394c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1395ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup");
1396c7eeac06SToby Isaac   if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity);
1397c7eeac06SToby Isaac   forest->weightCapacity = capacity;
1398c7eeac06SToby Isaac   PetscFunctionReturn(0);
1399c7eeac06SToby Isaac }
1400c7eeac06SToby Isaac 
14019be51f97SToby Isaac /*@
14029be51f97SToby Isaac   DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see
14039be51f97SToby Isaac   DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each process's cells to the
14049be51f97SToby Isaac   process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that the current process
14059be51f97SToby Isaac   should not have any cells after repartitioning.
14069be51f97SToby Isaac 
14079be51f97SToby Isaac   Not collective
14089be51f97SToby Isaac 
14099be51f97SToby Isaac   Input parameter:
14109be51f97SToby Isaac . dm - the forest
14119be51f97SToby Isaac 
14129be51f97SToby Isaac   Output parameter:
14139be51f97SToby Isaac . capacity - this process's capacity
14149be51f97SToby Isaac 
14159be51f97SToby Isaac   Level: advanced
14169be51f97SToby Isaac 
14179be51f97SToby Isaac .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
14189be51f97SToby Isaac @*/
1419c7eeac06SToby Isaac PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity)
1420c7eeac06SToby Isaac {
1421c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1422c7eeac06SToby Isaac 
1423c7eeac06SToby Isaac   PetscFunctionBegin;
1424c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1425c7eeac06SToby Isaac   PetscValidRealPointer(capacity,2);
1426c7eeac06SToby Isaac   *capacity = forest->weightCapacity;
1427c7eeac06SToby Isaac   PetscFunctionReturn(0);
1428c7eeac06SToby Isaac }
1429c7eeac06SToby Isaac 
14300709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm)
1431db4d5e8cSToby Isaac {
1432db4d5e8cSToby Isaac   DM_Forest                  *forest = (DM_Forest*) dm->data;
143356ba9f64SToby Isaac   PetscBool                  flg, flg1, flg2, flg3, flg4;
1434dd8e54a2SToby Isaac   DMForestTopology           oldTopo;
1435c7eeac06SToby Isaac   char                       stringBuffer[256];
1436dd8e54a2SToby Isaac   PetscViewer                viewer;
1437dd8e54a2SToby Isaac   PetscViewerFormat          format;
143856ba9f64SToby Isaac   PetscInt                   adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade;
1439c7eeac06SToby Isaac   PetscReal                  weightsFactor;
1440c7eeac06SToby Isaac   DMForestAdaptivityStrategy adaptStrategy;
1441db4d5e8cSToby Isaac   PetscErrorCode             ierr;
1442db4d5e8cSToby Isaac 
1443db4d5e8cSToby Isaac   PetscFunctionBegin;
1444db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
144558762b62SToby Isaac   forest->setfromoptionscalled = PETSC_TRUE;
1446dd8e54a2SToby Isaac   ierr                         = DMForestGetTopology(dm, &oldTopo);CHKERRQ(ierr);
1447a3eda1eaSToby Isaac   ierr                         = PetscOptionsHead(PetscOptionsObject,"DMForest Options");CHKERRQ(ierr);
144856ba9f64SToby Isaac   ierr                         = PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);CHKERRQ(ierr);
144956ba9f64SToby Isaac   ierr                         = PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);CHKERRQ(ierr);
145056ba9f64SToby Isaac   ierr                         = PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);CHKERRQ(ierr);
145156ba9f64SToby Isaac   ierr                         = PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);CHKERRQ(ierr);
1452f885a11aSToby 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}");
145356ba9f64SToby Isaac   if (flg1) {
145456ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,(DMForestTopology)stringBuffer);CHKERRQ(ierr);
145556ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
145620e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
145756ba9f64SToby Isaac   }
145856ba9f64SToby Isaac   if (flg2) {
1459dd8e54a2SToby Isaac     DM base;
1460dd8e54a2SToby Isaac 
1461dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&base);CHKERRQ(ierr);
1462dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1463dd8e54a2SToby Isaac     ierr = DMLoad(base,viewer);CHKERRQ(ierr);
1464dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1465dd8e54a2SToby Isaac     ierr = DMForestSetBaseDM(dm,base);CHKERRQ(ierr);
1466dd8e54a2SToby Isaac     ierr = DMDestroy(&base);CHKERRQ(ierr);
146756ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
146820e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
1469dd8e54a2SToby Isaac   }
147056ba9f64SToby Isaac   if (flg3) {
1471dd8e54a2SToby Isaac     DM coarse;
1472dd8e54a2SToby Isaac 
1473dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&coarse);CHKERRQ(ierr);
1474dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1475dd8e54a2SToby Isaac     ierr = DMLoad(coarse,viewer);CHKERRQ(ierr);
1476dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
147720e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,coarse);CHKERRQ(ierr);
1478dd8e54a2SToby Isaac     ierr = DMDestroy(&coarse);CHKERRQ(ierr);
147956ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
148056ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1481dd8e54a2SToby Isaac   }
148256ba9f64SToby Isaac   if (flg4) {
1483dd8e54a2SToby Isaac     DM fine;
1484dd8e54a2SToby Isaac 
1485dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&fine);CHKERRQ(ierr);
1486dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1487dd8e54a2SToby Isaac     ierr = DMLoad(fine,viewer);CHKERRQ(ierr);
1488dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
148920e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,fine);CHKERRQ(ierr);
1490dd8e54a2SToby Isaac     ierr = DMDestroy(&fine);CHKERRQ(ierr);
149156ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
149256ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1493dd8e54a2SToby Isaac   }
1494dd8e54a2SToby Isaac   ierr = DMForestGetAdjacencyDimension(dm,&adjDim);CHKERRQ(ierr);
1495dd8e54a2SToby Isaac   ierr = PetscOptionsInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg);CHKERRQ(ierr);
1496dd8e54a2SToby Isaac   if (flg) {
1497dd8e54a2SToby Isaac     ierr = DMForestSetAdjacencyDimension(dm,adjDim);CHKERRQ(ierr);
1498f885a11aSToby Isaac   } else {
1499dd8e54a2SToby Isaac     ierr = DMForestGetAdjacencyCodimension(dm,&adjCodim);CHKERRQ(ierr);
1500dd8e54a2SToby Isaac     ierr = PetscOptionsInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg);CHKERRQ(ierr);
1501dd8e54a2SToby Isaac     if (flg) {
1502dd8e54a2SToby Isaac       ierr = DMForestSetAdjacencyCodimension(dm,adjCodim);CHKERRQ(ierr);
1503dd8e54a2SToby Isaac     }
1504dd8e54a2SToby Isaac   }
1505dd8e54a2SToby Isaac   ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr);
1506dd8e54a2SToby Isaac   ierr = PetscOptionsInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg);CHKERRQ(ierr);
1507dd8e54a2SToby Isaac   if (flg) {
1508dd8e54a2SToby Isaac     ierr = DMForestSetPartitionOverlap(dm,overlap);CHKERRQ(ierr);
1509dd8e54a2SToby Isaac   }
1510a6121fbdSMatthew G. Knepley #if 0
1511a6121fbdSMatthew 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);
1512a6121fbdSMatthew G. Knepley   if (flg) {
1513a6121fbdSMatthew G. Knepley     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1514a6121fbdSMatthew G. Knepley     ierr = DMForestSetInitialRefinement(dm,minRefinement);CHKERRQ(ierr);
1515a6121fbdSMatthew G. Knepley   }
1516a6121fbdSMatthew 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);
1517a6121fbdSMatthew G. Knepley   if (flg) {
1518a6121fbdSMatthew G. Knepley     ierr = DMForestSetMinimumRefinement(dm,0);CHKERRQ(ierr);
1519a6121fbdSMatthew G. Knepley     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
1520a6121fbdSMatthew G. Knepley   }
1521a6121fbdSMatthew G. Knepley #endif
1522dd8e54a2SToby Isaac   ierr = DMForestGetMinimumRefinement(dm,&minRefinement);CHKERRQ(ierr);
1523dd8e54a2SToby Isaac   ierr = PetscOptionsInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg);CHKERRQ(ierr);
1524dd8e54a2SToby Isaac   if (flg) {
1525dd8e54a2SToby Isaac     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1526db4d5e8cSToby Isaac   }
152756ba9f64SToby Isaac   ierr = DMForestGetInitialRefinement(dm,&initRefinement);CHKERRQ(ierr);
152856ba9f64SToby Isaac   ierr = PetscOptionsInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg);CHKERRQ(ierr);
152956ba9f64SToby Isaac   if (flg) {
153056ba9f64SToby Isaac     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
153156ba9f64SToby Isaac   }
1532c7eeac06SToby Isaac   ierr = DMForestGetMaximumRefinement(dm,&maxRefinement);CHKERRQ(ierr);
1533c7eeac06SToby Isaac   ierr = PetscOptionsInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg);CHKERRQ(ierr);
1534c7eeac06SToby Isaac   if (flg) {
1535c7eeac06SToby Isaac     ierr = DMForestSetMaximumRefinement(dm,maxRefinement);CHKERRQ(ierr);
1536c7eeac06SToby Isaac   }
1537c7eeac06SToby Isaac   ierr = DMForestGetAdaptivityStrategy(dm,&adaptStrategy);CHKERRQ(ierr);
1538c7eeac06SToby Isaac   ierr = PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,256,&flg);CHKERRQ(ierr);
1539c7eeac06SToby Isaac   if (flg) {
1540c7eeac06SToby Isaac     ierr = DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);CHKERRQ(ierr);
1541c7eeac06SToby Isaac   }
1542c7eeac06SToby Isaac   ierr = DMForestGetGradeFactor(dm,&grade);CHKERRQ(ierr);
1543c7eeac06SToby Isaac   ierr = PetscOptionsInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg);CHKERRQ(ierr);
1544c7eeac06SToby Isaac   if (flg) {
1545c7eeac06SToby Isaac     ierr = DMForestSetGradeFactor(dm,grade);CHKERRQ(ierr);
1546c7eeac06SToby Isaac   }
1547c7eeac06SToby Isaac   ierr = DMForestGetCellWeightFactor(dm,&weightsFactor);CHKERRQ(ierr);
1548c7eeac06SToby Isaac   ierr = PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);CHKERRQ(ierr);
1549c7eeac06SToby Isaac   if (flg) {
1550c7eeac06SToby Isaac     ierr = DMForestSetCellWeightFactor(dm,weightsFactor);CHKERRQ(ierr);
1551c7eeac06SToby Isaac   }
1552db4d5e8cSToby Isaac   ierr = PetscOptionsTail();CHKERRQ(ierr);
1553db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1554db4d5e8cSToby Isaac }
1555db4d5e8cSToby Isaac 
1556276c5506SMatthew G. Knepley PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1557d8984e3bSMatthew G. Knepley {
1558d8984e3bSMatthew G. Knepley   PetscErrorCode ierr;
1559d8984e3bSMatthew G. Knepley 
1560d8984e3bSMatthew G. Knepley   PetscFunctionBegin;
1561d8984e3bSMatthew G. Knepley   if (subdm) {ierr = DMClone(dm, subdm);CHKERRQ(ierr);}
1562d8984e3bSMatthew G. Knepley   ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr);
1563d8984e3bSMatthew G. Knepley   PetscFunctionReturn(0);
1564d8984e3bSMatthew G. Knepley }
1565d8984e3bSMatthew G. Knepley 
15665421bac9SToby Isaac PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
15675421bac9SToby Isaac {
15685421bac9SToby Isaac   DMLabel        refine;
15695421bac9SToby Isaac   DM             fineDM;
15705421bac9SToby Isaac   PetscErrorCode ierr;
15715421bac9SToby Isaac 
15725421bac9SToby Isaac   PetscFunctionBegin;
15735421bac9SToby Isaac   ierr = DMGetFineDM(dm,&fineDM);CHKERRQ(ierr);
15745421bac9SToby Isaac   if (fineDM) {
15755421bac9SToby Isaac     ierr       = PetscObjectReference((PetscObject)fineDM);CHKERRQ(ierr);
15765421bac9SToby Isaac     *dmRefined = fineDM;
15775421bac9SToby Isaac     PetscFunctionReturn(0);
15785421bac9SToby Isaac   }
15795421bac9SToby Isaac   ierr = DMForestTemplate(dm,comm,dmRefined);CHKERRQ(ierr);
15805421bac9SToby Isaac   ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr);
15815421bac9SToby Isaac   if (!refine) {
1582a1b0c543SToby Isaac     ierr = DMLabelCreate("refine",&refine);CHKERRQ(ierr);
1583a1b0c543SToby Isaac     ierr = DMLabelSetDefaultValue(refine,DM_ADAPT_REFINE);CHKERRQ(ierr);
15845421bac9SToby Isaac   }
1585a1b0c543SToby Isaac   else {
1586a1b0c543SToby Isaac     refine->refct++;
1587a1b0c543SToby Isaac   }
1588a1b0c543SToby Isaac   ierr = DMForestSetAdaptivityLabel(*dmRefined,refine);CHKERRQ(ierr);
1589a1b0c543SToby Isaac   ierr = DMLabelDestroy(&refine);CHKERRQ(ierr);
15905421bac9SToby Isaac   PetscFunctionReturn(0);
15915421bac9SToby Isaac }
15925421bac9SToby Isaac 
15935421bac9SToby Isaac PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
15945421bac9SToby Isaac {
15955421bac9SToby Isaac   DMLabel        coarsen;
15965421bac9SToby Isaac   DM             coarseDM;
15975421bac9SToby Isaac   PetscErrorCode ierr;
15985421bac9SToby Isaac 
15995421bac9SToby Isaac   PetscFunctionBegin;
16004098eed7SToby Isaac   {
16014098eed7SToby Isaac     PetscMPIInt mpiComparison;
16024098eed7SToby Isaac     MPI_Comm    dmcomm = PetscObjectComm((PetscObject)dm);
16034098eed7SToby Isaac 
16044098eed7SToby Isaac     ierr = MPI_Comm_compare(comm, dmcomm, &mpiComparison);CHKERRQ(ierr);
1605f885a11aSToby Isaac     if (mpiComparison != MPI_IDENT && mpiComparison != MPI_CONGRUENT) SETERRQ(dmcomm,PETSC_ERR_SUP,"No support for different communicators yet");
16064098eed7SToby Isaac   }
16075421bac9SToby Isaac   ierr = DMGetCoarseDM(dm,&coarseDM);CHKERRQ(ierr);
16085421bac9SToby Isaac   if (coarseDM) {
16095421bac9SToby Isaac     ierr         = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr);
16105421bac9SToby Isaac     *dmCoarsened = coarseDM;
16115421bac9SToby Isaac     PetscFunctionReturn(0);
16125421bac9SToby Isaac   }
16135421bac9SToby Isaac   ierr = DMForestTemplate(dm,comm,dmCoarsened);CHKERRQ(ierr);
16147ee90a76SStefano Zampini   ierr = DMForestSetAdaptivityPurpose(*dmCoarsened,DM_ADAPT_COARSEN);CHKERRQ(ierr);
16155421bac9SToby Isaac   ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr);
16165421bac9SToby Isaac   if (!coarsen) {
1617a1b0c543SToby Isaac     ierr = DMLabelCreate("coarsen",&coarsen);CHKERRQ(ierr);
1618a1b0c543SToby Isaac     ierr = DMLabelSetDefaultValue(coarsen,DM_ADAPT_COARSEN);CHKERRQ(ierr);
1619a1b0c543SToby Isaac   } else {
1620a1b0c543SToby Isaac     coarsen->refct++;
16215421bac9SToby Isaac   }
1622a1b0c543SToby Isaac   ierr = DMForestSetAdaptivityLabel(*dmCoarsened,coarsen);CHKERRQ(ierr);
1623a1b0c543SToby Isaac   ierr = DMLabelDestroy(&coarsen);CHKERRQ(ierr);
16245421bac9SToby Isaac   PetscFunctionReturn(0);
16255421bac9SToby Isaac }
16265421bac9SToby Isaac 
1627a1b0c543SToby Isaac static PetscErrorCode DMAdaptLabel_Forest(DM dm, DMLabel label, DM *adaptedDM)
162809350103SToby Isaac {
162909350103SToby Isaac   PetscBool      success;
163009350103SToby Isaac   PetscErrorCode ierr;
163109350103SToby Isaac 
163209350103SToby Isaac   PetscFunctionBegin;
163309350103SToby Isaac   ierr = DMForestTemplate(dm,PetscObjectComm((PetscObject)dm),adaptedDM);CHKERRQ(ierr);
1634a1b0c543SToby Isaac   ierr = DMForestSetAdaptivityLabel(*adaptedDM,label);CHKERRQ(ierr);
163509350103SToby Isaac   ierr = DMSetUp(*adaptedDM);CHKERRQ(ierr);
163609350103SToby Isaac   ierr = DMForestGetAdaptivitySuccess(*adaptedDM,&success);CHKERRQ(ierr);
163709350103SToby Isaac   if (!success) {
163809350103SToby Isaac     ierr = DMDestroy(adaptedDM);CHKERRQ(ierr);
163909350103SToby Isaac     *adaptedDM = NULL;
164009350103SToby Isaac   }
164109350103SToby Isaac   PetscFunctionReturn(0);
164209350103SToby Isaac }
164309350103SToby Isaac 
1644d222f98bSToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm)
1645d222f98bSToby Isaac {
1646d222f98bSToby Isaac   PetscErrorCode ierr;
1647d222f98bSToby Isaac 
1648d222f98bSToby Isaac   PetscFunctionBegin;
1649d222f98bSToby Isaac   ierr = PetscMemzero(dm->ops,sizeof(*(dm->ops)));CHKERRQ(ierr);
1650d222f98bSToby Isaac 
1651d222f98bSToby Isaac   dm->ops->clone          = DMClone_Forest;
1652d222f98bSToby Isaac   dm->ops->setfromoptions = DMSetFromOptions_Forest;
1653d222f98bSToby Isaac   dm->ops->destroy        = DMDestroy_Forest;
1654d8984e3bSMatthew G. Knepley   dm->ops->createsubdm    = DMCreateSubDM_Forest;
16555421bac9SToby Isaac   dm->ops->refine         = DMRefine_Forest;
16565421bac9SToby Isaac   dm->ops->coarsen        = DMCoarsen_Forest;
16570d1cd5e0SMatthew G. Knepley   dm->ops->adaptlabel     = DMAdaptLabel_Forest;
1658d222f98bSToby Isaac   PetscFunctionReturn(0);
1659d222f98bSToby Isaac }
1660d222f98bSToby Isaac 
16619be51f97SToby Isaac /*MC
16629be51f97SToby Isaac 
1663bae1f979SBarry Smith      DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh.  Forests usually have a base DM
1664bae1f979SBarry Smith   (see DMForestGetBaseDM()), from which it is refined.  The refinement and partitioning of forests is considered
1665bae1f979SBarry Smith   immutable after DMSetUp() is called.  To adapt a mesh, one should call DMForestTemplate() to create a new mesh that
1666bae1f979SBarry Smith   will default to being identical to it, specify how that mesh should differ, and then calling DMSetUp() on the new
1667bae1f979SBarry Smith   mesh.
1668bae1f979SBarry Smith 
1669bae1f979SBarry Smith   To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the
1670bae1f979SBarry Smith   previous mesh whose values indicate which cells should be refined (DM_ADAPT_REFINE) or coarsened (DM_ADAPT_COARSEN)
1671bae1f979SBarry Smith   and how (subtypes are free to allow additional values for things like anisotropic refinement).  The label should be
1672bae1f979SBarry Smith   given to the *new* mesh with DMForestSetAdaptivityLabel().
16739be51f97SToby Isaac 
16749be51f97SToby Isaac   Level: advanced
16759be51f97SToby Isaac 
16769be51f97SToby Isaac .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel()
16779be51f97SToby Isaac M*/
16789be51f97SToby Isaac 
1679db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1680db4d5e8cSToby Isaac {
1681db4d5e8cSToby Isaac   DM_Forest      *forest;
1682db4d5e8cSToby Isaac   PetscErrorCode ierr;
1683db4d5e8cSToby Isaac 
1684db4d5e8cSToby Isaac   PetscFunctionBegin;
1685db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1686db4d5e8cSToby Isaac   ierr                         = PetscNewLog(dm,&forest);CHKERRQ(ierr);
1687db4d5e8cSToby Isaac   dm->dim                      = 0;
1688db4d5e8cSToby Isaac   dm->data                     = forest;
1689db4d5e8cSToby Isaac   forest->refct                = 1;
1690db4d5e8cSToby Isaac   forest->data                 = NULL;
169158762b62SToby Isaac   forest->setfromoptionscalled = PETSC_FALSE;
1692db4d5e8cSToby Isaac   forest->topology             = NULL;
1693cd3c525cSToby Isaac   forest->adapt                = NULL;
1694db4d5e8cSToby Isaac   forest->base                 = NULL;
16956a87ffbfSToby Isaac   forest->adaptPurpose         = DM_ADAPT_DETERMINE;
1696db4d5e8cSToby Isaac   forest->adjDim               = PETSC_DEFAULT;
1697db4d5e8cSToby Isaac   forest->overlap              = PETSC_DEFAULT;
1698db4d5e8cSToby Isaac   forest->minRefinement        = PETSC_DEFAULT;
1699db4d5e8cSToby Isaac   forest->maxRefinement        = PETSC_DEFAULT;
170056ba9f64SToby Isaac   forest->initRefinement       = PETSC_DEFAULT;
1701c7eeac06SToby Isaac   forest->cStart               = PETSC_DETERMINE;
1702c7eeac06SToby Isaac   forest->cEnd                 = PETSC_DETERMINE;
1703cd3c525cSToby Isaac   forest->cellSF               = NULL;
1704ebdf65a2SToby Isaac   forest->adaptLabel           = NULL;
1705db4d5e8cSToby Isaac   forest->gradeFactor          = 2;
1706db4d5e8cSToby Isaac   forest->cellWeights          = NULL;
1707db4d5e8cSToby Isaac   forest->cellWeightsCopyMode  = PETSC_USE_POINTER;
1708db4d5e8cSToby Isaac   forest->weightsFactor        = 1.;
1709db4d5e8cSToby Isaac   forest->weightCapacity       = 1.;
1710a73e2921SToby Isaac   ierr                         = DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);CHKERRQ(ierr);
1711d222f98bSToby Isaac   ierr                         = DMInitialize_Forest(dm);CHKERRQ(ierr);
1712db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1713db4d5e8cSToby Isaac }
1714