xref: /petsc/src/dm/impls/forest/forest.c (revision bf2d5fbb3a3ef6516951834f394af24d0a387002)
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) {
30728dfcf7cSStefano Zampini     PetscBool        isper;
30828dfcf7cSStefano Zampini     const PetscReal *maxCell, *L;
30928dfcf7cSStefano Zampini     const DMBoundaryType *bd;
31028dfcf7cSStefano 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);
31628dfcf7cSStefano Zampini     ierr = DMGetPeriodicity(base,&isper,&maxCell,&L,&bd);CHKERRQ(ierr);
31728dfcf7cSStefano Zampini     ierr = DMSetPeriodicity(dm,isper,maxCell,L,bd);CHKERRQ(ierr);
31828dfcf7cSStefano Zampini   } else {
31928dfcf7cSStefano 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:
428*bf2d5fbbSStefano Zampini   case DM_ADAPT_COARSEN_LAST:
42926d9498aSToby Isaac     ierr = DMSetFineDM(dm,adapt);CHKERRQ(ierr);
43026d9498aSToby Isaac     break;
43126d9498aSToby Isaac   default:
43226d9498aSToby Isaac     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
43326d9498aSToby Isaac   }
434dd8e54a2SToby Isaac   PetscFunctionReturn(0);
435dd8e54a2SToby Isaac }
436dd8e54a2SToby Isaac 
4379be51f97SToby Isaac /*@
4389be51f97SToby Isaac   DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted.
4399be51f97SToby Isaac 
4409be51f97SToby Isaac   Not collective
4419be51f97SToby Isaac 
4429be51f97SToby Isaac   Input Parameter:
4439be51f97SToby Isaac . dm - the forest
4449be51f97SToby Isaac 
4459be51f97SToby Isaac   Output Parameter:
4469be51f97SToby Isaac . adapt - the forest from which dm is/was adapted
4479be51f97SToby Isaac 
4489be51f97SToby Isaac   Level: intermediate
4499be51f97SToby Isaac 
4509be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose()
4519be51f97SToby Isaac @*/
452ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt)
453dd8e54a2SToby Isaac {
454ba936b91SToby Isaac   DM_Forest      *forest;
45526d9498aSToby Isaac   PetscErrorCode ierr;
456dd8e54a2SToby Isaac 
457dd8e54a2SToby Isaac   PetscFunctionBegin;
458dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
459ba936b91SToby Isaac   forest = (DM_Forest*) dm->data;
46026d9498aSToby Isaac   switch (forest->adaptPurpose) {
461cd3c525cSToby Isaac   case DM_ADAPT_DETERMINE:
462ba936b91SToby Isaac     *adapt = forest->adapt;
46326d9498aSToby Isaac     break;
464a1b0c543SToby Isaac   case DM_ADAPT_REFINE:
46526d9498aSToby Isaac     ierr = DMGetCoarseDM(dm,adapt);CHKERRQ(ierr);
46626d9498aSToby Isaac     break;
467a1b0c543SToby Isaac   case DM_ADAPT_COARSEN:
468*bf2d5fbbSStefano Zampini   case DM_ADAPT_COARSEN_LAST:
46926d9498aSToby Isaac     ierr = DMGetFineDM(dm,adapt);CHKERRQ(ierr);
47026d9498aSToby Isaac     break;
47126d9498aSToby Isaac   default:
47226d9498aSToby Isaac     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
47326d9498aSToby Isaac   }
47426d9498aSToby Isaac   PetscFunctionReturn(0);
47526d9498aSToby Isaac }
47626d9498aSToby Isaac 
4779be51f97SToby Isaac /*@
4789be51f97SToby Isaac   DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its
479a1b0c543SToby Isaac   source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening
480cd3c525cSToby Isaac   (DM_ADAPT_COARSEN), or undefined (DM_ADAPT_DETERMINE).  This only matters for the purposes of reference counting:
4819be51f97SToby Isaac   during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse
4829be51f97SToby Isaac   relationship (see DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does
4839be51f97SToby Isaac   not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can
4849be51f97SToby Isaac   cause a memory leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
4859be51f97SToby Isaac 
4869be51f97SToby Isaac   Logically collective on dm
4879be51f97SToby Isaac 
4889be51f97SToby Isaac   Input Parameters:
4899be51f97SToby Isaac + dm - the forest
490*bf2d5fbbSStefano Zampini - purpose - the adaptivity purpose
4919be51f97SToby Isaac 
4929be51f97SToby Isaac   Level: advanced
4939be51f97SToby Isaac 
494*bf2d5fbbSStefano Zampini .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest(), DMAdaptFlag
4959be51f97SToby Isaac @*/
496a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMAdaptFlag purpose)
49726d9498aSToby Isaac {
49826d9498aSToby Isaac   DM_Forest      *forest;
49926d9498aSToby Isaac   PetscErrorCode ierr;
50026d9498aSToby Isaac 
50126d9498aSToby Isaac   PetscFunctionBegin;
50226d9498aSToby Isaac   forest = (DM_Forest*) dm->data;
5039be51f97SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
50426d9498aSToby Isaac   if (purpose != forest->adaptPurpose) {
50526d9498aSToby Isaac     DM adapt;
50626d9498aSToby Isaac 
50726d9498aSToby Isaac     ierr = DMForestGetAdaptivityForest(dm,&adapt);CHKERRQ(ierr);
50826d9498aSToby Isaac     ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr);
50926d9498aSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
510f885a11aSToby Isaac 
51126d9498aSToby Isaac     forest->adaptPurpose = purpose;
512f885a11aSToby Isaac 
51326d9498aSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,adapt);CHKERRQ(ierr);
51426d9498aSToby Isaac     ierr = DMDestroy(&adapt);CHKERRQ(ierr);
51526d9498aSToby Isaac   }
516dd8e54a2SToby Isaac   PetscFunctionReturn(0);
517dd8e54a2SToby Isaac }
518dd8e54a2SToby Isaac 
51956c0450aSToby Isaac /*@
52056c0450aSToby Isaac   DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with
521*bf2d5fbbSStefano Zampini   DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening (DM_ADAPT_COARSEN),
522*bf2d5fbbSStefano Zampini   coarsening only the last level (DM_ADAPT_COARSEN_LAST) or undefined (DM_ADAPT_DETERMINE).
523*bf2d5fbbSStefano Zampini   This only matters for the purposes of reference counting: during DMDestroy(), cyclic
52456c0450aSToby Isaac   references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see
52556c0450aSToby Isaac   DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does not maintain a
52656c0450aSToby Isaac   reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory
52756c0450aSToby Isaac   leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
52856c0450aSToby Isaac 
52956c0450aSToby Isaac   Not collective
53056c0450aSToby Isaac 
53156c0450aSToby Isaac   Input Parameter:
53256c0450aSToby Isaac . dm - the forest
53356c0450aSToby Isaac 
53456c0450aSToby Isaac   Output Parameter:
535*bf2d5fbbSStefano Zampini . purpose - the adaptivity purpose
53656c0450aSToby Isaac 
53756c0450aSToby Isaac   Level: advanced
53856c0450aSToby Isaac 
539*bf2d5fbbSStefano Zampini .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest(), DMAdaptFlag
54056c0450aSToby Isaac @*/
541a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMAdaptFlag *purpose)
54256c0450aSToby Isaac {
54356c0450aSToby Isaac   DM_Forest *forest;
54456c0450aSToby Isaac 
54556c0450aSToby Isaac   PetscFunctionBegin;
54656c0450aSToby Isaac   forest   = (DM_Forest*) dm->data;
54756c0450aSToby Isaac   *purpose = forest->adaptPurpose;
54856c0450aSToby Isaac   PetscFunctionReturn(0);
54956c0450aSToby Isaac }
55056c0450aSToby Isaac 
5519be51f97SToby Isaac /*@
5529be51f97SToby Isaac   DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine
5539be51f97SToby Isaac   cell adjacency (for the purposes of partitioning and overlap).
5549be51f97SToby Isaac 
5559be51f97SToby Isaac   Logically collective on dm
5569be51f97SToby Isaac 
5579be51f97SToby Isaac   Input Parameters:
5589be51f97SToby Isaac + dm - the forest
5599be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency)
5609be51f97SToby Isaac 
5619be51f97SToby Isaac   Level: intermediate
5629be51f97SToby Isaac 
5639be51f97SToby Isaac .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap()
5649be51f97SToby Isaac @*/
565dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim)
566dd8e54a2SToby Isaac {
567dd8e54a2SToby Isaac   PetscInt       dim;
568dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
569dd8e54a2SToby Isaac   PetscErrorCode ierr;
570dd8e54a2SToby Isaac 
571dd8e54a2SToby Isaac   PetscFunctionBegin;
572dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
573ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup");
574dd8e54a2SToby Isaac   if (adjDim < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %d", adjDim);
575dd8e54a2SToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
576dd8e54a2SToby Isaac   if (adjDim > dim) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %d: %d", dim, adjDim);
577dd8e54a2SToby Isaac   forest->adjDim = adjDim;
578dd8e54a2SToby Isaac   PetscFunctionReturn(0);
579dd8e54a2SToby Isaac }
580dd8e54a2SToby Isaac 
5819be51f97SToby Isaac /*@
5829be51f97SToby Isaac   DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that,
5839be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
5849be51f97SToby Isaac 
5859be51f97SToby Isaac   Logically collective on dm
5869be51f97SToby Isaac 
5879be51f97SToby Isaac   Input Parameters:
5889be51f97SToby Isaac + dm - the forest
5899be51f97SToby Isaac - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
5909be51f97SToby Isaac 
5919be51f97SToby Isaac   Level: intermediate
5929be51f97SToby Isaac 
5939be51f97SToby Isaac .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension()
5949be51f97SToby Isaac @*/
595dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim)
596dd8e54a2SToby Isaac {
597dd8e54a2SToby Isaac   PetscInt       dim;
598dd8e54a2SToby Isaac   PetscErrorCode ierr;
599dd8e54a2SToby Isaac 
600dd8e54a2SToby Isaac   PetscFunctionBegin;
601dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
602dd8e54a2SToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
603dd8e54a2SToby Isaac   ierr = DMForestSetAdjacencyDimension(dm,dim-adjCodim);CHKERRQ(ierr);
604dd8e54a2SToby Isaac   PetscFunctionReturn(0);
605dd8e54a2SToby Isaac }
606dd8e54a2SToby Isaac 
6079be51f97SToby Isaac /*@
6089be51f97SToby Isaac   DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the
6099be51f97SToby Isaac   purposes of partitioning and overlap).
6109be51f97SToby Isaac 
6119be51f97SToby Isaac   Not collective
6129be51f97SToby Isaac 
6139be51f97SToby Isaac   Input Parameter:
6149be51f97SToby Isaac . dm - the forest
6159be51f97SToby Isaac 
6169be51f97SToby Isaac   Output Parameter:
6179be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency)
6189be51f97SToby Isaac 
6199be51f97SToby Isaac   Level: intermediate
6209be51f97SToby Isaac 
6219be51f97SToby Isaac .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap()
6229be51f97SToby Isaac @*/
623dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim)
624dd8e54a2SToby Isaac {
625dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
626dd8e54a2SToby Isaac 
627dd8e54a2SToby Isaac   PetscFunctionBegin;
628dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
629dd8e54a2SToby Isaac   PetscValidIntPointer(adjDim,2);
630dd8e54a2SToby Isaac   *adjDim = forest->adjDim;
631dd8e54a2SToby Isaac   PetscFunctionReturn(0);
632dd8e54a2SToby Isaac }
633dd8e54a2SToby Isaac 
6349be51f97SToby Isaac /*@
6359be51f97SToby Isaac   DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that,
6369be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
6379be51f97SToby Isaac 
6389be51f97SToby Isaac   Not collective
6399be51f97SToby Isaac 
6409be51f97SToby Isaac   Input Parameter:
6419be51f97SToby Isaac . dm - the forest
6429be51f97SToby Isaac 
6439be51f97SToby Isaac   Output Parameter:
6449be51f97SToby Isaac . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
6459be51f97SToby Isaac 
6469be51f97SToby Isaac   Level: intermediate
6479be51f97SToby Isaac 
6489be51f97SToby Isaac .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension()
6499be51f97SToby Isaac @*/
650dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim)
651dd8e54a2SToby Isaac {
652dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
653dd8e54a2SToby Isaac   PetscInt       dim;
654dd8e54a2SToby Isaac   PetscErrorCode ierr;
655dd8e54a2SToby Isaac 
656dd8e54a2SToby Isaac   PetscFunctionBegin;
657dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
658dd8e54a2SToby Isaac   PetscValidIntPointer(adjCodim,2);
659dd8e54a2SToby Isaac   ierr      = DMGetDimension(dm,&dim);CHKERRQ(ierr);
660dd8e54a2SToby Isaac   *adjCodim = dim - forest->adjDim;
661dd8e54a2SToby Isaac   PetscFunctionReturn(0);
662dd8e54a2SToby Isaac }
663dd8e54a2SToby Isaac 
6649be51f97SToby Isaac /*@
6659be51f97SToby Isaac   DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel
6669be51f97SToby Isaac   partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding
6679be51f97SToby Isaac   adjacent cells
6689be51f97SToby Isaac 
6699be51f97SToby Isaac   Logically collective on dm
6709be51f97SToby Isaac 
6719be51f97SToby Isaac   Input Parameters:
6729be51f97SToby Isaac + dm - the forest
6739be51f97SToby Isaac - overlap - default 0
6749be51f97SToby Isaac 
6759be51f97SToby Isaac   Level: intermediate
6769be51f97SToby Isaac 
6779be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
6789be51f97SToby Isaac @*/
679dd8e54a2SToby Isaac PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap)
680dd8e54a2SToby Isaac {
681dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
682dd8e54a2SToby Isaac 
683dd8e54a2SToby Isaac   PetscFunctionBegin;
684dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
685ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup");
686dd8e54a2SToby Isaac   if (overlap < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %d", overlap);
687dd8e54a2SToby Isaac   forest->overlap = overlap;
688dd8e54a2SToby Isaac   PetscFunctionReturn(0);
689dd8e54a2SToby Isaac }
690dd8e54a2SToby Isaac 
6919be51f97SToby Isaac /*@
6929be51f97SToby Isaac   DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values
6939be51f97SToby Isaac   > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells
6949be51f97SToby Isaac 
6959be51f97SToby Isaac   Not collective
6969be51f97SToby Isaac 
6979be51f97SToby Isaac   Input Parameter:
6989be51f97SToby Isaac . dm - the forest
6999be51f97SToby Isaac 
7009be51f97SToby Isaac   Output Parameter:
7019be51f97SToby Isaac . overlap - default 0
7029be51f97SToby Isaac 
7039be51f97SToby Isaac   Level: intermediate
7049be51f97SToby Isaac 
7059be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
7069be51f97SToby Isaac @*/
707dd8e54a2SToby Isaac PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap)
708dd8e54a2SToby Isaac {
709dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
710dd8e54a2SToby Isaac 
711dd8e54a2SToby Isaac   PetscFunctionBegin;
712dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
713dd8e54a2SToby Isaac   PetscValidIntPointer(overlap,2);
714dd8e54a2SToby Isaac   *overlap = forest->overlap;
715dd8e54a2SToby Isaac   PetscFunctionReturn(0);
716dd8e54a2SToby Isaac }
717dd8e54a2SToby Isaac 
7189be51f97SToby Isaac /*@
7199be51f97SToby Isaac   DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base
7209be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest
7219be51f97SToby Isaac   (see DMForestGetAdaptivityForest()) this limits the amount of coarsening.
7229be51f97SToby Isaac 
7239be51f97SToby Isaac   Logically collective on dm
7249be51f97SToby Isaac 
7259be51f97SToby Isaac   Input Parameters:
7269be51f97SToby Isaac + dm - the forest
7279be51f97SToby Isaac - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7289be51f97SToby Isaac 
7299be51f97SToby Isaac   Level: intermediate
7309be51f97SToby Isaac 
7319be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
7329be51f97SToby Isaac @*/
733dd8e54a2SToby Isaac PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement)
734dd8e54a2SToby Isaac {
735dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
736dd8e54a2SToby Isaac 
737dd8e54a2SToby Isaac   PetscFunctionBegin;
738dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
739ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup");
740dd8e54a2SToby Isaac   forest->minRefinement = minRefinement;
741dd8e54a2SToby Isaac   PetscFunctionReturn(0);
742dd8e54a2SToby Isaac }
743dd8e54a2SToby Isaac 
7449be51f97SToby Isaac /*@
7459be51f97SToby Isaac   DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see
7469be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest (see
7479be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of coarsening.
7489be51f97SToby Isaac 
7499be51f97SToby Isaac   Not collective
7509be51f97SToby Isaac 
7519be51f97SToby Isaac   Input Parameter:
7529be51f97SToby Isaac . dm - the forest
7539be51f97SToby Isaac 
7549be51f97SToby Isaac   Output Parameter:
7559be51f97SToby Isaac . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7569be51f97SToby Isaac 
7579be51f97SToby Isaac   Level: intermediate
7589be51f97SToby Isaac 
7599be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestGetMaximumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
7609be51f97SToby Isaac @*/
761dd8e54a2SToby Isaac PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement)
762dd8e54a2SToby Isaac {
763dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
764dd8e54a2SToby Isaac 
765dd8e54a2SToby Isaac   PetscFunctionBegin;
766dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
767dd8e54a2SToby Isaac   PetscValidIntPointer(minRefinement,2);
768dd8e54a2SToby Isaac   *minRefinement = forest->minRefinement;
769dd8e54a2SToby Isaac   PetscFunctionReturn(0);
770dd8e54a2SToby Isaac }
771dd8e54a2SToby Isaac 
7729be51f97SToby Isaac /*@
7739be51f97SToby Isaac   DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base
7749be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.
7759be51f97SToby Isaac 
7769be51f97SToby Isaac   Logically collective on dm
7779be51f97SToby Isaac 
7789be51f97SToby Isaac   Input Parameters:
7799be51f97SToby Isaac + dm - the forest
7809be51f97SToby Isaac - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
7819be51f97SToby Isaac 
7829be51f97SToby Isaac   Level: intermediate
7839be51f97SToby Isaac 
7849be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
7859be51f97SToby Isaac @*/
78656ba9f64SToby Isaac PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement)
78756ba9f64SToby Isaac {
78856ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
78956ba9f64SToby Isaac 
79056ba9f64SToby Isaac   PetscFunctionBegin;
79156ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
79256ba9f64SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup");
79356ba9f64SToby Isaac   forest->initRefinement = initRefinement;
79456ba9f64SToby Isaac   PetscFunctionReturn(0);
79556ba9f64SToby Isaac }
79656ba9f64SToby Isaac 
7979be51f97SToby Isaac /*@
7989be51f97SToby Isaac   DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see
7999be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.
8009be51f97SToby Isaac 
8019be51f97SToby Isaac   Not collective
8029be51f97SToby Isaac 
8039be51f97SToby Isaac   Input Parameter:
8049be51f97SToby Isaac . dm - the forest
8059be51f97SToby Isaac 
8069be51f97SToby Isaac   Output Paramater:
807*bf2d5fbbSStefano Zampini . initRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8089be51f97SToby Isaac 
8099be51f97SToby Isaac   Level: intermediate
8109be51f97SToby Isaac 
8119be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
8129be51f97SToby Isaac @*/
81356ba9f64SToby Isaac PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement)
81456ba9f64SToby Isaac {
81556ba9f64SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
81656ba9f64SToby Isaac 
81756ba9f64SToby Isaac   PetscFunctionBegin;
81856ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
81956ba9f64SToby Isaac   PetscValidIntPointer(initRefinement,2);
82056ba9f64SToby Isaac   *initRefinement = forest->initRefinement;
82156ba9f64SToby Isaac   PetscFunctionReturn(0);
82256ba9f64SToby Isaac }
82356ba9f64SToby Isaac 
8249be51f97SToby Isaac /*@
8259be51f97SToby Isaac   DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base
8269be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest
8279be51f97SToby Isaac   (see DMForestGetAdaptivityForest()), this limits the amount of refinement.
8289be51f97SToby Isaac 
8299be51f97SToby Isaac   Logically collective on dm
8309be51f97SToby Isaac 
8319be51f97SToby Isaac   Input Parameters:
8329be51f97SToby Isaac + dm - the forest
8339be51f97SToby Isaac - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8349be51f97SToby Isaac 
8359be51f97SToby Isaac   Level: intermediate
8369be51f97SToby Isaac 
8379be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityDM()
8389be51f97SToby Isaac @*/
839c7eeac06SToby Isaac PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement)
840dd8e54a2SToby Isaac {
841dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
842dd8e54a2SToby Isaac 
843dd8e54a2SToby Isaac   PetscFunctionBegin;
844dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
845ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup");
846c7eeac06SToby Isaac   forest->maxRefinement = maxRefinement;
847dd8e54a2SToby Isaac   PetscFunctionReturn(0);
848dd8e54a2SToby Isaac }
849dd8e54a2SToby Isaac 
8509be51f97SToby Isaac /*@
8519be51f97SToby Isaac   DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see
8529be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest (see
8539be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of refinement.
8549be51f97SToby Isaac 
8559be51f97SToby Isaac   Not collective
8569be51f97SToby Isaac 
8579be51f97SToby Isaac   Input Parameter:
8589be51f97SToby Isaac . dm - the forest
8599be51f97SToby Isaac 
8609be51f97SToby Isaac   Output Parameter:
8619be51f97SToby Isaac . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
8629be51f97SToby Isaac 
8639be51f97SToby Isaac   Level: intermediate
8649be51f97SToby Isaac 
8659be51f97SToby Isaac .seealso: DMForestSetMaximumRefinement(), DMForestGetMinimumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
8669be51f97SToby Isaac @*/
867c7eeac06SToby Isaac PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement)
868dd8e54a2SToby Isaac {
869dd8e54a2SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
870dd8e54a2SToby Isaac 
871dd8e54a2SToby Isaac   PetscFunctionBegin;
872dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
873c7eeac06SToby Isaac   PetscValidIntPointer(maxRefinement,2);
874c7eeac06SToby Isaac   *maxRefinement = forest->maxRefinement;
875dd8e54a2SToby Isaac   PetscFunctionReturn(0);
876dd8e54a2SToby Isaac }
877c7eeac06SToby Isaac 
8789be51f97SToby Isaac /*@C
8799be51f97SToby Isaac   DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes.
8809be51f97SToby Isaac   Subtypes of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree
8819be51f97SToby Isaac   for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to
8829be51f97SToby Isaac   specify refinement/coarsening.
8839be51f97SToby Isaac 
8849be51f97SToby Isaac   Logically collective on dm
8859be51f97SToby Isaac 
8869be51f97SToby Isaac   Input Parameters:
8879be51f97SToby Isaac + dm - the forest
8889be51f97SToby Isaac - adaptStrategy - default DMFORESTADAPTALL
8899be51f97SToby Isaac 
8909be51f97SToby Isaac   Level: advanced
8919be51f97SToby Isaac 
8929be51f97SToby Isaac .seealso: DMForestGetAdaptivityStrategy()
8939be51f97SToby Isaac @*/
894c7eeac06SToby Isaac PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy)
895c7eeac06SToby Isaac {
896c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
897c7eeac06SToby Isaac   PetscErrorCode ierr;
898c7eeac06SToby Isaac 
899c7eeac06SToby Isaac   PetscFunctionBegin;
900c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
901c7eeac06SToby Isaac   ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr);
902a73e2921SToby Isaac   ierr = PetscStrallocpy((const char*) adaptStrategy,(char**)&forest->adaptStrategy);CHKERRQ(ierr);
903c7eeac06SToby Isaac   PetscFunctionReturn(0);
904c7eeac06SToby Isaac }
905c7eeac06SToby Isaac 
9069be51f97SToby Isaac /*@C
9079be51f97SToby Isaac   DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes.  Subtypes
9089be51f97SToby Isaac   of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all
9099be51f97SToby Isaac   processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only
9109be51f97SToby Isaac   one process needs to specify refinement/coarsening.
9119be51f97SToby Isaac 
9129be51f97SToby Isaac   Not collective
9139be51f97SToby Isaac 
9149be51f97SToby Isaac   Input Parameter:
9159be51f97SToby Isaac . dm - the forest
9169be51f97SToby Isaac 
9179be51f97SToby Isaac   Output Parameter:
9189be51f97SToby Isaac . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL)
9199be51f97SToby Isaac 
9209be51f97SToby Isaac   Level: advanced
9219be51f97SToby Isaac 
9229be51f97SToby Isaac .seealso: DMForestSetAdaptivityStrategy()
9239be51f97SToby Isaac @*/
924c7eeac06SToby Isaac PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy)
925c7eeac06SToby Isaac {
926c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
927c7eeac06SToby Isaac 
928c7eeac06SToby Isaac   PetscFunctionBegin;
929c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
930c7eeac06SToby Isaac   PetscValidPointer(adaptStrategy,2);
931c7eeac06SToby Isaac   *adaptStrategy = forest->adaptStrategy;
932c7eeac06SToby Isaac   PetscFunctionReturn(0);
933c7eeac06SToby Isaac }
934c7eeac06SToby Isaac 
9352a133e43SToby Isaac /*@
9362a133e43SToby Isaac   DMForestGetAdaptivitySuccess - Return whether the requested adaptation (refinement, coarsening, repartitioning,
9372a133e43SToby Isaac   etc.) was successful.  PETSC_FALSE indicates that the post-adaptation forest is the same as the pre-adpatation
9382a133e43SToby Isaac   forest.  A requested adaptation may have been unsuccessful if, for example, the requested refinement would have
9392a133e43SToby Isaac   exceeded the maximum refinement level.
9402a133e43SToby Isaac 
9412a133e43SToby Isaac   Collective on dm
9422a133e43SToby Isaac 
9432a133e43SToby Isaac   Input Parameter:
9442a133e43SToby Isaac 
9452a133e43SToby Isaac . dm - the post-adaptation forest
9462a133e43SToby Isaac 
9472a133e43SToby Isaac   Output Parameter:
9482a133e43SToby Isaac 
9492a133e43SToby Isaac . success - PETSC_TRUE if the post-adaptation forest is different from the pre-adaptation forest.
9502a133e43SToby Isaac 
9512a133e43SToby Isaac   Level: intermediate
9522a133e43SToby Isaac 
9532a133e43SToby Isaac .see
9542a133e43SToby Isaac @*/
9552a133e43SToby Isaac PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success)
9562a133e43SToby Isaac {
9572a133e43SToby Isaac   DM_Forest      *forest;
9582a133e43SToby Isaac   PetscErrorCode ierr;
9592a133e43SToby Isaac 
9602a133e43SToby Isaac   PetscFunctionBegin;
9612a133e43SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
9622a133e43SToby Isaac   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() has not been called yet.");
9632a133e43SToby Isaac   forest = (DM_Forest *) dm->data;
9642a133e43SToby Isaac   ierr = (forest->getadaptivitysuccess)(dm,success);CHKERRQ(ierr);
9652a133e43SToby Isaac   PetscFunctionReturn(0);
9662a133e43SToby Isaac }
9672a133e43SToby Isaac 
968bf9b5d84SToby Isaac /*@
969bf9b5d84SToby Isaac   DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed
970bf9b5d84SToby 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().
971bf9b5d84SToby Isaac 
972bf9b5d84SToby Isaac   Logically collective on dm
973bf9b5d84SToby Isaac 
974bf9b5d84SToby Isaac   Input Parameters:
975bf9b5d84SToby Isaac + dm - the post-adaptation forest
976bf9b5d84SToby Isaac - computeSF - default PETSC_TRUE
977bf9b5d84SToby Isaac 
978bf9b5d84SToby Isaac   Level: advanced
979bf9b5d84SToby Isaac 
980bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
981bf9b5d84SToby Isaac @*/
982bf9b5d84SToby Isaac PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF)
983bf9b5d84SToby Isaac {
984bf9b5d84SToby Isaac   DM_Forest *forest;
985bf9b5d84SToby Isaac 
986bf9b5d84SToby Isaac   PetscFunctionBegin;
987bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
988bf9b5d84SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called");
989bf9b5d84SToby Isaac   forest                 = (DM_Forest*) dm->data;
990bf9b5d84SToby Isaac   forest->computeAdaptSF = computeSF;
991bf9b5d84SToby Isaac   PetscFunctionReturn(0);
992bf9b5d84SToby Isaac }
993bf9b5d84SToby Isaac 
9940eb7e1eaSToby Isaac PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time)
99580b27e07SToby Isaac {
99680b27e07SToby Isaac   DM_Forest      *forest;
99780b27e07SToby Isaac   PetscErrorCode ierr;
99880b27e07SToby Isaac 
99980b27e07SToby Isaac   PetscFunctionBegin;
100080b27e07SToby Isaac   PetscValidHeaderSpecific(dmIn   ,DM_CLASSID  ,1);
100180b27e07SToby Isaac   PetscValidHeaderSpecific(vecIn  ,VEC_CLASSID ,2);
100280b27e07SToby Isaac   PetscValidHeaderSpecific(dmOut  ,DM_CLASSID  ,3);
100380b27e07SToby Isaac   PetscValidHeaderSpecific(vecOut ,VEC_CLASSID ,4);
100480b27e07SToby Isaac   forest = (DM_Forest *) dmIn->data;
100580b27e07SToby Isaac   if (!forest->transfervec) SETERRQ(PetscObjectComm((PetscObject)dmIn),PETSC_ERR_SUP,"DMForestTransferVec() not implemented");
10060eb7e1eaSToby Isaac   ierr = (forest->transfervec)(dmIn,vecIn,dmOut,vecOut,useBCs,time);CHKERRQ(ierr);
100780b27e07SToby Isaac   PetscFunctionReturn(0);
100880b27e07SToby Isaac }
100980b27e07SToby Isaac 
1010bf9b5d84SToby Isaac /*@
1011bf9b5d84SToby Isaac   DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the
1012bf9b5d84SToby Isaac   pre-adaptation forest to the post-adaptiation forest.  After DMSetUp() is called, these transfer PetscSFs can be
1013bf9b5d84SToby Isaac   accessed with DMForestGetAdaptivitySF().
1014bf9b5d84SToby Isaac 
1015bf9b5d84SToby Isaac   Not collective
1016bf9b5d84SToby Isaac 
1017bf9b5d84SToby Isaac   Input Parameter:
1018bf9b5d84SToby Isaac . dm - the post-adaptation forest
1019bf9b5d84SToby Isaac 
1020bf9b5d84SToby Isaac   Output Parameter:
1021bf9b5d84SToby Isaac . computeSF - default PETSC_TRUE
1022bf9b5d84SToby Isaac 
1023bf9b5d84SToby Isaac   Level: advanced
1024bf9b5d84SToby Isaac 
1025bf9b5d84SToby Isaac .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
1026bf9b5d84SToby Isaac @*/
1027bf9b5d84SToby Isaac PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF)
1028bf9b5d84SToby Isaac {
1029bf9b5d84SToby Isaac   DM_Forest *forest;
1030bf9b5d84SToby Isaac 
1031bf9b5d84SToby Isaac   PetscFunctionBegin;
1032bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1033bf9b5d84SToby Isaac   forest     = (DM_Forest*) dm->data;
1034bf9b5d84SToby Isaac   *computeSF = forest->computeAdaptSF;
1035bf9b5d84SToby Isaac   PetscFunctionReturn(0);
1036bf9b5d84SToby Isaac }
1037bf9b5d84SToby Isaac 
1038bf9b5d84SToby Isaac /*@
1039bf9b5d84SToby Isaac   DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest.
1040bf9b5d84SToby Isaac   Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be
1041bf9b5d84SToby Isaac   some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa.  Therefore there are two
1042bf9b5d84SToby Isaac   PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates
1043bf9b5d84SToby Isaac   pre-adaptation fine cells to post-adaptation coarse cells.
1044bf9b5d84SToby Isaac 
1045bf9b5d84SToby Isaac   Not collective
1046bf9b5d84SToby Isaac 
1047bf9b5d84SToby Isaac   Input Parameter:
1048bf9b5d84SToby Isaac   dm - the post-adaptation forest
1049bf9b5d84SToby Isaac 
1050bf9b5d84SToby Isaac   Output Parameter:
10510f17b9e3SToby Isaac   preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post-
10520f17b9e3SToby Isaac   coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre-
1053bf9b5d84SToby Isaac 
1054bf9b5d84SToby Isaac   Level: advanced
1055bf9b5d84SToby Isaac 
1056bf9b5d84SToby Isaac .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF()
1057bf9b5d84SToby Isaac @*/
10580f17b9e3SToby Isaac PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine)
1059bf9b5d84SToby Isaac {
1060bf9b5d84SToby Isaac   DM_Forest      *forest;
1061bf9b5d84SToby Isaac   PetscErrorCode ierr;
1062bf9b5d84SToby Isaac 
1063bf9b5d84SToby Isaac   PetscFunctionBegin;
1064bf9b5d84SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1065bf9b5d84SToby Isaac   ierr   = DMSetUp(dm);CHKERRQ(ierr);
1066bf9b5d84SToby Isaac   forest = (DM_Forest*) dm->data;
1067f885a11aSToby Isaac   if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine;
1068f885a11aSToby Isaac   if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine;
1069bf9b5d84SToby Isaac   PetscFunctionReturn(0);
1070bf9b5d84SToby Isaac }
1071bf9b5d84SToby Isaac 
10729be51f97SToby Isaac /*@
10739be51f97SToby Isaac   DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to
10749be51f97SToby Isaac   indicate that the diameter of neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may
10759be51f97SToby Isaac   only support one particular choice of grading factor.
10769be51f97SToby Isaac 
10779be51f97SToby Isaac   Logically collective on dm
10789be51f97SToby Isaac 
10799be51f97SToby Isaac   Input Parameters:
10809be51f97SToby Isaac + dm - the forest
10819be51f97SToby Isaac - grade - the grading factor
10829be51f97SToby Isaac 
10839be51f97SToby Isaac   Level: advanced
10849be51f97SToby Isaac 
10859be51f97SToby Isaac .seealso: DMForestGetGradeFactor()
10869be51f97SToby Isaac @*/
1087c7eeac06SToby Isaac PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade)
1088c7eeac06SToby Isaac {
1089c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1090c7eeac06SToby Isaac 
1091c7eeac06SToby Isaac   PetscFunctionBegin;
1092c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1093ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup");
1094c7eeac06SToby Isaac   forest->gradeFactor = grade;
1095c7eeac06SToby Isaac   PetscFunctionReturn(0);
1096c7eeac06SToby Isaac }
1097c7eeac06SToby Isaac 
10989be51f97SToby Isaac /*@
10999be51f97SToby Isaac   DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of
11009be51f97SToby Isaac   neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may only support one particular
11019be51f97SToby Isaac   choice of grading factor.
11029be51f97SToby Isaac 
11039be51f97SToby Isaac   Not collective
11049be51f97SToby Isaac 
11059be51f97SToby Isaac   Input Parameter:
11069be51f97SToby Isaac . dm - the forest
11079be51f97SToby Isaac 
11089be51f97SToby Isaac   Output Parameter:
11099be51f97SToby Isaac . grade - the grading factor
11109be51f97SToby Isaac 
11119be51f97SToby Isaac   Level: advanced
11129be51f97SToby Isaac 
11139be51f97SToby Isaac .seealso: DMForestSetGradeFactor()
11149be51f97SToby Isaac @*/
1115c7eeac06SToby Isaac PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade)
1116c7eeac06SToby Isaac {
1117c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1118c7eeac06SToby Isaac 
1119c7eeac06SToby Isaac   PetscFunctionBegin;
1120c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1121c7eeac06SToby Isaac   PetscValidIntPointer(grade,2);
1122c7eeac06SToby Isaac   *grade = forest->gradeFactor;
1123c7eeac06SToby Isaac   PetscFunctionReturn(0);
1124c7eeac06SToby Isaac }
1125c7eeac06SToby Isaac 
11269be51f97SToby Isaac /*@
11279be51f97SToby Isaac   DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes
11289be51f97SToby Isaac   the cell weight (see DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be
11299be51f97SToby Isaac   (cellWeight) * (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on
11309be51f97SToby Isaac   its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the
11319be51f97SToby Isaac   computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
11329be51f97SToby Isaac 
11339be51f97SToby Isaac   Logically collective on dm
11349be51f97SToby Isaac 
11359be51f97SToby Isaac   Input Parameters:
11369be51f97SToby Isaac + dm - the forest
11379be51f97SToby Isaac - weightsFactors - default 1.
11389be51f97SToby Isaac 
11399be51f97SToby Isaac   Level: advanced
11409be51f97SToby Isaac 
11419be51f97SToby Isaac .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights()
11429be51f97SToby Isaac @*/
1143ef51cf95SToby Isaac PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor)
1144c7eeac06SToby Isaac {
1145c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1146c7eeac06SToby Isaac 
1147c7eeac06SToby Isaac   PetscFunctionBegin;
1148c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1149ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup");
1150c7eeac06SToby Isaac   forest->weightsFactor = weightsFactor;
1151c7eeac06SToby Isaac   PetscFunctionReturn(0);
1152c7eeac06SToby Isaac }
1153c7eeac06SToby Isaac 
11549be51f97SToby Isaac /*@
11559be51f97SToby Isaac   DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see
11569be51f97SToby Isaac   DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be (cellWeight) *
11579be51f97SToby Isaac   (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on its level; a
11589be51f97SToby Isaac   factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation
11599be51f97SToby Isaac   associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
11609be51f97SToby Isaac 
11619be51f97SToby Isaac   Not collective
11629be51f97SToby Isaac 
11639be51f97SToby Isaac   Input Parameter:
11649be51f97SToby Isaac . dm - the forest
11659be51f97SToby Isaac 
11669be51f97SToby Isaac   Output Parameter:
11679be51f97SToby Isaac . weightsFactors - default 1.
11689be51f97SToby Isaac 
11699be51f97SToby Isaac   Level: advanced
11709be51f97SToby Isaac 
11719be51f97SToby Isaac .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights()
11729be51f97SToby Isaac @*/
1173ef51cf95SToby Isaac PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor)
1174c7eeac06SToby Isaac {
1175c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1176c7eeac06SToby Isaac 
1177c7eeac06SToby Isaac   PetscFunctionBegin;
1178c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1179c7eeac06SToby Isaac   PetscValidRealPointer(weightsFactor,2);
1180c7eeac06SToby Isaac   *weightsFactor = forest->weightsFactor;
1181c7eeac06SToby Isaac   PetscFunctionReturn(0);
1182c7eeac06SToby Isaac }
1183c7eeac06SToby Isaac 
11849be51f97SToby Isaac /*@
11859be51f97SToby Isaac   DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process
11869be51f97SToby Isaac 
11879be51f97SToby Isaac   Not collective
11889be51f97SToby Isaac 
11899be51f97SToby Isaac   Input Parameter:
11909be51f97SToby Isaac . dm - the forest
11919be51f97SToby Isaac 
11929be51f97SToby Isaac   Output Parameters:
11939be51f97SToby Isaac + cStart - the first cell on this process
11949be51f97SToby Isaac - cEnd - one after the final cell on this process
11959be51f97SToby Isaac 
11961a244344SSatish Balay   Level: intermediate
11979be51f97SToby Isaac 
11989be51f97SToby Isaac .seealso: DMForestGetCellSF()
11999be51f97SToby Isaac @*/
1200c7eeac06SToby Isaac PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd)
1201c7eeac06SToby Isaac {
1202c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1203c7eeac06SToby Isaac   PetscErrorCode ierr;
1204c7eeac06SToby Isaac 
1205c7eeac06SToby Isaac   PetscFunctionBegin;
1206c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1207c7eeac06SToby Isaac   PetscValidIntPointer(cStart,2);
1208c7eeac06SToby Isaac   PetscValidIntPointer(cEnd,2);
1209c7eeac06SToby Isaac   if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) {
1210c7eeac06SToby Isaac     ierr = forest->createcellchart(dm,&forest->cStart,&forest->cEnd);CHKERRQ(ierr);
1211c7eeac06SToby Isaac   }
1212c7eeac06SToby Isaac   *cStart =  forest->cStart;
1213c7eeac06SToby Isaac   *cEnd   =  forest->cEnd;
1214c7eeac06SToby Isaac   PetscFunctionReturn(0);
1215c7eeac06SToby Isaac }
1216c7eeac06SToby Isaac 
12179be51f97SToby Isaac /*@
12189be51f97SToby Isaac   DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes
12199be51f97SToby Isaac 
12209be51f97SToby Isaac   Not collective
12219be51f97SToby Isaac 
12229be51f97SToby Isaac   Input Parameter:
12239be51f97SToby Isaac . dm - the forest
12249be51f97SToby Isaac 
12259be51f97SToby Isaac   Output Parameter:
12269be51f97SToby Isaac . cellSF - the PetscSF
12279be51f97SToby Isaac 
12281a244344SSatish Balay   Level: intermediate
12299be51f97SToby Isaac 
12309be51f97SToby Isaac .seealso: DMForestGetCellChart()
12319be51f97SToby Isaac @*/
1232c7eeac06SToby Isaac PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF)
1233c7eeac06SToby Isaac {
1234c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1235c7eeac06SToby Isaac   PetscErrorCode ierr;
1236c7eeac06SToby Isaac 
1237c7eeac06SToby Isaac   PetscFunctionBegin;
1238c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1239c7eeac06SToby Isaac   PetscValidPointer(cellSF,2);
1240c7eeac06SToby Isaac   if ((!forest->cellSF) && forest->createcellsf) {
1241c7eeac06SToby Isaac     ierr = forest->createcellsf(dm,&forest->cellSF);CHKERRQ(ierr);
1242c7eeac06SToby Isaac   }
1243c7eeac06SToby Isaac   *cellSF = forest->cellSF;
1244c7eeac06SToby Isaac   PetscFunctionReturn(0);
1245c7eeac06SToby Isaac }
1246c7eeac06SToby Isaac 
12479be51f97SToby Isaac /*@C
12489be51f97SToby Isaac   DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see
12499be51f97SToby Isaac   DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination).  The
1250cd3c525cSToby Isaac   interpretation of the label values is up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP,
1251cd3c525cSToby Isaac   DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have been reserved as choices that should be accepted by all subtypes.
12529be51f97SToby Isaac 
12539be51f97SToby Isaac   Logically collective on dm
12549be51f97SToby Isaac 
12559be51f97SToby Isaac   Input Parameters:
12569be51f97SToby Isaac - dm - the forest
1257a1b0c543SToby Isaac + adaptLabel - the label in the pre-adaptation forest
12589be51f97SToby Isaac 
12599be51f97SToby Isaac   Level: intermediate
12609be51f97SToby Isaac 
12619be51f97SToby Isaac .seealso DMForestGetAdaptivityLabel()
12629be51f97SToby Isaac @*/
1263a1b0c543SToby Isaac PetscErrorCode DMForestSetAdaptivityLabel(DM dm, DMLabel adaptLabel)
1264c7eeac06SToby Isaac {
1265c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1266c7eeac06SToby Isaac   PetscErrorCode ierr;
1267c7eeac06SToby Isaac 
1268c7eeac06SToby Isaac   PetscFunctionBegin;
1269c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1270a1b0c543SToby Isaac   adaptLabel->refct++;
1271a1b0c543SToby Isaac   if (forest->adaptLabel) {ierr = DMLabelDestroy(&forest->adaptLabel);CHKERRQ(ierr);}
1272a1b0c543SToby Isaac   forest->adaptLabel = adaptLabel;
1273c7eeac06SToby Isaac   PetscFunctionReturn(0);
1274c7eeac06SToby Isaac }
1275c7eeac06SToby Isaac 
12769be51f97SToby Isaac /*@C
12779be51f97SToby Isaac   DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that
12789be51f97SToby Isaac   holds the adaptation flags (refinement, coarsening, or some combination).  The interpretation of the label values is
1279cd3c525cSToby Isaac   up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have
1280cd3c525cSToby Isaac   been reserved as choices that should be accepted by all subtypes.
12819be51f97SToby Isaac 
12829be51f97SToby Isaac   Not collective
12839be51f97SToby Isaac 
12849be51f97SToby Isaac   Input Parameter:
12859be51f97SToby Isaac . dm - the forest
12869be51f97SToby Isaac 
12879be51f97SToby Isaac   Output Parameter:
12889be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest
12899be51f97SToby Isaac 
12909be51f97SToby Isaac   Level: intermediate
12919be51f97SToby Isaac 
12929be51f97SToby Isaac .seealso DMForestSetAdaptivityLabel()
12939be51f97SToby Isaac @*/
1294a1b0c543SToby Isaac PetscErrorCode DMForestGetAdaptivityLabel(DM dm, DMLabel *adaptLabel)
1295c7eeac06SToby Isaac {
1296c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1297c7eeac06SToby Isaac 
1298c7eeac06SToby Isaac   PetscFunctionBegin;
1299c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1300ba936b91SToby Isaac   *adaptLabel = forest->adaptLabel;
1301c7eeac06SToby Isaac   PetscFunctionReturn(0);
1302c7eeac06SToby Isaac }
1303c7eeac06SToby Isaac 
13049be51f97SToby Isaac /*@
13059be51f97SToby Isaac   DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
13069be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
13079be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
13089be51f97SToby Isaac   of 1.
13099be51f97SToby Isaac 
13109be51f97SToby Isaac   Logically collective on dm
13119be51f97SToby Isaac 
13129be51f97SToby Isaac   Input Parameters:
13139be51f97SToby Isaac + dm - the forest
13149be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
13159be51f97SToby Isaac - copyMode - how weights should reference weights
13169be51f97SToby Isaac 
13179be51f97SToby Isaac   Level: advanced
13189be51f97SToby Isaac 
13199be51f97SToby Isaac .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity()
13209be51f97SToby Isaac @*/
1321c7eeac06SToby Isaac PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode)
1322c7eeac06SToby Isaac {
1323c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest*) dm->data;
1324c7eeac06SToby Isaac   PetscInt       cStart, cEnd;
1325c7eeac06SToby Isaac   PetscErrorCode ierr;
1326c7eeac06SToby Isaac 
1327c7eeac06SToby Isaac   PetscFunctionBegin;
1328c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1329c7eeac06SToby Isaac   ierr = DMForestGetCellChart(dm,&cStart,&cEnd);CHKERRQ(ierr);
1330c7eeac06SToby Isaac   if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd);
1331c7eeac06SToby Isaac   if (copyMode == PETSC_COPY_VALUES) {
1332c7eeac06SToby Isaac     if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) {
1333c7eeac06SToby Isaac       ierr = PetscMalloc1(cEnd-cStart,&forest->cellWeights);CHKERRQ(ierr);
1334c7eeac06SToby Isaac     }
1335c7eeac06SToby Isaac     ierr                        = PetscMemcpy(forest->cellWeights,weights,(cEnd-cStart)*sizeof(*weights));CHKERRQ(ierr);
1336c7eeac06SToby Isaac     forest->cellWeightsCopyMode = PETSC_OWN_POINTER;
1337c7eeac06SToby Isaac     PetscFunctionReturn(0);
1338c7eeac06SToby Isaac   }
1339c7eeac06SToby Isaac   if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) {
1340c7eeac06SToby Isaac     ierr = PetscFree(forest->cellWeights);CHKERRQ(ierr);
1341c7eeac06SToby Isaac   }
1342c7eeac06SToby Isaac   forest->cellWeights         = weights;
1343c7eeac06SToby Isaac   forest->cellWeightsCopyMode = copyMode;
1344c7eeac06SToby Isaac   PetscFunctionReturn(0);
1345c7eeac06SToby Isaac }
1346c7eeac06SToby Isaac 
13479be51f97SToby Isaac /*@
13489be51f97SToby Isaac   DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
13499be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
13509be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
13519be51f97SToby Isaac   of 1.
13529be51f97SToby Isaac 
13539be51f97SToby Isaac   Not collective
13549be51f97SToby Isaac 
13559be51f97SToby Isaac   Input Parameter:
13569be51f97SToby Isaac . dm - the forest
13579be51f97SToby Isaac 
13589be51f97SToby Isaac   Output Parameter:
13599be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
13609be51f97SToby Isaac 
13619be51f97SToby Isaac   Level: advanced
13629be51f97SToby Isaac 
13639be51f97SToby Isaac .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity()
13649be51f97SToby Isaac @*/
1365c7eeac06SToby Isaac PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights)
1366c7eeac06SToby Isaac {
1367c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1368c7eeac06SToby Isaac 
1369c7eeac06SToby Isaac   PetscFunctionBegin;
1370c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1371c7eeac06SToby Isaac   PetscValidPointer(weights,2);
1372c7eeac06SToby Isaac   *weights = forest->cellWeights;
1373c7eeac06SToby Isaac   PetscFunctionReturn(0);
1374c7eeac06SToby Isaac }
1375c7eeac06SToby Isaac 
13769be51f97SToby Isaac /*@
13779be51f97SToby Isaac   DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning
13789be51f97SToby Isaac   a pre-adaptation forest (see DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each
13799be51f97SToby Isaac   process's cells to the process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that
13809be51f97SToby Isaac   the current process should not have any cells after repartitioning.
13819be51f97SToby Isaac 
13829be51f97SToby Isaac   Logically Collective on dm
13839be51f97SToby Isaac 
13849be51f97SToby Isaac   Input parameters:
13859be51f97SToby Isaac + dm - the forest
13869be51f97SToby Isaac - capacity - this process's capacity
13879be51f97SToby Isaac 
13889be51f97SToby Isaac   Level: advanced
13899be51f97SToby Isaac 
13909be51f97SToby Isaac .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
13919be51f97SToby Isaac @*/
1392c7eeac06SToby Isaac PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity)
1393c7eeac06SToby Isaac {
1394c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1395c7eeac06SToby Isaac 
1396c7eeac06SToby Isaac   PetscFunctionBegin;
1397c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1398ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup");
1399c7eeac06SToby Isaac   if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity);
1400c7eeac06SToby Isaac   forest->weightCapacity = capacity;
1401c7eeac06SToby Isaac   PetscFunctionReturn(0);
1402c7eeac06SToby Isaac }
1403c7eeac06SToby Isaac 
14049be51f97SToby Isaac /*@
14059be51f97SToby Isaac   DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see
14069be51f97SToby Isaac   DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each process's cells to the
14079be51f97SToby Isaac   process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that the current process
14089be51f97SToby Isaac   should not have any cells after repartitioning.
14099be51f97SToby Isaac 
14109be51f97SToby Isaac   Not collective
14119be51f97SToby Isaac 
14129be51f97SToby Isaac   Input parameter:
14139be51f97SToby Isaac . dm - the forest
14149be51f97SToby Isaac 
14159be51f97SToby Isaac   Output parameter:
14169be51f97SToby Isaac . capacity - this process's capacity
14179be51f97SToby Isaac 
14189be51f97SToby Isaac   Level: advanced
14199be51f97SToby Isaac 
14209be51f97SToby Isaac .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
14219be51f97SToby Isaac @*/
1422c7eeac06SToby Isaac PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity)
1423c7eeac06SToby Isaac {
1424c7eeac06SToby Isaac   DM_Forest *forest = (DM_Forest*) dm->data;
1425c7eeac06SToby Isaac 
1426c7eeac06SToby Isaac   PetscFunctionBegin;
1427c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1428c7eeac06SToby Isaac   PetscValidRealPointer(capacity,2);
1429c7eeac06SToby Isaac   *capacity = forest->weightCapacity;
1430c7eeac06SToby Isaac   PetscFunctionReturn(0);
1431c7eeac06SToby Isaac }
1432c7eeac06SToby Isaac 
14330709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm)
1434db4d5e8cSToby Isaac {
1435db4d5e8cSToby Isaac   DM_Forest                  *forest = (DM_Forest*) dm->data;
143656ba9f64SToby Isaac   PetscBool                  flg, flg1, flg2, flg3, flg4;
1437dd8e54a2SToby Isaac   DMForestTopology           oldTopo;
1438c7eeac06SToby Isaac   char                       stringBuffer[256];
1439dd8e54a2SToby Isaac   PetscViewer                viewer;
1440dd8e54a2SToby Isaac   PetscViewerFormat          format;
144156ba9f64SToby Isaac   PetscInt                   adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade;
1442c7eeac06SToby Isaac   PetscReal                  weightsFactor;
1443c7eeac06SToby Isaac   DMForestAdaptivityStrategy adaptStrategy;
1444db4d5e8cSToby Isaac   PetscErrorCode             ierr;
1445db4d5e8cSToby Isaac 
1446db4d5e8cSToby Isaac   PetscFunctionBegin;
1447db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
144858762b62SToby Isaac   forest->setfromoptionscalled = PETSC_TRUE;
1449dd8e54a2SToby Isaac   ierr                         = DMForestGetTopology(dm, &oldTopo);CHKERRQ(ierr);
1450a3eda1eaSToby Isaac   ierr                         = PetscOptionsHead(PetscOptionsObject,"DMForest Options");CHKERRQ(ierr);
145156ba9f64SToby Isaac   ierr                         = PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);CHKERRQ(ierr);
145256ba9f64SToby Isaac   ierr                         = PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);CHKERRQ(ierr);
145356ba9f64SToby Isaac   ierr                         = PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);CHKERRQ(ierr);
145456ba9f64SToby Isaac   ierr                         = PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);CHKERRQ(ierr);
1455f885a11aSToby 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}");
145656ba9f64SToby Isaac   if (flg1) {
145756ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,(DMForestTopology)stringBuffer);CHKERRQ(ierr);
145856ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
145920e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
146056ba9f64SToby Isaac   }
146156ba9f64SToby Isaac   if (flg2) {
1462dd8e54a2SToby Isaac     DM base;
1463dd8e54a2SToby Isaac 
1464dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&base);CHKERRQ(ierr);
1465dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1466dd8e54a2SToby Isaac     ierr = DMLoad(base,viewer);CHKERRQ(ierr);
1467dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1468dd8e54a2SToby Isaac     ierr = DMForestSetBaseDM(dm,base);CHKERRQ(ierr);
1469dd8e54a2SToby Isaac     ierr = DMDestroy(&base);CHKERRQ(ierr);
147056ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
147120e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
1472dd8e54a2SToby Isaac   }
147356ba9f64SToby Isaac   if (flg3) {
1474dd8e54a2SToby Isaac     DM coarse;
1475dd8e54a2SToby Isaac 
1476dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&coarse);CHKERRQ(ierr);
1477dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1478dd8e54a2SToby Isaac     ierr = DMLoad(coarse,viewer);CHKERRQ(ierr);
1479dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
148020e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,coarse);CHKERRQ(ierr);
1481dd8e54a2SToby Isaac     ierr = DMDestroy(&coarse);CHKERRQ(ierr);
148256ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
148356ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1484dd8e54a2SToby Isaac   }
148556ba9f64SToby Isaac   if (flg4) {
1486dd8e54a2SToby Isaac     DM fine;
1487dd8e54a2SToby Isaac 
1488dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&fine);CHKERRQ(ierr);
1489dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1490dd8e54a2SToby Isaac     ierr = DMLoad(fine,viewer);CHKERRQ(ierr);
1491dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
149220e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,fine);CHKERRQ(ierr);
1493dd8e54a2SToby Isaac     ierr = DMDestroy(&fine);CHKERRQ(ierr);
149456ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
149556ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1496dd8e54a2SToby Isaac   }
1497dd8e54a2SToby Isaac   ierr = DMForestGetAdjacencyDimension(dm,&adjDim);CHKERRQ(ierr);
1498dd8e54a2SToby Isaac   ierr = PetscOptionsInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg);CHKERRQ(ierr);
1499dd8e54a2SToby Isaac   if (flg) {
1500dd8e54a2SToby Isaac     ierr = DMForestSetAdjacencyDimension(dm,adjDim);CHKERRQ(ierr);
1501f885a11aSToby Isaac   } else {
1502dd8e54a2SToby Isaac     ierr = DMForestGetAdjacencyCodimension(dm,&adjCodim);CHKERRQ(ierr);
1503dd8e54a2SToby Isaac     ierr = PetscOptionsInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg);CHKERRQ(ierr);
1504dd8e54a2SToby Isaac     if (flg) {
1505dd8e54a2SToby Isaac       ierr = DMForestSetAdjacencyCodimension(dm,adjCodim);CHKERRQ(ierr);
1506dd8e54a2SToby Isaac     }
1507dd8e54a2SToby Isaac   }
1508dd8e54a2SToby Isaac   ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr);
1509dd8e54a2SToby Isaac   ierr = PetscOptionsInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg);CHKERRQ(ierr);
1510dd8e54a2SToby Isaac   if (flg) {
1511dd8e54a2SToby Isaac     ierr = DMForestSetPartitionOverlap(dm,overlap);CHKERRQ(ierr);
1512dd8e54a2SToby Isaac   }
1513a6121fbdSMatthew G. Knepley #if 0
1514a6121fbdSMatthew 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);
1515a6121fbdSMatthew G. Knepley   if (flg) {
1516a6121fbdSMatthew G. Knepley     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1517a6121fbdSMatthew G. Knepley     ierr = DMForestSetInitialRefinement(dm,minRefinement);CHKERRQ(ierr);
1518a6121fbdSMatthew G. Knepley   }
1519a6121fbdSMatthew 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);
1520a6121fbdSMatthew G. Knepley   if (flg) {
1521a6121fbdSMatthew G. Knepley     ierr = DMForestSetMinimumRefinement(dm,0);CHKERRQ(ierr);
1522a6121fbdSMatthew G. Knepley     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
1523a6121fbdSMatthew G. Knepley   }
1524a6121fbdSMatthew G. Knepley #endif
1525dd8e54a2SToby Isaac   ierr = DMForestGetMinimumRefinement(dm,&minRefinement);CHKERRQ(ierr);
1526dd8e54a2SToby Isaac   ierr = PetscOptionsInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg);CHKERRQ(ierr);
1527dd8e54a2SToby Isaac   if (flg) {
1528dd8e54a2SToby Isaac     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1529db4d5e8cSToby Isaac   }
153056ba9f64SToby Isaac   ierr = DMForestGetInitialRefinement(dm,&initRefinement);CHKERRQ(ierr);
153156ba9f64SToby Isaac   ierr = PetscOptionsInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg);CHKERRQ(ierr);
153256ba9f64SToby Isaac   if (flg) {
153356ba9f64SToby Isaac     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
153456ba9f64SToby Isaac   }
1535c7eeac06SToby Isaac   ierr = DMForestGetMaximumRefinement(dm,&maxRefinement);CHKERRQ(ierr);
1536c7eeac06SToby Isaac   ierr = PetscOptionsInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg);CHKERRQ(ierr);
1537c7eeac06SToby Isaac   if (flg) {
1538c7eeac06SToby Isaac     ierr = DMForestSetMaximumRefinement(dm,maxRefinement);CHKERRQ(ierr);
1539c7eeac06SToby Isaac   }
1540c7eeac06SToby Isaac   ierr = DMForestGetAdaptivityStrategy(dm,&adaptStrategy);CHKERRQ(ierr);
1541c7eeac06SToby Isaac   ierr = PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,256,&flg);CHKERRQ(ierr);
1542c7eeac06SToby Isaac   if (flg) {
1543c7eeac06SToby Isaac     ierr = DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);CHKERRQ(ierr);
1544c7eeac06SToby Isaac   }
1545c7eeac06SToby Isaac   ierr = DMForestGetGradeFactor(dm,&grade);CHKERRQ(ierr);
1546c7eeac06SToby Isaac   ierr = PetscOptionsInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg);CHKERRQ(ierr);
1547c7eeac06SToby Isaac   if (flg) {
1548c7eeac06SToby Isaac     ierr = DMForestSetGradeFactor(dm,grade);CHKERRQ(ierr);
1549c7eeac06SToby Isaac   }
1550c7eeac06SToby Isaac   ierr = DMForestGetCellWeightFactor(dm,&weightsFactor);CHKERRQ(ierr);
1551c7eeac06SToby Isaac   ierr = PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);CHKERRQ(ierr);
1552c7eeac06SToby Isaac   if (flg) {
1553c7eeac06SToby Isaac     ierr = DMForestSetCellWeightFactor(dm,weightsFactor);CHKERRQ(ierr);
1554c7eeac06SToby Isaac   }
1555db4d5e8cSToby Isaac   ierr = PetscOptionsTail();CHKERRQ(ierr);
1556db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1557db4d5e8cSToby Isaac }
1558db4d5e8cSToby Isaac 
1559276c5506SMatthew G. Knepley PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, const PetscInt fields[], IS *is, DM *subdm)
1560d8984e3bSMatthew G. Knepley {
1561d8984e3bSMatthew G. Knepley   PetscErrorCode ierr;
1562d8984e3bSMatthew G. Knepley 
1563d8984e3bSMatthew G. Knepley   PetscFunctionBegin;
1564d8984e3bSMatthew G. Knepley   if (subdm) {ierr = DMClone(dm, subdm);CHKERRQ(ierr);}
1565d8984e3bSMatthew G. Knepley   ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr);
1566d8984e3bSMatthew G. Knepley   PetscFunctionReturn(0);
1567d8984e3bSMatthew G. Knepley }
1568d8984e3bSMatthew G. Knepley 
15695421bac9SToby Isaac PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
15705421bac9SToby Isaac {
15715421bac9SToby Isaac   DMLabel        refine;
15725421bac9SToby Isaac   DM             fineDM;
15735421bac9SToby Isaac   PetscErrorCode ierr;
15745421bac9SToby Isaac 
15755421bac9SToby Isaac   PetscFunctionBegin;
15765421bac9SToby Isaac   ierr = DMGetFineDM(dm,&fineDM);CHKERRQ(ierr);
15775421bac9SToby Isaac   if (fineDM) {
15785421bac9SToby Isaac     ierr       = PetscObjectReference((PetscObject)fineDM);CHKERRQ(ierr);
15795421bac9SToby Isaac     *dmRefined = fineDM;
15805421bac9SToby Isaac     PetscFunctionReturn(0);
15815421bac9SToby Isaac   }
15825421bac9SToby Isaac   ierr = DMForestTemplate(dm,comm,dmRefined);CHKERRQ(ierr);
15835421bac9SToby Isaac   ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr);
15845421bac9SToby Isaac   if (!refine) {
1585a1b0c543SToby Isaac     ierr = DMLabelCreate("refine",&refine);CHKERRQ(ierr);
1586a1b0c543SToby Isaac     ierr = DMLabelSetDefaultValue(refine,DM_ADAPT_REFINE);CHKERRQ(ierr);
15875421bac9SToby Isaac   }
1588a1b0c543SToby Isaac   else {
1589a1b0c543SToby Isaac     refine->refct++;
1590a1b0c543SToby Isaac   }
1591a1b0c543SToby Isaac   ierr = DMForestSetAdaptivityLabel(*dmRefined,refine);CHKERRQ(ierr);
1592a1b0c543SToby Isaac   ierr = DMLabelDestroy(&refine);CHKERRQ(ierr);
15935421bac9SToby Isaac   PetscFunctionReturn(0);
15945421bac9SToby Isaac }
15955421bac9SToby Isaac 
15965421bac9SToby Isaac PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
15975421bac9SToby Isaac {
15985421bac9SToby Isaac   DMLabel        coarsen;
15995421bac9SToby Isaac   DM             coarseDM;
16005421bac9SToby Isaac   PetscErrorCode ierr;
16015421bac9SToby Isaac 
16025421bac9SToby Isaac   PetscFunctionBegin;
16034098eed7SToby Isaac   {
16044098eed7SToby Isaac     PetscMPIInt mpiComparison;
16054098eed7SToby Isaac     MPI_Comm    dmcomm = PetscObjectComm((PetscObject)dm);
16064098eed7SToby Isaac 
16074098eed7SToby Isaac     ierr = MPI_Comm_compare(comm, dmcomm, &mpiComparison);CHKERRQ(ierr);
1608f885a11aSToby Isaac     if (mpiComparison != MPI_IDENT && mpiComparison != MPI_CONGRUENT) SETERRQ(dmcomm,PETSC_ERR_SUP,"No support for different communicators yet");
16094098eed7SToby Isaac   }
16105421bac9SToby Isaac   ierr = DMGetCoarseDM(dm,&coarseDM);CHKERRQ(ierr);
16115421bac9SToby Isaac   if (coarseDM) {
16125421bac9SToby Isaac     ierr         = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr);
16135421bac9SToby Isaac     *dmCoarsened = coarseDM;
16145421bac9SToby Isaac     PetscFunctionReturn(0);
16155421bac9SToby Isaac   }
16165421bac9SToby Isaac   ierr = DMForestTemplate(dm,comm,dmCoarsened);CHKERRQ(ierr);
16177ee90a76SStefano Zampini   ierr = DMForestSetAdaptivityPurpose(*dmCoarsened,DM_ADAPT_COARSEN);CHKERRQ(ierr);
16185421bac9SToby Isaac   ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr);
16195421bac9SToby Isaac   if (!coarsen) {
1620a1b0c543SToby Isaac     ierr = DMLabelCreate("coarsen",&coarsen);CHKERRQ(ierr);
1621a1b0c543SToby Isaac     ierr = DMLabelSetDefaultValue(coarsen,DM_ADAPT_COARSEN);CHKERRQ(ierr);
1622a1b0c543SToby Isaac   } else {
1623a1b0c543SToby Isaac     coarsen->refct++;
16245421bac9SToby Isaac   }
1625a1b0c543SToby Isaac   ierr = DMForestSetAdaptivityLabel(*dmCoarsened,coarsen);CHKERRQ(ierr);
1626a1b0c543SToby Isaac   ierr = DMLabelDestroy(&coarsen);CHKERRQ(ierr);
16275421bac9SToby Isaac   PetscFunctionReturn(0);
16285421bac9SToby Isaac }
16295421bac9SToby Isaac 
1630a1b0c543SToby Isaac static PetscErrorCode DMAdaptLabel_Forest(DM dm, DMLabel label, DM *adaptedDM)
163109350103SToby Isaac {
163209350103SToby Isaac   PetscBool      success;
163309350103SToby Isaac   PetscErrorCode ierr;
163409350103SToby Isaac 
163509350103SToby Isaac   PetscFunctionBegin;
163609350103SToby Isaac   ierr = DMForestTemplate(dm,PetscObjectComm((PetscObject)dm),adaptedDM);CHKERRQ(ierr);
1637a1b0c543SToby Isaac   ierr = DMForestSetAdaptivityLabel(*adaptedDM,label);CHKERRQ(ierr);
163809350103SToby Isaac   ierr = DMSetUp(*adaptedDM);CHKERRQ(ierr);
163909350103SToby Isaac   ierr = DMForestGetAdaptivitySuccess(*adaptedDM,&success);CHKERRQ(ierr);
164009350103SToby Isaac   if (!success) {
164109350103SToby Isaac     ierr = DMDestroy(adaptedDM);CHKERRQ(ierr);
164209350103SToby Isaac     *adaptedDM = NULL;
164309350103SToby Isaac   }
164409350103SToby Isaac   PetscFunctionReturn(0);
164509350103SToby Isaac }
164609350103SToby Isaac 
1647d222f98bSToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm)
1648d222f98bSToby Isaac {
1649d222f98bSToby Isaac   PetscErrorCode ierr;
1650d222f98bSToby Isaac 
1651d222f98bSToby Isaac   PetscFunctionBegin;
1652d222f98bSToby Isaac   ierr = PetscMemzero(dm->ops,sizeof(*(dm->ops)));CHKERRQ(ierr);
1653d222f98bSToby Isaac 
1654d222f98bSToby Isaac   dm->ops->clone          = DMClone_Forest;
1655d222f98bSToby Isaac   dm->ops->setfromoptions = DMSetFromOptions_Forest;
1656d222f98bSToby Isaac   dm->ops->destroy        = DMDestroy_Forest;
1657d8984e3bSMatthew G. Knepley   dm->ops->createsubdm    = DMCreateSubDM_Forest;
16585421bac9SToby Isaac   dm->ops->refine         = DMRefine_Forest;
16595421bac9SToby Isaac   dm->ops->coarsen        = DMCoarsen_Forest;
16600d1cd5e0SMatthew G. Knepley   dm->ops->adaptlabel     = DMAdaptLabel_Forest;
1661d222f98bSToby Isaac   PetscFunctionReturn(0);
1662d222f98bSToby Isaac }
1663d222f98bSToby Isaac 
16649be51f97SToby Isaac /*MC
16659be51f97SToby Isaac 
1666bae1f979SBarry Smith      DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh.  Forests usually have a base DM
1667bae1f979SBarry Smith   (see DMForestGetBaseDM()), from which it is refined.  The refinement and partitioning of forests is considered
1668bae1f979SBarry Smith   immutable after DMSetUp() is called.  To adapt a mesh, one should call DMForestTemplate() to create a new mesh that
1669bae1f979SBarry Smith   will default to being identical to it, specify how that mesh should differ, and then calling DMSetUp() on the new
1670bae1f979SBarry Smith   mesh.
1671bae1f979SBarry Smith 
1672bae1f979SBarry Smith   To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the
1673bae1f979SBarry Smith   previous mesh whose values indicate which cells should be refined (DM_ADAPT_REFINE) or coarsened (DM_ADAPT_COARSEN)
1674bae1f979SBarry Smith   and how (subtypes are free to allow additional values for things like anisotropic refinement).  The label should be
1675bae1f979SBarry Smith   given to the *new* mesh with DMForestSetAdaptivityLabel().
16769be51f97SToby Isaac 
16779be51f97SToby Isaac   Level: advanced
16789be51f97SToby Isaac 
16799be51f97SToby Isaac .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel()
16809be51f97SToby Isaac M*/
16819be51f97SToby Isaac 
1682db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1683db4d5e8cSToby Isaac {
1684db4d5e8cSToby Isaac   DM_Forest      *forest;
1685db4d5e8cSToby Isaac   PetscErrorCode ierr;
1686db4d5e8cSToby Isaac 
1687db4d5e8cSToby Isaac   PetscFunctionBegin;
1688db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1689db4d5e8cSToby Isaac   ierr                         = PetscNewLog(dm,&forest);CHKERRQ(ierr);
1690db4d5e8cSToby Isaac   dm->dim                      = 0;
1691db4d5e8cSToby Isaac   dm->data                     = forest;
1692db4d5e8cSToby Isaac   forest->refct                = 1;
1693db4d5e8cSToby Isaac   forest->data                 = NULL;
169458762b62SToby Isaac   forest->setfromoptionscalled = PETSC_FALSE;
1695db4d5e8cSToby Isaac   forest->topology             = NULL;
1696cd3c525cSToby Isaac   forest->adapt                = NULL;
1697db4d5e8cSToby Isaac   forest->base                 = NULL;
16986a87ffbfSToby Isaac   forest->adaptPurpose         = DM_ADAPT_DETERMINE;
1699db4d5e8cSToby Isaac   forest->adjDim               = PETSC_DEFAULT;
1700db4d5e8cSToby Isaac   forest->overlap              = PETSC_DEFAULT;
1701db4d5e8cSToby Isaac   forest->minRefinement        = PETSC_DEFAULT;
1702db4d5e8cSToby Isaac   forest->maxRefinement        = PETSC_DEFAULT;
170356ba9f64SToby Isaac   forest->initRefinement       = PETSC_DEFAULT;
1704c7eeac06SToby Isaac   forest->cStart               = PETSC_DETERMINE;
1705c7eeac06SToby Isaac   forest->cEnd                 = PETSC_DETERMINE;
1706cd3c525cSToby Isaac   forest->cellSF               = NULL;
1707ebdf65a2SToby Isaac   forest->adaptLabel           = NULL;
1708db4d5e8cSToby Isaac   forest->gradeFactor          = 2;
1709db4d5e8cSToby Isaac   forest->cellWeights          = NULL;
1710db4d5e8cSToby Isaac   forest->cellWeightsCopyMode  = PETSC_USE_POINTER;
1711db4d5e8cSToby Isaac   forest->weightsFactor        = 1.;
1712db4d5e8cSToby Isaac   forest->weightCapacity       = 1.;
1713a73e2921SToby Isaac   ierr                         = DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);CHKERRQ(ierr);
1714d222f98bSToby Isaac   ierr                         = DMInitialize_Forest(dm);CHKERRQ(ierr);
1715db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1716db4d5e8cSToby Isaac }
1717