xref: /petsc/src/dm/impls/forest/forest.c (revision 9be51f9731d718151c5cb83f1f553d7ad20b1d74)
1*9be51f97SToby Isaac #include <petsc/private/dmforestimpl.h> /*I "petscdmforest.h" I*/
2*9be51f97SToby Isaac #include <petsc/private/dmimpl.h>       /*I "petscdm.h" I*/
3ef19d27cSToby Isaac #include <petscsf.h>
4db4d5e8cSToby Isaac 
527d4645fSToby Isaac PetscBool DMForestPackageInitialized = PETSC_FALSE;
627d4645fSToby Isaac 
727d4645fSToby Isaac typedef struct _DMForestTypeLink *DMForestTypeLink;
827d4645fSToby Isaac 
927d4645fSToby Isaac struct _DMForestTypeLink
1027d4645fSToby Isaac {
1127d4645fSToby Isaac   char *name;
1227d4645fSToby Isaac   DMForestTypeLink next;
1327d4645fSToby Isaac };
1427d4645fSToby Isaac 
1527d4645fSToby Isaac DMForestTypeLink DMForestTypeList;
1627d4645fSToby Isaac 
1727d4645fSToby Isaac #undef __FUNCT__
1827d4645fSToby Isaac #define __FUNCT__ "DMForestPackageFinalize"
1927d4645fSToby Isaac static PetscErrorCode DMForestPackageFinalize(void)
2027d4645fSToby Isaac {
2127d4645fSToby Isaac   DMForestTypeLink oldLink, link = DMForestTypeList;
2227d4645fSToby Isaac   PetscErrorCode ierr;
2327d4645fSToby Isaac 
2427d4645fSToby Isaac   PetscFunctionBegin;
2527d4645fSToby Isaac   while (link) {
2627d4645fSToby Isaac     oldLink = link;
2727d4645fSToby Isaac     ierr = PetscFree(oldLink->name);
2827d4645fSToby Isaac     link = oldLink->next;
2927d4645fSToby Isaac     ierr = PetscFree(oldLink);CHKERRQ(ierr);
3027d4645fSToby Isaac   }
3127d4645fSToby Isaac   PetscFunctionReturn(0);
3227d4645fSToby Isaac }
3327d4645fSToby Isaac 
3427d4645fSToby Isaac #undef __FUNCT__
3527d4645fSToby Isaac #define __FUNCT__ "DMForestPackageInitialize"
3627d4645fSToby Isaac static PetscErrorCode DMForestPackageInitialize(void)
3727d4645fSToby Isaac {
3827d4645fSToby Isaac   PetscErrorCode ierr;
3927d4645fSToby Isaac 
4027d4645fSToby Isaac   PetscFunctionBegin;
4127d4645fSToby Isaac   if (DMForestPackageInitialized) PetscFunctionReturn(0);
4227d4645fSToby Isaac   DMForestPackageInitialized = PETSC_TRUE;
4327d4645fSToby Isaac   ierr = DMForestRegisterType(DMFOREST);CHKERRQ(ierr);
4427d4645fSToby Isaac   ierr = PetscRegisterFinalize(DMForestPackageFinalize);CHKERRQ(ierr);
4527d4645fSToby Isaac   PetscFunctionReturn(0);
4627d4645fSToby Isaac }
4727d4645fSToby Isaac 
4827d4645fSToby Isaac #undef __FUNCT__
4927d4645fSToby Isaac #define __FUNCT__ "DMForestRegisterType"
50*9be51f97SToby Isaac /*@C
51*9be51f97SToby Isaac   DMForestRegisterType - Registers a DMType as a subtype of DMFOREST (so that DMIsForest() will be correct)
52*9be51f97SToby Isaac 
53*9be51f97SToby Isaac   Not Collective
54*9be51f97SToby Isaac 
55*9be51f97SToby Isaac   Input parameter:
56*9be51f97SToby Isaac . name - the name of the type
57*9be51f97SToby Isaac 
58*9be51f97SToby Isaac   Level: advanced
59*9be51f97SToby Isaac 
60*9be51f97SToby Isaac .seealso: DMFOREST, DMIsForest()
61*9be51f97SToby Isaac @*/
6227d4645fSToby Isaac PetscErrorCode DMForestRegisterType(DMType name)
6327d4645fSToby Isaac {
6427d4645fSToby Isaac   DMForestTypeLink link;
6527d4645fSToby Isaac   PetscErrorCode ierr;
6627d4645fSToby Isaac 
6727d4645fSToby Isaac   PetscFunctionBegin;
6827d4645fSToby Isaac   ierr = DMForestPackageInitialize();CHKERRQ(ierr);
6927d4645fSToby Isaac   ierr = PetscNew(&link);CHKERRQ(ierr);
7027d4645fSToby Isaac   ierr = PetscStrallocpy(name,&link->name);CHKERRQ(ierr);
7127d4645fSToby Isaac   link->next = DMForestTypeList;
7227d4645fSToby Isaac   DMForestTypeList = link;
7327d4645fSToby Isaac   PetscFunctionReturn(0);
7427d4645fSToby Isaac }
7527d4645fSToby Isaac 
7627d4645fSToby Isaac #undef __FUNCT__
7727d4645fSToby Isaac #define __FUNCT__ "DMIsForest"
78*9be51f97SToby Isaac /*@
79*9be51f97SToby Isaac   DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes
80*9be51f97SToby Isaac 
81*9be51f97SToby Isaac   Not Collective
82*9be51f97SToby Isaac 
83*9be51f97SToby Isaac   Input parameter:
84*9be51f97SToby Isaac . dm - the DM object
85*9be51f97SToby Isaac 
86*9be51f97SToby Isaac   Output parameter:
87*9be51f97SToby Isaac . isForest - whether dm is a subtype of DMFOREST
88*9be51f97SToby Isaac 
89*9be51f97SToby Isaac   Level: intermediate
90*9be51f97SToby Isaac 
91*9be51f97SToby Isaac .seealso: DMFOREST, DMForestRegisterType()
92*9be51f97SToby Isaac @*/
9327d4645fSToby Isaac PetscErrorCode DMIsForest(DM dm, PetscBool *isForest)
9427d4645fSToby Isaac {
9527d4645fSToby Isaac   DMForestTypeLink link = DMForestTypeList;
9627d4645fSToby Isaac   PetscErrorCode ierr;
9727d4645fSToby Isaac 
9827d4645fSToby Isaac   PetscFunctionBegin;
9927d4645fSToby Isaac   while (link) {
10027d4645fSToby Isaac     PetscBool sameType;
10127d4645fSToby Isaac     ierr = PetscObjectTypeCompare((PetscObject)dm,link->name,&sameType);CHKERRQ(ierr);
10227d4645fSToby Isaac     if (sameType) {
10327d4645fSToby Isaac       *isForest = PETSC_TRUE;
10427d4645fSToby Isaac       PetscFunctionReturn(0);
10527d4645fSToby Isaac     }
10627d4645fSToby Isaac     link = link->next;
10727d4645fSToby Isaac   }
10827d4645fSToby Isaac   *isForest = PETSC_FALSE;
10927d4645fSToby Isaac   PetscFunctionReturn(0);
11027d4645fSToby Isaac }
11127d4645fSToby Isaac 
112db4d5e8cSToby Isaac #undef __FUNCT__
113a0452a8eSToby Isaac #define __FUNCT__ "DMForestTemplate"
114*9be51f97SToby Isaac /*@
115*9be51f97SToby Isaac   DMForestTemplate - Create a new DM that will be adapted from a source DM.  The new DM reproduces the configuration
116*9be51f97SToby 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
117*9be51f97SToby Isaac   (by, e.g., refinement or repartitioning).  The source DM is also set as the adaptivity source DM of the new DM (see
118*9be51f97SToby Isaac   DMForestSetAdaptivityForest()).
119*9be51f97SToby Isaac 
120*9be51f97SToby Isaac   Collective on dm
121*9be51f97SToby Isaac 
122*9be51f97SToby Isaac   Input Parameters:
123*9be51f97SToby Isaac + dm - the source DM object
124*9be51f97SToby Isaac - comm - the communicator for the new DM (this communicator is currently ignored, but is present so that DMForestTemplate() can be used within DMCoarsen())
125*9be51f97SToby Isaac 
126*9be51f97SToby Isaac   Output Parameter:
127*9be51f97SToby Isaac . tdm - the new DM object
128*9be51f97SToby Isaac 
129*9be51f97SToby Isaac   Level: intermediate
130*9be51f97SToby Isaac 
131*9be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest()
132*9be51f97SToby Isaac @*/
133*9be51f97SToby Isaac PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm)
134a0452a8eSToby Isaac {
135a0452a8eSToby Isaac   DM_Forest        *forest = (DM_Forest *) dm->data;
13620e8089bSToby Isaac   DMType           type;
137a0452a8eSToby Isaac   DM               base;
138a0452a8eSToby Isaac   DMForestTopology topology;
139a0452a8eSToby Isaac   PetscInt         dim, overlap, ref, factor;
140a0452a8eSToby Isaac   DMForestAdaptivityStrategy strat;
141795844e7SToby Isaac   PetscDS          ds;
142795844e7SToby Isaac   void             *ctx;
143a0452a8eSToby Isaac   PetscErrorCode   ierr;
144a0452a8eSToby Isaac 
145a0452a8eSToby Isaac   PetscFunctionBegin;
146a0452a8eSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
14720e8089bSToby Isaac   ierr = DMCreate(PetscObjectComm((PetscObject)dm),tdm);CHKERRQ(ierr);
14820e8089bSToby Isaac   ierr = DMGetType(dm,&type);CHKERRQ(ierr);
14920e8089bSToby Isaac   ierr = DMSetType(*tdm,type);CHKERRQ(ierr);
150a0452a8eSToby Isaac   ierr = DMForestGetBaseDM(dm,&base);CHKERRQ(ierr);
15120e8089bSToby Isaac   ierr = DMForestSetBaseDM(*tdm,base);CHKERRQ(ierr);
152a0452a8eSToby Isaac   ierr = DMForestGetTopology(dm,&topology);CHKERRQ(ierr);
15320e8089bSToby Isaac   ierr = DMForestSetTopology(*tdm,topology);CHKERRQ(ierr);
154a0452a8eSToby Isaac   ierr = DMForestGetAdjacencyDimension(dm,&dim);CHKERRQ(ierr);
15520e8089bSToby Isaac   ierr = DMForestSetAdjacencyDimension(*tdm,dim);CHKERRQ(ierr);
156a0452a8eSToby Isaac   ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr);
15720e8089bSToby Isaac   ierr = DMForestSetPartitionOverlap(*tdm,overlap);CHKERRQ(ierr);
158a0452a8eSToby Isaac   ierr = DMForestGetMinimumRefinement(dm,&ref);CHKERRQ(ierr);
15920e8089bSToby Isaac   ierr = DMForestSetMinimumRefinement(*tdm,ref);CHKERRQ(ierr);
160a0452a8eSToby Isaac   ierr = DMForestGetMaximumRefinement(dm,&ref);CHKERRQ(ierr);
16120e8089bSToby Isaac   ierr = DMForestSetMaximumRefinement(*tdm,ref);CHKERRQ(ierr);
162a0452a8eSToby Isaac   ierr = DMForestGetAdaptivityStrategy(dm,&strat);CHKERRQ(ierr);
16320e8089bSToby Isaac   ierr = DMForestSetAdaptivityStrategy(*tdm,strat);CHKERRQ(ierr);
164a0452a8eSToby Isaac   ierr = DMForestGetGradeFactor(dm,&factor);CHKERRQ(ierr);
16520e8089bSToby Isaac   ierr = DMForestSetGradeFactor(*tdm,factor);CHKERRQ(ierr);
166a0452a8eSToby Isaac   if (forest->ftemplate) {
16720e8089bSToby Isaac     ierr = (forest->ftemplate) (dm, *tdm);CHKERRQ(ierr);
168a0452a8eSToby Isaac   }
16920e8089bSToby Isaac   ierr = DMForestSetAdaptivityForest(*tdm,dm);CHKERRQ(ierr);
170795844e7SToby Isaac   ierr = DMGetDS(dm,&ds);CHKERRQ(ierr);
171795844e7SToby Isaac   ierr = DMSetDS(*tdm,ds);CHKERRQ(ierr);
172795844e7SToby Isaac   ierr = DMGetApplicationContext(dm,&ctx);CHKERRQ(ierr);
173795844e7SToby Isaac   ierr = DMSetApplicationContext(*tdm,&ctx);CHKERRQ(ierr);
174795844e7SToby Isaac   if (dm->maxCell) {
175795844e7SToby Isaac     const PetscReal *maxCell, *L;
176795844e7SToby Isaac     const DMBoundaryType *bd;
177795844e7SToby Isaac 
178795844e7SToby Isaac     ierr = DMGetPeriodicity(dm,&maxCell,&L,&bd);CHKERRQ(ierr);
179795844e7SToby Isaac     ierr = DMSetPeriodicity(*tdm,maxCell,L,bd);CHKERRQ(ierr);
180795844e7SToby Isaac   }
181a0452a8eSToby Isaac   PetscFunctionReturn(0);
182a0452a8eSToby Isaac }
183a0452a8eSToby Isaac 
18401d9d024SToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm);
18501d9d024SToby Isaac 
186a0452a8eSToby Isaac #undef __FUNCT__
187db4d5e8cSToby Isaac #define __FUNCT__ "DMClone_Forest"
188db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm)
189db4d5e8cSToby Isaac {
190db4d5e8cSToby Isaac   DM_Forest        *forest = (DM_Forest *) dm->data;
191db4d5e8cSToby Isaac   const char       *type;
192db4d5e8cSToby Isaac   PetscErrorCode ierr;
193db4d5e8cSToby Isaac 
194db4d5e8cSToby Isaac   PetscFunctionBegin;
195db4d5e8cSToby Isaac   forest->refct++;
196db4d5e8cSToby Isaac   (*newdm)->data = forest;
197db4d5e8cSToby Isaac   ierr = PetscObjectGetType((PetscObject) dm, &type);CHKERRQ(ierr);
198db4d5e8cSToby Isaac   ierr = PetscObjectChangeTypeName((PetscObject) *newdm, type);CHKERRQ(ierr);
19901d9d024SToby Isaac   ierr = DMInitialize_Forest(*newdm);CHKERRQ(ierr);
200db4d5e8cSToby Isaac   PetscFunctionReturn(0);
201db4d5e8cSToby Isaac }
202db4d5e8cSToby Isaac 
203db4d5e8cSToby Isaac #undef __FUNCT__
204db4d5e8cSToby Isaac #define __FUNCT__ "DMDestroy_Forest"
205d222f98bSToby Isaac static PetscErrorCode DMDestroy_Forest(DM dm)
206db4d5e8cSToby Isaac {
207db4d5e8cSToby Isaac   DM_Forest     *forest = (DM_Forest*) dm->data;
208db4d5e8cSToby Isaac   PetscErrorCode ierr;
209db4d5e8cSToby Isaac 
210db4d5e8cSToby Isaac   PetscFunctionBegin;
211db4d5e8cSToby Isaac   if (--forest->refct > 0) PetscFunctionReturn(0);
212d222f98bSToby Isaac   if (forest->destroy) {ierr = forest->destroy(dm);CHKERRQ(ierr);}
213db4d5e8cSToby Isaac   ierr = PetscSFDestroy(&forest->cellSF);CHKERRQ(ierr);
214ebdf65a2SToby Isaac   ierr = PetscFree(forest->adaptLabel);CHKERRQ(ierr);
2159a81d013SToby Isaac   ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr);
21656ba9f64SToby Isaac   ierr = DMDestroy(&forest->base);CHKERRQ(ierr);
217ba936b91SToby Isaac   ierr = DMDestroy(&forest->adapt);CHKERRQ(ierr);
21830f902e7SToby Isaac   ierr = PetscFree(forest->topology);CHKERRQ(ierr);
21930f902e7SToby Isaac   ierr = PetscFree(forest);CHKERRQ(ierr);
220db4d5e8cSToby Isaac   PetscFunctionReturn(0);
221db4d5e8cSToby Isaac }
222db4d5e8cSToby Isaac 
223db4d5e8cSToby Isaac #undef __FUNCT__
224dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetTopology"
225*9be51f97SToby Isaac /*@C
226*9be51f97SToby Isaac   DMForestSetTopology - Set the topology of a DMForest during the pre-setup phase.  The topology is a string (e.g.
227*9be51f97SToby Isaac   "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest durint
228*9be51f97SToby Isaac   DMSetUp().
229*9be51f97SToby Isaac 
230*9be51f97SToby Isaac   Logically collective on dm
231*9be51f97SToby Isaac 
232*9be51f97SToby Isaac   Input parameters:
233*9be51f97SToby Isaac + dm - the forest
234*9be51f97SToby Isaac - topology - the topology of the forest
235*9be51f97SToby Isaac 
236*9be51f97SToby Isaac   Level: intermediate
237*9be51f97SToby Isaac 
238*9be51f97SToby Isaac .seealso(): DMForestGetTopology(), DMForestSetBaseDM()
239*9be51f97SToby Isaac @*/
240dd8e54a2SToby Isaac PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology)
241db4d5e8cSToby Isaac {
242db4d5e8cSToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
243db4d5e8cSToby Isaac   PetscErrorCode ierr;
244db4d5e8cSToby Isaac 
245db4d5e8cSToby Isaac   PetscFunctionBegin;
246db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
247ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the topology after setup");
248dd8e54a2SToby Isaac   ierr = PetscFree(forest->topology);CHKERRQ(ierr);
249dd8e54a2SToby Isaac   ierr = PetscStrallocpy((const char *)topology,(char **) &forest->topology);CHKERRQ(ierr);
250db4d5e8cSToby Isaac   PetscFunctionReturn(0);
251db4d5e8cSToby Isaac }
252db4d5e8cSToby Isaac 
253db4d5e8cSToby Isaac #undef __FUNCT__
254dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetTopology"
255*9be51f97SToby Isaac /*@C
256*9be51f97SToby Isaac   DMForestGetTopology - Get a string describing the topology of a DMForest.
257*9be51f97SToby Isaac 
258*9be51f97SToby Isaac   Not collective
259*9be51f97SToby Isaac 
260*9be51f97SToby Isaac   Input parameter:
261*9be51f97SToby Isaac . dm - the forest
262*9be51f97SToby Isaac 
263*9be51f97SToby Isaac   Output parameter:
264*9be51f97SToby Isaac . topology - the topology of the forest (e.g., 'cube', 'shell')
265*9be51f97SToby Isaac 
266*9be51f97SToby Isaac   Level: intermediate
267*9be51f97SToby Isaac 
268*9be51f97SToby Isaac .seealso: DMForestSetTopology()
269*9be51f97SToby Isaac @*/
270dd8e54a2SToby Isaac PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology)
271dd8e54a2SToby Isaac {
272dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
273dd8e54a2SToby Isaac 
274dd8e54a2SToby Isaac   PetscFunctionBegin;
275dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
276dd8e54a2SToby Isaac   PetscValidPointer(topology,2);
277dd8e54a2SToby Isaac   *topology = forest->topology;
278dd8e54a2SToby Isaac   PetscFunctionReturn(0);
279dd8e54a2SToby Isaac }
280dd8e54a2SToby Isaac 
281dd8e54a2SToby Isaac #undef __FUNCT__
282dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetBaseDM"
283*9be51f97SToby Isaac /*@
284*9be51f97SToby Isaac   DMForestSetBaseDM - During the pre-setup phase, set the DM that defines the base mesh of a DMForest forest.  The
285*9be51f97SToby Isaac   forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its
286*9be51f97SToby Isaac   base.  In general, two forest must share a bse to be comparable, to do things like construct interpolators.
287*9be51f97SToby Isaac 
288*9be51f97SToby Isaac   Logically collective on dm
289*9be51f97SToby Isaac 
290*9be51f97SToby Isaac   Input Parameters:
291*9be51f97SToby Isaac + dm - the forest
292*9be51f97SToby Isaac - base - the base DM of the forest
293*9be51f97SToby Isaac 
294*9be51f97SToby Isaac   Level: intermediate
295*9be51f97SToby Isaac 
296*9be51f97SToby Isaac .seealso(): DMForestGetBaseDM()
297*9be51f97SToby Isaac @*/
298dd8e54a2SToby Isaac PetscErrorCode DMForestSetBaseDM(DM dm, DM base)
299dd8e54a2SToby Isaac {
300dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
301dd8e54a2SToby Isaac   PetscInt       dim, dimEmbed;
302dd8e54a2SToby Isaac   PetscErrorCode ierr;
303dd8e54a2SToby Isaac 
304dd8e54a2SToby Isaac   PetscFunctionBegin;
305dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
306ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the base after setup");
307dd8e54a2SToby Isaac   ierr = PetscObjectReference((PetscObject)base);CHKERRQ(ierr);
308dd8e54a2SToby Isaac   ierr = DMDestroy(&forest->base);CHKERRQ(ierr);
309dd8e54a2SToby Isaac   forest->base = base;
310a0452a8eSToby Isaac   if (base) {
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);
316a0452a8eSToby Isaac   }
317dd8e54a2SToby Isaac   PetscFunctionReturn(0);
318dd8e54a2SToby Isaac }
319dd8e54a2SToby Isaac 
320dd8e54a2SToby Isaac #undef __FUNCT__
321dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetBaseDM"
322*9be51f97SToby Isaac /*@
323*9be51f97SToby Isaac   DMForestGetBaseDM - Get the base DM of a DMForest forest.  The forest will be hierarchically refined from the base,
324*9be51f97SToby Isaac   and all refinements/coarsenings of the forest will share its base.  In general, two forest must share a bse to be
325*9be51f97SToby Isaac   comparable, to do things like construct interpolators.
326*9be51f97SToby Isaac 
327*9be51f97SToby Isaac   Not collective
328*9be51f97SToby Isaac 
329*9be51f97SToby Isaac   Input Parameter:
330*9be51f97SToby Isaac . dm - the forest
331*9be51f97SToby Isaac 
332*9be51f97SToby Isaac   Output Parameter:
333*9be51f97SToby Isaac . base - the base DM of the forest
334*9be51f97SToby Isaac 
335*9be51f97SToby Isaac   Level: intermediate
336*9be51f97SToby Isaac 
337*9be51f97SToby Isaac .seealso(); DMForestSetBaseDM()
338*9be51f97SToby Isaac @*/
339dd8e54a2SToby Isaac PetscErrorCode DMForestGetBaseDM(DM dm, DM *base)
340dd8e54a2SToby Isaac {
341dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
342dd8e54a2SToby Isaac 
343dd8e54a2SToby Isaac   PetscFunctionBegin;
344dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
345dd8e54a2SToby Isaac   PetscValidPointer(base, 2);
346dd8e54a2SToby Isaac   *base = forest->base;
347dd8e54a2SToby Isaac   PetscFunctionReturn(0);
348dd8e54a2SToby Isaac }
349dd8e54a2SToby Isaac 
350dd8e54a2SToby Isaac #undef __FUNCT__
351ba936b91SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityForest"
352*9be51f97SToby Isaac /*@
353*9be51f97SToby Isaac   DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be
354*9be51f97SToby Isaac   adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) im DMSetUp().  Usually not needed
355*9be51f97SToby Isaac   by users directly: DMForestTemplate() constructs a new forest to be adapted from an old forest and calls this
356*9be51f97SToby Isaac   routine.
357*9be51f97SToby Isaac 
358*9be51f97SToby Isaac   Logically collective on dm
359*9be51f97SToby Isaac 
360*9be51f97SToby Isaac   Input Parameter:
361*9be51f97SToby Isaac + dm - the new forest, which will be constructed from adapt
362*9be51f97SToby Isaac - adapt - the old forest
363*9be51f97SToby Isaac 
364*9be51f97SToby Isaac   Level: intermediate
365*9be51f97SToby Isaac 
366*9be51f97SToby Isaac .seealso: DMForestGetAdaptivityForest(), DMForestSetAdaptivityPurpose()
367*9be51f97SToby Isaac @*/
368ba936b91SToby Isaac PetscErrorCode DMForestSetAdaptivityForest(DM dm,DM adapt)
369dd8e54a2SToby Isaac {
370ba936b91SToby Isaac   DM_Forest        *forest;
371dd8e54a2SToby Isaac   PetscErrorCode   ierr;
372dd8e54a2SToby Isaac 
373dd8e54a2SToby Isaac   PetscFunctionBegin;
374dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
375ba936b91SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
376ba936b91SToby Isaac   forest = (DM_Forest *) dm->data;
377ba936b91SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
37826d9498aSToby Isaac   switch (forest->adaptPurpose) {
37926d9498aSToby Isaac   case DM_FOREST_KEEP:
380ba936b91SToby Isaac     ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr);
381ba936b91SToby Isaac     ierr = DMDestroy(&(forest->adapt));CHKERRQ(ierr);
382ba936b91SToby Isaac     forest->adapt = adapt;
38326d9498aSToby Isaac     break;
38426d9498aSToby Isaac   case DM_FOREST_REFINE:
38526d9498aSToby Isaac     ierr = DMSetCoarseDM(dm,adapt);CHKERRQ(ierr);
38626d9498aSToby Isaac     break;
38726d9498aSToby Isaac   case DM_FOREST_COARSEN:
38826d9498aSToby Isaac     ierr = DMSetFineDM(dm,adapt);CHKERRQ(ierr);
38926d9498aSToby Isaac     break;
39026d9498aSToby Isaac   default:
39126d9498aSToby Isaac     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
39226d9498aSToby Isaac   }
393dd8e54a2SToby Isaac   PetscFunctionReturn(0);
394dd8e54a2SToby Isaac }
395dd8e54a2SToby Isaac 
396dd8e54a2SToby Isaac #undef __FUNCT__
397ba936b91SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityForest"
398*9be51f97SToby Isaac /*@
399*9be51f97SToby Isaac   DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted.
400*9be51f97SToby Isaac 
401*9be51f97SToby Isaac   Not collective
402*9be51f97SToby Isaac 
403*9be51f97SToby Isaac   Input Parameter:
404*9be51f97SToby Isaac . dm - the forest
405*9be51f97SToby Isaac 
406*9be51f97SToby Isaac   Output Parameter:
407*9be51f97SToby Isaac . adapt - the forest from which dm is/was adapted
408*9be51f97SToby Isaac 
409*9be51f97SToby Isaac   Level: intermediate
410*9be51f97SToby Isaac 
411*9be51f97SToby Isaac .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose()
412*9be51f97SToby Isaac @*/
413ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt)
414dd8e54a2SToby Isaac {
415ba936b91SToby Isaac   DM_Forest        *forest;
41626d9498aSToby Isaac   PetscErrorCode   ierr;
417dd8e54a2SToby Isaac 
418dd8e54a2SToby Isaac   PetscFunctionBegin;
419dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
420ba936b91SToby Isaac   forest = (DM_Forest *) dm->data;
42126d9498aSToby Isaac   switch (forest->adaptPurpose) {
42226d9498aSToby Isaac   case DM_FOREST_KEEP:
423ba936b91SToby Isaac     *adapt = forest->adapt;
42426d9498aSToby Isaac     break;
42526d9498aSToby Isaac   case DM_FOREST_REFINE:
42626d9498aSToby Isaac     ierr = DMGetCoarseDM(dm,adapt);CHKERRQ(ierr);
42726d9498aSToby Isaac     break;
42826d9498aSToby Isaac   case DM_FOREST_COARSEN:
42926d9498aSToby Isaac     ierr = DMGetFineDM(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   }
43426d9498aSToby Isaac   PetscFunctionReturn(0);
43526d9498aSToby Isaac }
43626d9498aSToby Isaac 
43726d9498aSToby Isaac #undef __FUNCT__
43826d9498aSToby Isaac #define __FUNCT__ "DMForestSetAdaptivityPurpose"
439*9be51f97SToby Isaac /*@
440*9be51f97SToby Isaac   DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its
441*9be51f97SToby Isaac   source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_FOREST_REFINE), coarsening
442*9be51f97SToby Isaac   (DM_FOREST_COARSEN), or undefined (DM_FOREST_NONE).  This only matters for the purposes of reference counting:
443*9be51f97SToby Isaac   during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse
444*9be51f97SToby Isaac   relationship (see DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does
445*9be51f97SToby Isaac   not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can
446*9be51f97SToby Isaac   cause a memory leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
447*9be51f97SToby Isaac 
448*9be51f97SToby Isaac   Logically collective on dm
449*9be51f97SToby Isaac 
450*9be51f97SToby Isaac   Input Parameters:
451*9be51f97SToby Isaac + dm - the forest
452*9be51f97SToby Isaac - purpose - the adaptivity purpose (DM_FOREST_NONE/DM_FOREST_REFINE/DM_FOREST_COARSEN)
453*9be51f97SToby Isaac 
454*9be51f97SToby Isaac   Level: advanced
455*9be51f97SToby Isaac 
456*9be51f97SToby Isaac .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest()
457*9be51f97SToby Isaac @*/
45826d9498aSToby Isaac PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMForestAdaptivityPurpose purpose)
45926d9498aSToby Isaac {
46026d9498aSToby Isaac   DM_Forest      *forest;
46126d9498aSToby Isaac   PetscErrorCode ierr;
46226d9498aSToby Isaac 
46326d9498aSToby Isaac   PetscFunctionBegin;
46426d9498aSToby Isaac   forest = (DM_Forest *) dm->data;
465*9be51f97SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
46626d9498aSToby Isaac   if (purpose != forest->adaptPurpose) {
46726d9498aSToby Isaac     DM adapt;
46826d9498aSToby Isaac 
46926d9498aSToby Isaac     ierr = DMForestGetAdaptivityForest(dm,&adapt);CHKERRQ(ierr);
47026d9498aSToby Isaac     ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr);
47126d9498aSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
47226d9498aSToby Isaac     forest->adaptPurpose = purpose;
47326d9498aSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,adapt);CHKERRQ(ierr);
47426d9498aSToby Isaac     ierr = DMDestroy(&adapt);CHKERRQ(ierr);
47526d9498aSToby Isaac   }
476dd8e54a2SToby Isaac   PetscFunctionReturn(0);
477dd8e54a2SToby Isaac }
478dd8e54a2SToby Isaac 
479dd8e54a2SToby Isaac #undef __FUNCT__
480dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetAdjacencyDimension"
481*9be51f97SToby Isaac /*@
482*9be51f97SToby Isaac   DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine
483*9be51f97SToby Isaac   cell adjacency (for the purposes of partitioning and overlap).
484*9be51f97SToby Isaac 
485*9be51f97SToby Isaac   Logically collective on dm
486*9be51f97SToby Isaac 
487*9be51f97SToby Isaac   Input Parameters:
488*9be51f97SToby Isaac + dm - the forest
489*9be51f97SToby Isaac - adjDim - default 0 (i.e., vertices determine adjacency)
490*9be51f97SToby Isaac 
491*9be51f97SToby Isaac   Level: intermediate
492*9be51f97SToby Isaac 
493*9be51f97SToby Isaac .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap()
494*9be51f97SToby Isaac @*/
495dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim)
496dd8e54a2SToby Isaac {
497dd8e54a2SToby Isaac   PetscInt        dim;
498dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
499dd8e54a2SToby Isaac   PetscErrorCode  ierr;
500dd8e54a2SToby Isaac 
501dd8e54a2SToby Isaac   PetscFunctionBegin;
502dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
503ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup");
504dd8e54a2SToby Isaac   if (adjDim < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %d", adjDim);
505dd8e54a2SToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
506dd8e54a2SToby Isaac   if (adjDim > dim) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %d: %d", dim, adjDim);
507dd8e54a2SToby Isaac   forest->adjDim = adjDim;
508dd8e54a2SToby Isaac   PetscFunctionReturn(0);
509dd8e54a2SToby Isaac }
510dd8e54a2SToby Isaac 
511dd8e54a2SToby Isaac #undef __FUNCT__
512dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetAdjacencyCodimension"
513*9be51f97SToby Isaac /*@
514*9be51f97SToby Isaac   DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that,
515*9be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
516*9be51f97SToby Isaac 
517*9be51f97SToby Isaac   Logically collective on dm
518*9be51f97SToby Isaac 
519*9be51f97SToby Isaac   Input Parameters:
520*9be51f97SToby Isaac + dm - the forest
521*9be51f97SToby Isaac - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
522*9be51f97SToby Isaac 
523*9be51f97SToby Isaac   Level: intermediate
524*9be51f97SToby Isaac 
525*9be51f97SToby Isaac .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension()
526*9be51f97SToby Isaac @*/
527dd8e54a2SToby Isaac PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim)
528dd8e54a2SToby Isaac {
529dd8e54a2SToby Isaac   PetscInt        dim;
530dd8e54a2SToby Isaac   PetscErrorCode  ierr;
531dd8e54a2SToby Isaac 
532dd8e54a2SToby Isaac   PetscFunctionBegin;
533dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
534dd8e54a2SToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
535dd8e54a2SToby Isaac   ierr = DMForestSetAdjacencyDimension(dm,dim-adjCodim);CHKERRQ(ierr);
536dd8e54a2SToby Isaac   PetscFunctionReturn(0);
537dd8e54a2SToby Isaac }
538dd8e54a2SToby Isaac 
539dd8e54a2SToby Isaac #undef __FUNCT__
540dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetAdjacencyDimension"
541*9be51f97SToby Isaac /*@
542*9be51f97SToby Isaac   DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the
543*9be51f97SToby Isaac   purposes of partitioning and overlap).
544*9be51f97SToby Isaac 
545*9be51f97SToby Isaac   Not collective
546*9be51f97SToby Isaac 
547*9be51f97SToby Isaac   Input Parameter:
548*9be51f97SToby Isaac . dm - the forest
549*9be51f97SToby Isaac 
550*9be51f97SToby Isaac   Output Parameter:
551*9be51f97SToby Isaac . adjDim - default 0 (i.e., vertices determine adjacency)
552*9be51f97SToby Isaac 
553*9be51f97SToby Isaac   Level: intermediate
554*9be51f97SToby Isaac 
555*9be51f97SToby Isaac .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap()
556*9be51f97SToby Isaac @*/
557dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim)
558dd8e54a2SToby Isaac {
559dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
560dd8e54a2SToby Isaac 
561dd8e54a2SToby Isaac   PetscFunctionBegin;
562dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
563dd8e54a2SToby Isaac   PetscValidIntPointer(adjDim,2);
564dd8e54a2SToby Isaac   *adjDim = forest->adjDim;
565dd8e54a2SToby Isaac   PetscFunctionReturn(0);
566dd8e54a2SToby Isaac }
567dd8e54a2SToby Isaac 
568dd8e54a2SToby Isaac #undef __FUNCT__
569dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetAdjacencyCodimension"
570*9be51f97SToby Isaac /*@
571*9be51f97SToby Isaac   DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that,
572*9be51f97SToby Isaac   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
573*9be51f97SToby Isaac 
574*9be51f97SToby Isaac   Not collective
575*9be51f97SToby Isaac 
576*9be51f97SToby Isaac   Input Parameter:
577*9be51f97SToby Isaac . dm - the forest
578*9be51f97SToby Isaac 
579*9be51f97SToby Isaac   Output Parameter:
580*9be51f97SToby Isaac . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
581*9be51f97SToby Isaac 
582*9be51f97SToby Isaac   Level: intermediate
583*9be51f97SToby Isaac 
584*9be51f97SToby Isaac .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension()
585*9be51f97SToby Isaac @*/
586dd8e54a2SToby Isaac PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim)
587dd8e54a2SToby Isaac {
588dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
589dd8e54a2SToby Isaac   PetscInt       dim;
590dd8e54a2SToby Isaac   PetscErrorCode ierr;
591dd8e54a2SToby Isaac 
592dd8e54a2SToby Isaac   PetscFunctionBegin;
593dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
594dd8e54a2SToby Isaac   PetscValidIntPointer(adjCodim,2);
595dd8e54a2SToby Isaac   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
596dd8e54a2SToby Isaac   *adjCodim = dim - forest->adjDim;
597dd8e54a2SToby Isaac   PetscFunctionReturn(0);
598dd8e54a2SToby Isaac }
599dd8e54a2SToby Isaac 
600dd8e54a2SToby Isaac #undef __FUNCT__
601ef51cf95SToby Isaac #define __FUNCT__ "DMForestSetPartitionOverlap"
602*9be51f97SToby Isaac /*@
603*9be51f97SToby Isaac   DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel
604*9be51f97SToby Isaac   partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding
605*9be51f97SToby Isaac   adjacent cells
606*9be51f97SToby Isaac 
607*9be51f97SToby Isaac   Logically collective on dm
608*9be51f97SToby Isaac 
609*9be51f97SToby Isaac   Input Parameters:
610*9be51f97SToby Isaac + dm - the forest
611*9be51f97SToby Isaac - overlap - default 0
612*9be51f97SToby Isaac 
613*9be51f97SToby Isaac   Level: intermediate
614*9be51f97SToby Isaac 
615*9be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
616*9be51f97SToby Isaac @*/
617dd8e54a2SToby Isaac PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap)
618dd8e54a2SToby Isaac {
619dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
620dd8e54a2SToby Isaac 
621dd8e54a2SToby Isaac   PetscFunctionBegin;
622dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
623ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup");
624dd8e54a2SToby Isaac   if (overlap < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %d", overlap);
625dd8e54a2SToby Isaac   forest->overlap = overlap;
626dd8e54a2SToby Isaac   PetscFunctionReturn(0);
627dd8e54a2SToby Isaac }
628dd8e54a2SToby Isaac 
629dd8e54a2SToby Isaac #undef __FUNCT__
630dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetPartitionOverlap"
631*9be51f97SToby Isaac /*@
632*9be51f97SToby Isaac   DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values
633*9be51f97SToby Isaac   > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells
634*9be51f97SToby Isaac 
635*9be51f97SToby Isaac   Not collective
636*9be51f97SToby Isaac 
637*9be51f97SToby Isaac   Input Parameter:
638*9be51f97SToby Isaac . dm - the forest
639*9be51f97SToby Isaac 
640*9be51f97SToby Isaac   Output Parameter:
641*9be51f97SToby Isaac . overlap - default 0
642*9be51f97SToby Isaac 
643*9be51f97SToby Isaac   Level: intermediate
644*9be51f97SToby Isaac 
645*9be51f97SToby Isaac .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
646*9be51f97SToby Isaac @*/
647dd8e54a2SToby Isaac PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap)
648dd8e54a2SToby Isaac {
649dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
650dd8e54a2SToby Isaac 
651dd8e54a2SToby Isaac   PetscFunctionBegin;
652dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
653dd8e54a2SToby Isaac   PetscValidIntPointer(overlap,2);
654dd8e54a2SToby Isaac   *overlap = forest->overlap;
655dd8e54a2SToby Isaac   PetscFunctionReturn(0);
656dd8e54a2SToby Isaac }
657dd8e54a2SToby Isaac 
658dd8e54a2SToby Isaac #undef __FUNCT__
659dd8e54a2SToby Isaac #define __FUNCT__ "DMForestSetMinimumRefinement"
660*9be51f97SToby Isaac /*@
661*9be51f97SToby Isaac   DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base
662*9be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest
663*9be51f97SToby Isaac   (see DMForestGetAdaptivityForest()) this limits the amount of coarsening.
664*9be51f97SToby Isaac 
665*9be51f97SToby Isaac   Logically collective on dm
666*9be51f97SToby Isaac 
667*9be51f97SToby Isaac   Input Parameters:
668*9be51f97SToby Isaac + dm - the forest
669*9be51f97SToby Isaac - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
670*9be51f97SToby Isaac 
671*9be51f97SToby Isaac   Level: intermediate
672*9be51f97SToby Isaac 
673*9be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
674*9be51f97SToby Isaac @*/
675dd8e54a2SToby Isaac PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement)
676dd8e54a2SToby Isaac {
677dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
678dd8e54a2SToby Isaac 
679dd8e54a2SToby Isaac   PetscFunctionBegin;
680dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
681ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup");
682dd8e54a2SToby Isaac   forest->minRefinement = minRefinement;
683dd8e54a2SToby Isaac   PetscFunctionReturn(0);
684dd8e54a2SToby Isaac }
685dd8e54a2SToby Isaac 
686dd8e54a2SToby Isaac #undef __FUNCT__
687dd8e54a2SToby Isaac #define __FUNCT__ "DMForestGetMinimumRefinement"
688*9be51f97SToby Isaac /*@
689*9be51f97SToby Isaac   DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see
690*9be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest (see
691*9be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of coarsening.
692*9be51f97SToby Isaac 
693*9be51f97SToby Isaac   Not collective
694*9be51f97SToby Isaac 
695*9be51f97SToby Isaac   Input Parameter:
696*9be51f97SToby Isaac . dm - the forest
697*9be51f97SToby Isaac 
698*9be51f97SToby Isaac   Output Parameter:
699*9be51f97SToby Isaac . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
700*9be51f97SToby Isaac 
701*9be51f97SToby Isaac   Level: intermediate
702*9be51f97SToby Isaac 
703*9be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestGetMaximumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
704*9be51f97SToby Isaac @*/
705dd8e54a2SToby Isaac PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement)
706dd8e54a2SToby Isaac {
707dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
708dd8e54a2SToby Isaac 
709dd8e54a2SToby Isaac   PetscFunctionBegin;
710dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
711dd8e54a2SToby Isaac   PetscValidIntPointer(minRefinement,2);
712dd8e54a2SToby Isaac   *minRefinement = forest->minRefinement;
713dd8e54a2SToby Isaac   PetscFunctionReturn(0);
714dd8e54a2SToby Isaac }
715dd8e54a2SToby Isaac 
716dd8e54a2SToby Isaac #undef __FUNCT__
71756ba9f64SToby Isaac #define __FUNCT__ "DMForestSetInitialRefinement"
718*9be51f97SToby Isaac /*@
719*9be51f97SToby Isaac   DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base
720*9be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.
721*9be51f97SToby Isaac 
722*9be51f97SToby Isaac   Logically collective on dm
723*9be51f97SToby Isaac 
724*9be51f97SToby Isaac   Input Parameters:
725*9be51f97SToby Isaac + dm - the forest
726*9be51f97SToby Isaac - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
727*9be51f97SToby Isaac 
728*9be51f97SToby Isaac   Level: intermediate
729*9be51f97SToby Isaac 
730*9be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
731*9be51f97SToby Isaac @*/
73256ba9f64SToby Isaac PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement)
73356ba9f64SToby Isaac {
73456ba9f64SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
73556ba9f64SToby Isaac 
73656ba9f64SToby Isaac   PetscFunctionBegin;
73756ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
73856ba9f64SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup");
73956ba9f64SToby Isaac   forest->initRefinement = initRefinement;
74056ba9f64SToby Isaac   PetscFunctionReturn(0);
74156ba9f64SToby Isaac }
74256ba9f64SToby Isaac 
74356ba9f64SToby Isaac #undef __FUNCT__
74456ba9f64SToby Isaac #define __FUNCT__ "DMForestGetInitialRefinement"
745*9be51f97SToby Isaac /*@
746*9be51f97SToby Isaac   DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see
747*9be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.
748*9be51f97SToby Isaac 
749*9be51f97SToby Isaac   Not collective
750*9be51f97SToby Isaac 
751*9be51f97SToby Isaac   Input Parameter:
752*9be51f97SToby Isaac . dm - the forest
753*9be51f97SToby Isaac 
754*9be51f97SToby Isaac   Output Paramater:
755*9be51f97SToby Isaac . initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
756*9be51f97SToby Isaac 
757*9be51f97SToby Isaac   Level: intermediate
758*9be51f97SToby Isaac 
759*9be51f97SToby Isaac .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
760*9be51f97SToby Isaac @*/
76156ba9f64SToby Isaac PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement)
76256ba9f64SToby Isaac {
76356ba9f64SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
76456ba9f64SToby Isaac 
76556ba9f64SToby Isaac   PetscFunctionBegin;
76656ba9f64SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
76756ba9f64SToby Isaac   PetscValidIntPointer(initRefinement,2);
76856ba9f64SToby Isaac   *initRefinement = forest->initRefinement;
76956ba9f64SToby Isaac   PetscFunctionReturn(0);
77056ba9f64SToby Isaac }
77156ba9f64SToby Isaac 
77256ba9f64SToby Isaac #undef __FUNCT__
773c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetMaximumRefinement"
774*9be51f97SToby Isaac /*@
775*9be51f97SToby Isaac   DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base
776*9be51f97SToby Isaac   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest
777*9be51f97SToby Isaac   (see DMForestGetAdaptivityForest()), this limits the amount of refinement.
778*9be51f97SToby Isaac 
779*9be51f97SToby Isaac   Logically collective on dm
780*9be51f97SToby Isaac 
781*9be51f97SToby Isaac   Input Parameters:
782*9be51f97SToby Isaac + dm - the forest
783*9be51f97SToby Isaac - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
784*9be51f97SToby Isaac 
785*9be51f97SToby Isaac   Level: intermediate
786*9be51f97SToby Isaac 
787*9be51f97SToby Isaac .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityDM()
788*9be51f97SToby Isaac @*/
789c7eeac06SToby Isaac PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement)
790dd8e54a2SToby Isaac {
791dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
792dd8e54a2SToby Isaac 
793dd8e54a2SToby Isaac   PetscFunctionBegin;
794dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
795ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup");
796c7eeac06SToby Isaac   forest->maxRefinement = maxRefinement;
797dd8e54a2SToby Isaac   PetscFunctionReturn(0);
798dd8e54a2SToby Isaac }
799dd8e54a2SToby Isaac 
800dd8e54a2SToby Isaac #undef __FUNCT__
801c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetMaximumRefinement"
802*9be51f97SToby Isaac /*@
803*9be51f97SToby Isaac   DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see
804*9be51f97SToby Isaac   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest (see
805*9be51f97SToby Isaac   DMForestGetAdaptivityForest()), this limits the amount of refinement.
806*9be51f97SToby Isaac 
807*9be51f97SToby Isaac   Not collective
808*9be51f97SToby Isaac 
809*9be51f97SToby Isaac   Input Parameter:
810*9be51f97SToby Isaac . dm - the forest
811*9be51f97SToby Isaac 
812*9be51f97SToby Isaac   Output Parameter:
813*9be51f97SToby Isaac . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
814*9be51f97SToby Isaac 
815*9be51f97SToby Isaac   Level: intermediate
816*9be51f97SToby Isaac 
817*9be51f97SToby Isaac .seealso: DMForestSetMaximumRefinement(), DMForestGetMinimumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
818*9be51f97SToby Isaac @*/
819c7eeac06SToby Isaac PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement)
820dd8e54a2SToby Isaac {
821dd8e54a2SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
822dd8e54a2SToby Isaac 
823dd8e54a2SToby Isaac   PetscFunctionBegin;
824dd8e54a2SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
825c7eeac06SToby Isaac   PetscValidIntPointer(maxRefinement,2);
826c7eeac06SToby Isaac   *maxRefinement = forest->maxRefinement;
827dd8e54a2SToby Isaac   PetscFunctionReturn(0);
828dd8e54a2SToby Isaac }
829c7eeac06SToby Isaac 
830c7eeac06SToby Isaac #undef __FUNCT__
831c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityStrategy"
832*9be51f97SToby Isaac /*@C
833*9be51f97SToby Isaac   DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes.
834*9be51f97SToby Isaac   Subtypes of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree
835*9be51f97SToby Isaac   for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to
836*9be51f97SToby Isaac   specify refinement/coarsening.
837*9be51f97SToby Isaac 
838*9be51f97SToby Isaac   Logically collective on dm
839*9be51f97SToby Isaac 
840*9be51f97SToby Isaac   Input Parameters:
841*9be51f97SToby Isaac + dm - the forest
842*9be51f97SToby Isaac - adaptStrategy - default DMFORESTADAPTALL
843*9be51f97SToby Isaac 
844*9be51f97SToby Isaac   Level: advanced
845*9be51f97SToby Isaac 
846*9be51f97SToby Isaac .seealso: DMForestGetAdaptivityStrategy()
847*9be51f97SToby Isaac @*/
848c7eeac06SToby Isaac PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy)
849c7eeac06SToby Isaac {
850c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
851c7eeac06SToby Isaac   PetscErrorCode ierr;
852c7eeac06SToby Isaac 
853c7eeac06SToby Isaac   PetscFunctionBegin;
854c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
855c7eeac06SToby Isaac   ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr);
856a73e2921SToby Isaac   ierr = PetscStrallocpy((const char *) adaptStrategy,(char **)&forest->adaptStrategy);CHKERRQ(ierr);
857c7eeac06SToby Isaac   PetscFunctionReturn(0);
858c7eeac06SToby Isaac }
859c7eeac06SToby Isaac 
860c7eeac06SToby Isaac #undef __FUNCT__
861c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityStrategy"
862*9be51f97SToby Isaac /*@C
863*9be51f97SToby Isaac   DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes.  Subtypes
864*9be51f97SToby Isaac   of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all
865*9be51f97SToby Isaac   processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only
866*9be51f97SToby Isaac   one process needs to specify refinement/coarsening.
867*9be51f97SToby Isaac 
868*9be51f97SToby Isaac   Not collective
869*9be51f97SToby Isaac 
870*9be51f97SToby Isaac   Input Parameter:
871*9be51f97SToby Isaac . dm - the forest
872*9be51f97SToby Isaac 
873*9be51f97SToby Isaac   Output Parameter:
874*9be51f97SToby Isaac . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL)
875*9be51f97SToby Isaac 
876*9be51f97SToby Isaac   Level: advanced
877*9be51f97SToby Isaac 
878*9be51f97SToby Isaac .seealso: DMForestSetAdaptivityStrategy()
879*9be51f97SToby Isaac @*/
880c7eeac06SToby Isaac PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy)
881c7eeac06SToby Isaac {
882c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
883c7eeac06SToby Isaac 
884c7eeac06SToby Isaac   PetscFunctionBegin;
885c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
886c7eeac06SToby Isaac   PetscValidPointer(adaptStrategy,2);
887c7eeac06SToby Isaac   *adaptStrategy = forest->adaptStrategy;
888c7eeac06SToby Isaac   PetscFunctionReturn(0);
889c7eeac06SToby Isaac }
890c7eeac06SToby Isaac 
891c7eeac06SToby Isaac #undef __FUNCT__
892c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetGradeFactor"
893*9be51f97SToby Isaac /*@
894*9be51f97SToby Isaac   DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to
895*9be51f97SToby Isaac   indicate that the diameter of neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may
896*9be51f97SToby Isaac   only support one particular choice of grading factor.
897*9be51f97SToby Isaac 
898*9be51f97SToby Isaac   Logically collective on dm
899*9be51f97SToby Isaac 
900*9be51f97SToby Isaac   Input Parameters:
901*9be51f97SToby Isaac + dm - the forest
902*9be51f97SToby Isaac - grade - the grading factor
903*9be51f97SToby Isaac 
904*9be51f97SToby Isaac   Level: advanced
905*9be51f97SToby Isaac 
906*9be51f97SToby Isaac .seealso: DMForestGetGradeFactor()
907*9be51f97SToby Isaac @*/
908c7eeac06SToby Isaac PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade)
909c7eeac06SToby Isaac {
910c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
911c7eeac06SToby Isaac 
912c7eeac06SToby Isaac   PetscFunctionBegin;
913c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
914ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup");
915c7eeac06SToby Isaac   forest->gradeFactor = grade;
916c7eeac06SToby Isaac   PetscFunctionReturn(0);
917c7eeac06SToby Isaac }
918c7eeac06SToby Isaac 
919c7eeac06SToby Isaac #undef __FUNCT__
920c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetGradeFactor"
921*9be51f97SToby Isaac /*@
922*9be51f97SToby Isaac   DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of
923*9be51f97SToby Isaac   neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may only support one particular
924*9be51f97SToby Isaac   choice of grading factor.
925*9be51f97SToby Isaac 
926*9be51f97SToby Isaac   Not collective
927*9be51f97SToby Isaac 
928*9be51f97SToby Isaac   Input Parameter:
929*9be51f97SToby Isaac . dm - the forest
930*9be51f97SToby Isaac 
931*9be51f97SToby Isaac   Output Parameter:
932*9be51f97SToby Isaac . grade - the grading factor
933*9be51f97SToby Isaac 
934*9be51f97SToby Isaac   Level: advanced
935*9be51f97SToby Isaac 
936*9be51f97SToby Isaac .seealso: DMForestSetGradeFactor()
937*9be51f97SToby Isaac @*/
938c7eeac06SToby Isaac PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade)
939c7eeac06SToby Isaac {
940c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
941c7eeac06SToby Isaac 
942c7eeac06SToby Isaac   PetscFunctionBegin;
943c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
944c7eeac06SToby Isaac   PetscValidIntPointer(grade,2);
945c7eeac06SToby Isaac   *grade = forest->gradeFactor;
946c7eeac06SToby Isaac   PetscFunctionReturn(0);
947c7eeac06SToby Isaac }
948c7eeac06SToby Isaac 
949c7eeac06SToby Isaac #undef __FUNCT__
950ef51cf95SToby Isaac #define __FUNCT__ "DMForestSetCellWeightFactor"
951*9be51f97SToby Isaac /*@
952*9be51f97SToby Isaac   DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes
953*9be51f97SToby Isaac   the cell weight (see DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be
954*9be51f97SToby Isaac   (cellWeight) * (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on
955*9be51f97SToby Isaac   its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the
956*9be51f97SToby Isaac   computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
957*9be51f97SToby Isaac 
958*9be51f97SToby Isaac   Logically collective on dm
959*9be51f97SToby Isaac 
960*9be51f97SToby Isaac   Input Parameters:
961*9be51f97SToby Isaac + dm - the forest
962*9be51f97SToby Isaac - weightsFactors - default 1.
963*9be51f97SToby Isaac 
964*9be51f97SToby Isaac   Level: advanced
965*9be51f97SToby Isaac 
966*9be51f97SToby Isaac .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights()
967*9be51f97SToby Isaac @*/
968ef51cf95SToby Isaac PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor)
969c7eeac06SToby Isaac {
970c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
971c7eeac06SToby Isaac 
972c7eeac06SToby Isaac   PetscFunctionBegin;
973c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
974ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup");
975c7eeac06SToby Isaac   forest->weightsFactor = weightsFactor;
976c7eeac06SToby Isaac   PetscFunctionReturn(0);
977c7eeac06SToby Isaac }
978c7eeac06SToby Isaac 
979c7eeac06SToby Isaac #undef __FUNCT__
980ef51cf95SToby Isaac #define __FUNCT__ "DMForestGetCellWeightFactor"
981*9be51f97SToby Isaac /*@
982*9be51f97SToby Isaac   DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see
983*9be51f97SToby Isaac   DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be (cellWeight) *
984*9be51f97SToby Isaac   (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on its level; a
985*9be51f97SToby Isaac   factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation
986*9be51f97SToby Isaac   associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
987*9be51f97SToby Isaac 
988*9be51f97SToby Isaac   Not collective
989*9be51f97SToby Isaac 
990*9be51f97SToby Isaac   Input Parameter:
991*9be51f97SToby Isaac . dm - the forest
992*9be51f97SToby Isaac 
993*9be51f97SToby Isaac   Output Parameter:
994*9be51f97SToby Isaac . weightsFactors - default 1.
995*9be51f97SToby Isaac 
996*9be51f97SToby Isaac   Level: advanced
997*9be51f97SToby Isaac 
998*9be51f97SToby Isaac .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights()
999*9be51f97SToby Isaac @*/
1000ef51cf95SToby Isaac PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor)
1001c7eeac06SToby Isaac {
1002c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
1003c7eeac06SToby Isaac 
1004c7eeac06SToby Isaac   PetscFunctionBegin;
1005c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1006c7eeac06SToby Isaac   PetscValidRealPointer(weightsFactor,2);
1007c7eeac06SToby Isaac   *weightsFactor = forest->weightsFactor;
1008c7eeac06SToby Isaac   PetscFunctionReturn(0);
1009c7eeac06SToby Isaac }
1010c7eeac06SToby Isaac 
1011c7eeac06SToby Isaac #undef __FUNCT__
1012c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellChart"
1013*9be51f97SToby Isaac /*@
1014*9be51f97SToby Isaac   DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process
1015*9be51f97SToby Isaac 
1016*9be51f97SToby Isaac   Not collective
1017*9be51f97SToby Isaac 
1018*9be51f97SToby Isaac   Input Parameter:
1019*9be51f97SToby Isaac . dm - the forest
1020*9be51f97SToby Isaac 
1021*9be51f97SToby Isaac   Output Parameters:
1022*9be51f97SToby Isaac + cStart - the first cell on this process
1023*9be51f97SToby Isaac - cEnd - one after the final cell on this process
1024*9be51f97SToby Isaac 
1025*9be51f97SToby Isaac   level: intermediate
1026*9be51f97SToby Isaac 
1027*9be51f97SToby Isaac .seealso: DMForestGetCellSF()
1028*9be51f97SToby Isaac @*/
1029c7eeac06SToby Isaac PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd)
1030c7eeac06SToby Isaac {
1031c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
1032c7eeac06SToby Isaac   PetscErrorCode ierr;
1033c7eeac06SToby Isaac 
1034c7eeac06SToby Isaac   PetscFunctionBegin;
1035c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1036c7eeac06SToby Isaac   PetscValidIntPointer(cStart,2);
1037c7eeac06SToby Isaac   PetscValidIntPointer(cEnd,2);
1038c7eeac06SToby Isaac   if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) {
1039c7eeac06SToby Isaac     ierr = forest->createcellchart(dm,&forest->cStart,&forest->cEnd);CHKERRQ(ierr);
1040c7eeac06SToby Isaac   }
1041c7eeac06SToby Isaac   *cStart =  forest->cStart;
1042c7eeac06SToby Isaac   *cEnd   =  forest->cEnd;
1043c7eeac06SToby Isaac   PetscFunctionReturn(0);
1044c7eeac06SToby Isaac }
1045c7eeac06SToby Isaac 
1046c7eeac06SToby Isaac #undef __FUNCT__
1047c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellSF"
1048*9be51f97SToby Isaac /*@
1049*9be51f97SToby Isaac   DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes
1050*9be51f97SToby Isaac 
1051*9be51f97SToby Isaac   Not collective
1052*9be51f97SToby Isaac 
1053*9be51f97SToby Isaac   Input Parameter:
1054*9be51f97SToby Isaac . dm - the forest
1055*9be51f97SToby Isaac 
1056*9be51f97SToby Isaac   Output Parameter:
1057*9be51f97SToby Isaac . cellSF - the PetscSF
1058*9be51f97SToby Isaac 
1059*9be51f97SToby Isaac   level: intermediate
1060*9be51f97SToby Isaac 
1061*9be51f97SToby Isaac .seealso: DMForestGetCellChart()
1062*9be51f97SToby Isaac @*/
1063c7eeac06SToby Isaac PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF)
1064c7eeac06SToby Isaac {
1065c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
1066c7eeac06SToby Isaac   PetscErrorCode ierr;
1067c7eeac06SToby Isaac 
1068c7eeac06SToby Isaac   PetscFunctionBegin;
1069c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1070c7eeac06SToby Isaac   PetscValidPointer(cellSF,2);
1071c7eeac06SToby Isaac   if ((!forest->cellSF) && forest->createcellsf) {
1072c7eeac06SToby Isaac     ierr = forest->createcellsf(dm,&forest->cellSF);CHKERRQ(ierr);
1073c7eeac06SToby Isaac   }
1074c7eeac06SToby Isaac   *cellSF = forest->cellSF;
1075c7eeac06SToby Isaac   PetscFunctionReturn(0);
1076c7eeac06SToby Isaac }
1077c7eeac06SToby Isaac 
1078c7eeac06SToby Isaac #undef __FUNCT__
1079ebdf65a2SToby Isaac #define __FUNCT__ "DMForestSetAdaptivityLabel"
1080*9be51f97SToby Isaac /*@C
1081*9be51f97SToby Isaac   DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see
1082*9be51f97SToby Isaac   DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination).  The
1083*9be51f97SToby Isaac   interpretation of the label values is up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and
1084*9be51f97SToby Isaac   DM_FOREST_COARSEN have been reserved as choices that should be accepted by all subtypes.
1085*9be51f97SToby Isaac 
1086*9be51f97SToby Isaac   Logically collective on dm
1087*9be51f97SToby Isaac 
1088*9be51f97SToby Isaac   Input Parameters:
1089*9be51f97SToby Isaac - dm - the forest
1090*9be51f97SToby Isaac + adaptLabel - the name of the label in the pre-adaptation forest
1091*9be51f97SToby Isaac 
1092*9be51f97SToby Isaac   Level: intermediate
1093*9be51f97SToby Isaac 
1094*9be51f97SToby Isaac .seealso DMForestGetAdaptivityLabel()
1095*9be51f97SToby Isaac @*/
1096ebdf65a2SToby Isaac PetscErrorCode DMForestSetAdaptivityLabel(DM dm, const char * adaptLabel)
1097c7eeac06SToby Isaac {
1098c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
1099c7eeac06SToby Isaac   PetscErrorCode ierr;
1100c7eeac06SToby Isaac 
1101c7eeac06SToby Isaac   PetscFunctionBegin;
1102c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1103ebdf65a2SToby Isaac   ierr = PetscFree(forest->adaptLabel);CHKERRQ(ierr);
1104ebdf65a2SToby Isaac   ierr = PetscStrallocpy(adaptLabel,&forest->adaptLabel);CHKERRQ(ierr);
1105c7eeac06SToby Isaac   PetscFunctionReturn(0);
1106c7eeac06SToby Isaac }
1107c7eeac06SToby Isaac 
1108c7eeac06SToby Isaac #undef __FUNCT__
1109ebdf65a2SToby Isaac #define __FUNCT__ "DMForestGetAdaptivityLabel"
1110*9be51f97SToby Isaac /*@C
1111*9be51f97SToby Isaac   DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that
1112*9be51f97SToby Isaac   holds the adaptation flags (refinement, coarsening, or some combination).  The interpretation of the label values is
1113*9be51f97SToby Isaac   up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and DM_FOREST_COARSEN have been reserved as
1114*9be51f97SToby Isaac   choices that should be accepted by all subtypes.
1115*9be51f97SToby Isaac 
1116*9be51f97SToby Isaac   Not collective
1117*9be51f97SToby Isaac 
1118*9be51f97SToby Isaac   Input Parameter:
1119*9be51f97SToby Isaac . dm - the forest
1120*9be51f97SToby Isaac 
1121*9be51f97SToby Isaac   Output Parameter:
1122*9be51f97SToby Isaac . adaptLabel - the name of the label in the pre-adaptation forest
1123*9be51f97SToby Isaac 
1124*9be51f97SToby Isaac   Level: intermediate
1125*9be51f97SToby Isaac 
1126*9be51f97SToby Isaac .seealso DMForestSetAdaptivityLabel()
1127*9be51f97SToby Isaac @*/
1128ba936b91SToby Isaac PetscErrorCode DMForestGetAdaptivityLabel(DM dm, const char ** adaptLabel)
1129c7eeac06SToby Isaac {
1130c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
1131c7eeac06SToby Isaac 
1132c7eeac06SToby Isaac   PetscFunctionBegin;
1133c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1134ba936b91SToby Isaac   *adaptLabel = forest->adaptLabel;
1135c7eeac06SToby Isaac   PetscFunctionReturn(0);
1136c7eeac06SToby Isaac }
1137c7eeac06SToby Isaac 
1138c7eeac06SToby Isaac #undef __FUNCT__
1139c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetCellWeights"
1140*9be51f97SToby Isaac /*@
1141*9be51f97SToby Isaac   DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
1142*9be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
1143*9be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
1144*9be51f97SToby Isaac   of 1.
1145*9be51f97SToby Isaac 
1146*9be51f97SToby Isaac   Logically collective on dm
1147*9be51f97SToby Isaac 
1148*9be51f97SToby Isaac   Input Parameters:
1149*9be51f97SToby Isaac + dm - the forest
1150*9be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
1151*9be51f97SToby Isaac - copyMode - how weights should reference weights
1152*9be51f97SToby Isaac 
1153*9be51f97SToby Isaac   Level: advanced
1154*9be51f97SToby Isaac 
1155*9be51f97SToby Isaac .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity()
1156*9be51f97SToby Isaac @*/
1157c7eeac06SToby Isaac PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode)
1158c7eeac06SToby Isaac {
1159c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
1160c7eeac06SToby Isaac   PetscInt       cStart, cEnd;
1161c7eeac06SToby Isaac   PetscErrorCode ierr;
1162c7eeac06SToby Isaac 
1163c7eeac06SToby Isaac   PetscFunctionBegin;
1164c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1165c7eeac06SToby Isaac   ierr = DMForestGetCellChart(dm,&cStart,&cEnd);CHKERRQ(ierr);
1166c7eeac06SToby Isaac   if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd);
1167c7eeac06SToby Isaac   if (copyMode == PETSC_COPY_VALUES) {
1168c7eeac06SToby Isaac     if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) {
1169c7eeac06SToby Isaac       ierr = PetscMalloc1(cEnd-cStart,&forest->cellWeights);CHKERRQ(ierr);
1170c7eeac06SToby Isaac     }
1171c7eeac06SToby Isaac     ierr = PetscMemcpy(forest->cellWeights,weights,(cEnd-cStart)*sizeof(*weights));CHKERRQ(ierr);
1172c7eeac06SToby Isaac     forest->cellWeightsCopyMode = PETSC_OWN_POINTER;
1173c7eeac06SToby Isaac     PetscFunctionReturn(0);
1174c7eeac06SToby Isaac   }
1175c7eeac06SToby Isaac   if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) {
1176c7eeac06SToby Isaac     ierr = PetscFree(forest->cellWeights);CHKERRQ(ierr);
1177c7eeac06SToby Isaac   }
1178c7eeac06SToby Isaac   forest->cellWeights  = weights;
1179c7eeac06SToby Isaac   forest->cellWeightsCopyMode = copyMode;
1180c7eeac06SToby Isaac   PetscFunctionReturn(0);
1181c7eeac06SToby Isaac }
1182c7eeac06SToby Isaac 
1183c7eeac06SToby Isaac #undef __FUNCT__
1184c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetCellWeights"
1185*9be51f97SToby Isaac /*@
1186*9be51f97SToby Isaac   DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
1187*9be51f97SToby Isaac   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
1188*9be51f97SToby Isaac   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
1189*9be51f97SToby Isaac   of 1.
1190*9be51f97SToby Isaac 
1191*9be51f97SToby Isaac   Not collective
1192*9be51f97SToby Isaac 
1193*9be51f97SToby Isaac   Input Parameter:
1194*9be51f97SToby Isaac . dm - the forest
1195*9be51f97SToby Isaac 
1196*9be51f97SToby Isaac   Output Parameter:
1197*9be51f97SToby Isaac . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
1198*9be51f97SToby Isaac 
1199*9be51f97SToby Isaac   Level: advanced
1200*9be51f97SToby Isaac 
1201*9be51f97SToby Isaac .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity()
1202*9be51f97SToby Isaac @*/
1203c7eeac06SToby Isaac PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights)
1204c7eeac06SToby Isaac {
1205c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
1206c7eeac06SToby Isaac 
1207c7eeac06SToby Isaac   PetscFunctionBegin;
1208c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1209c7eeac06SToby Isaac   PetscValidPointer(weights,2);
1210c7eeac06SToby Isaac   *weights = forest->cellWeights;
1211c7eeac06SToby Isaac   PetscFunctionReturn(0);
1212c7eeac06SToby Isaac }
1213c7eeac06SToby Isaac 
1214c7eeac06SToby Isaac #undef __FUNCT__
1215c7eeac06SToby Isaac #define __FUNCT__ "DMForestSetWeightCapacity"
1216*9be51f97SToby Isaac /*@
1217*9be51f97SToby Isaac   DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning
1218*9be51f97SToby Isaac   a pre-adaptation forest (see DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each
1219*9be51f97SToby Isaac   process's cells to the process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that
1220*9be51f97SToby Isaac   the current process should not have any cells after repartitioning.
1221*9be51f97SToby Isaac 
1222*9be51f97SToby Isaac   Logically Collective on dm
1223*9be51f97SToby Isaac 
1224*9be51f97SToby Isaac   Input parameters:
1225*9be51f97SToby Isaac + dm - the forest
1226*9be51f97SToby Isaac - capacity - this process's capacity
1227*9be51f97SToby Isaac 
1228*9be51f97SToby Isaac   Level: advanced
1229*9be51f97SToby Isaac 
1230*9be51f97SToby Isaac .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
1231*9be51f97SToby Isaac @*/
1232c7eeac06SToby Isaac PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity)
1233c7eeac06SToby Isaac {
1234c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
1235c7eeac06SToby Isaac 
1236c7eeac06SToby Isaac   PetscFunctionBegin;
1237c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1238ef51cf95SToby Isaac   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup");
1239c7eeac06SToby Isaac   if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity);
1240c7eeac06SToby Isaac   forest->weightCapacity = capacity;
1241c7eeac06SToby Isaac   PetscFunctionReturn(0);
1242c7eeac06SToby Isaac }
1243c7eeac06SToby Isaac 
1244c7eeac06SToby Isaac #undef __FUNCT__
1245c7eeac06SToby Isaac #define __FUNCT__ "DMForestGetWeightCapacity"
1246*9be51f97SToby Isaac /*@
1247*9be51f97SToby Isaac   DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see
1248*9be51f97SToby Isaac   DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each process's cells to the
1249*9be51f97SToby Isaac   process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that the current process
1250*9be51f97SToby Isaac   should not have any cells after repartitioning.
1251*9be51f97SToby Isaac 
1252*9be51f97SToby Isaac   Not collective
1253*9be51f97SToby Isaac 
1254*9be51f97SToby Isaac   Input parameter:
1255*9be51f97SToby Isaac . dm - the forest
1256*9be51f97SToby Isaac 
1257*9be51f97SToby Isaac   Output parameter:
1258*9be51f97SToby Isaac . capacity - this process's capacity
1259*9be51f97SToby Isaac 
1260*9be51f97SToby Isaac   Level: advanced
1261*9be51f97SToby Isaac 
1262*9be51f97SToby Isaac .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
1263*9be51f97SToby Isaac @*/
1264c7eeac06SToby Isaac PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity)
1265c7eeac06SToby Isaac {
1266c7eeac06SToby Isaac   DM_Forest      *forest = (DM_Forest *) dm->data;
1267c7eeac06SToby Isaac 
1268c7eeac06SToby Isaac   PetscFunctionBegin;
1269c7eeac06SToby Isaac   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1270c7eeac06SToby Isaac   PetscValidRealPointer(capacity,2);
1271c7eeac06SToby Isaac   *capacity = forest->weightCapacity;
1272c7eeac06SToby Isaac   PetscFunctionReturn(0);
1273c7eeac06SToby Isaac }
1274c7eeac06SToby Isaac 
1275dd8e54a2SToby Isaac #undef __FUNCT__
1276db4d5e8cSToby Isaac #define __FUNCT__ "DMSetFromOptions_Forest"
12770709b2feSToby Isaac PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm)
1278db4d5e8cSToby Isaac {
1279db4d5e8cSToby Isaac   DM_Forest                  *forest = (DM_Forest *) dm->data;
128056ba9f64SToby Isaac   PetscBool                  flg, flg1, flg2, flg3, flg4;
1281dd8e54a2SToby Isaac   DMForestTopology           oldTopo;
1282c7eeac06SToby Isaac   char                       stringBuffer[256];
1283dd8e54a2SToby Isaac   PetscViewer                viewer;
1284dd8e54a2SToby Isaac   PetscViewerFormat          format;
128556ba9f64SToby Isaac   PetscInt                   adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade;
1286c7eeac06SToby Isaac   PetscReal                  weightsFactor;
1287c7eeac06SToby Isaac   DMForestAdaptivityStrategy adaptStrategy;
1288db4d5e8cSToby Isaac   PetscErrorCode             ierr;
1289db4d5e8cSToby Isaac 
1290db4d5e8cSToby Isaac   PetscFunctionBegin;
1291db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
129258762b62SToby Isaac   forest->setfromoptionscalled = PETSC_TRUE;
1293dd8e54a2SToby Isaac   ierr = DMForestGetTopology(dm, &oldTopo);CHKERRQ(ierr);
1294a3eda1eaSToby Isaac   ierr = PetscOptionsHead(PetscOptionsObject,"DMForest Options");CHKERRQ(ierr);
129556ba9f64SToby Isaac   ierr = PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);CHKERRQ(ierr);
129656ba9f64SToby Isaac   ierr = PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);CHKERRQ(ierr);
129756ba9f64SToby Isaac   ierr = PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);CHKERRQ(ierr);
129856ba9f64SToby Isaac   ierr = PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);CHKERRQ(ierr);
129956ba9f64SToby Isaac   if ((PetscInt) flg1 + (PetscInt) flg2 + (PetscInt) flg3 + (PetscInt) flg4 > 1) {
130056ba9f64SToby Isaac     SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Specify only one of -dm_forest_{topology,base_dm,coarse_forest,fine_forest}");
1301dd8e54a2SToby Isaac   }
130256ba9f64SToby Isaac   if (flg1) {
130356ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,(DMForestTopology)stringBuffer);CHKERRQ(ierr);
130456ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
130520e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
130656ba9f64SToby Isaac   }
130756ba9f64SToby Isaac   if (flg2) {
1308dd8e54a2SToby Isaac     DM         base;
1309dd8e54a2SToby Isaac 
1310dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&base);CHKERRQ(ierr);
1311dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1312dd8e54a2SToby Isaac     ierr = DMLoad(base,viewer);CHKERRQ(ierr);
1313dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1314dd8e54a2SToby Isaac     ierr = DMForestSetBaseDM(dm,base);CHKERRQ(ierr);
1315dd8e54a2SToby Isaac     ierr = DMDestroy(&base);CHKERRQ(ierr);
131656ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
131720e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
1318dd8e54a2SToby Isaac   }
131956ba9f64SToby Isaac   if (flg3) {
1320dd8e54a2SToby Isaac     DM         coarse;
1321dd8e54a2SToby Isaac 
1322dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&coarse);CHKERRQ(ierr);
1323dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1324dd8e54a2SToby Isaac     ierr = DMLoad(coarse,viewer);CHKERRQ(ierr);
1325dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
132620e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,coarse);CHKERRQ(ierr);
1327dd8e54a2SToby Isaac     ierr = DMDestroy(&coarse);CHKERRQ(ierr);
132856ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
132956ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1330dd8e54a2SToby Isaac   }
133156ba9f64SToby Isaac   if (flg4) {
1332dd8e54a2SToby Isaac     DM         fine;
1333dd8e54a2SToby Isaac 
1334dd8e54a2SToby Isaac     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&fine);CHKERRQ(ierr);
1335dd8e54a2SToby Isaac     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1336dd8e54a2SToby Isaac     ierr = DMLoad(fine,viewer);CHKERRQ(ierr);
1337dd8e54a2SToby Isaac     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
133820e8089bSToby Isaac     ierr = DMForestSetAdaptivityForest(dm,fine);CHKERRQ(ierr);
1339dd8e54a2SToby Isaac     ierr = DMDestroy(&fine);CHKERRQ(ierr);
134056ba9f64SToby Isaac     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
134156ba9f64SToby Isaac     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1342dd8e54a2SToby Isaac   }
1343dd8e54a2SToby Isaac   ierr = DMForestGetAdjacencyDimension(dm,&adjDim);CHKERRQ(ierr);
1344dd8e54a2SToby Isaac   ierr = PetscOptionsInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg);CHKERRQ(ierr);
1345dd8e54a2SToby Isaac   if (flg) {
1346dd8e54a2SToby Isaac     ierr = DMForestSetAdjacencyDimension(dm,adjDim);CHKERRQ(ierr);
1347dd8e54a2SToby Isaac   }
1348dd8e54a2SToby Isaac   else {
1349dd8e54a2SToby Isaac     ierr = DMForestGetAdjacencyCodimension(dm,&adjCodim);CHKERRQ(ierr);
1350dd8e54a2SToby Isaac     ierr = PetscOptionsInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg);CHKERRQ(ierr);
1351dd8e54a2SToby Isaac     if (flg) {
1352dd8e54a2SToby Isaac       ierr = DMForestSetAdjacencyCodimension(dm,adjCodim);CHKERRQ(ierr);
1353dd8e54a2SToby Isaac     }
1354dd8e54a2SToby Isaac   }
1355dd8e54a2SToby Isaac   ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr);
1356dd8e54a2SToby Isaac   ierr = PetscOptionsInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg);CHKERRQ(ierr);
1357dd8e54a2SToby Isaac   if (flg) {
1358dd8e54a2SToby Isaac     ierr = DMForestSetPartitionOverlap(dm,overlap);CHKERRQ(ierr);
1359dd8e54a2SToby Isaac   }
1360a6121fbdSMatthew G. Knepley #if 0
1361a6121fbdSMatthew 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);
1362a6121fbdSMatthew G. Knepley   if (flg) {
1363a6121fbdSMatthew G. Knepley     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1364a6121fbdSMatthew G. Knepley     ierr = DMForestSetInitialRefinement(dm,minRefinement);CHKERRQ(ierr);
1365a6121fbdSMatthew G. Knepley   }
1366a6121fbdSMatthew 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);
1367a6121fbdSMatthew G. Knepley   if (flg) {
1368a6121fbdSMatthew G. Knepley     ierr = DMForestSetMinimumRefinement(dm,0);CHKERRQ(ierr);
1369a6121fbdSMatthew G. Knepley     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
1370a6121fbdSMatthew G. Knepley   }
1371a6121fbdSMatthew G. Knepley #endif
1372dd8e54a2SToby Isaac   ierr = DMForestGetMinimumRefinement(dm,&minRefinement);CHKERRQ(ierr);
1373dd8e54a2SToby Isaac   ierr = PetscOptionsInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg);CHKERRQ(ierr);
1374dd8e54a2SToby Isaac   if (flg) {
1375dd8e54a2SToby Isaac     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1376db4d5e8cSToby Isaac   }
137756ba9f64SToby Isaac   ierr = DMForestGetInitialRefinement(dm,&initRefinement);CHKERRQ(ierr);
137856ba9f64SToby Isaac   ierr = PetscOptionsInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg);CHKERRQ(ierr);
137956ba9f64SToby Isaac   if (flg) {
138056ba9f64SToby Isaac     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
138156ba9f64SToby Isaac   }
1382c7eeac06SToby Isaac   ierr = DMForestGetMaximumRefinement(dm,&maxRefinement);CHKERRQ(ierr);
1383c7eeac06SToby Isaac   ierr = PetscOptionsInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg);CHKERRQ(ierr);
1384c7eeac06SToby Isaac   if (flg) {
1385c7eeac06SToby Isaac     ierr = DMForestSetMaximumRefinement(dm,maxRefinement);CHKERRQ(ierr);
1386c7eeac06SToby Isaac   }
1387c7eeac06SToby Isaac   ierr = DMForestGetAdaptivityStrategy(dm,&adaptStrategy);CHKERRQ(ierr);
1388c7eeac06SToby Isaac   ierr = PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,256,&flg);CHKERRQ(ierr);
1389c7eeac06SToby Isaac   if (flg) {
1390c7eeac06SToby Isaac     ierr = DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);CHKERRQ(ierr);
1391c7eeac06SToby Isaac   }
1392c7eeac06SToby Isaac   ierr = DMForestGetGradeFactor(dm,&grade);CHKERRQ(ierr);
1393c7eeac06SToby Isaac   ierr = PetscOptionsInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg);CHKERRQ(ierr);
1394c7eeac06SToby Isaac   if (flg) {
1395c7eeac06SToby Isaac     ierr = DMForestSetGradeFactor(dm,grade);CHKERRQ(ierr);
1396c7eeac06SToby Isaac   }
1397c7eeac06SToby Isaac   ierr = DMForestGetCellWeightFactor(dm,&weightsFactor);CHKERRQ(ierr);
1398c7eeac06SToby Isaac   ierr = PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);CHKERRQ(ierr);
1399c7eeac06SToby Isaac   if (flg) {
1400c7eeac06SToby Isaac     ierr = DMForestSetCellWeightFactor(dm,weightsFactor);CHKERRQ(ierr);
1401c7eeac06SToby Isaac   }
1402db4d5e8cSToby Isaac   ierr = PetscOptionsTail();CHKERRQ(ierr);
1403db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1404db4d5e8cSToby Isaac }
1405db4d5e8cSToby Isaac 
1406db4d5e8cSToby Isaac #undef __FUNCT__
1407d8984e3bSMatthew G. Knepley #define __FUNCT__ "DMCreateSubDM_Forest"
1408d8984e3bSMatthew G. Knepley PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
1409d8984e3bSMatthew G. Knepley {
1410d8984e3bSMatthew G. Knepley   PetscErrorCode ierr;
1411d8984e3bSMatthew G. Knepley 
1412d8984e3bSMatthew G. Knepley   PetscFunctionBegin;
1413d8984e3bSMatthew G. Knepley   if (subdm) {ierr = DMClone(dm, subdm);CHKERRQ(ierr);}
1414d8984e3bSMatthew G. Knepley   ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr);
1415d8984e3bSMatthew G. Knepley   PetscFunctionReturn(0);
1416d8984e3bSMatthew G. Knepley }
1417d8984e3bSMatthew G. Knepley 
1418d8984e3bSMatthew G. Knepley #undef __FUNCT__
14195421bac9SToby Isaac #define __FUNCT__ "DMRefine_Forest"
14205421bac9SToby Isaac PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
14215421bac9SToby Isaac {
14225421bac9SToby Isaac   DMLabel        refine;
14235421bac9SToby Isaac   DM             fineDM;
14245421bac9SToby Isaac   PetscErrorCode ierr;
14255421bac9SToby Isaac 
14265421bac9SToby Isaac   PetscFunctionBegin;
14275421bac9SToby Isaac   ierr = DMGetFineDM(dm,&fineDM);CHKERRQ(ierr);
14285421bac9SToby Isaac   if (fineDM) {
14295421bac9SToby Isaac     ierr = PetscObjectReference((PetscObject)fineDM);CHKERRQ(ierr);
14305421bac9SToby Isaac     *dmRefined = fineDM;
14315421bac9SToby Isaac     PetscFunctionReturn(0);
14325421bac9SToby Isaac   }
14335421bac9SToby Isaac   ierr = DMForestTemplate(dm,comm,dmRefined);CHKERRQ(ierr);
14345421bac9SToby Isaac   ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr);
14355421bac9SToby Isaac   if (!refine) {
14365421bac9SToby Isaac     ierr = DMCreateLabel(dm,"refine");CHKERRQ(ierr);
14375421bac9SToby Isaac     ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr);
14385421bac9SToby Isaac     ierr = DMLabelSetDefaultValue(refine,DM_FOREST_REFINE);CHKERRQ(ierr);
14395421bac9SToby Isaac   }
14405421bac9SToby Isaac   ierr = DMForestSetAdaptivityLabel(*dmRefined,"refine");CHKERRQ(ierr);
14415421bac9SToby Isaac   PetscFunctionReturn(0);
14425421bac9SToby Isaac }
14435421bac9SToby Isaac 
14445421bac9SToby Isaac #undef __FUNCT__
14455421bac9SToby Isaac #define __FUNCT__ "DMCoarsen_Forest"
14465421bac9SToby Isaac PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
14475421bac9SToby Isaac {
14485421bac9SToby Isaac   DMLabel        coarsen;
14495421bac9SToby Isaac   DM             coarseDM;
14505421bac9SToby Isaac   PetscErrorCode ierr;
14515421bac9SToby Isaac 
14525421bac9SToby Isaac   PetscFunctionBegin;
14535421bac9SToby Isaac   ierr = DMGetCoarseDM(dm,&coarseDM);CHKERRQ(ierr);
14545421bac9SToby Isaac   if (coarseDM) {
14555421bac9SToby Isaac     ierr = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr);
14565421bac9SToby Isaac     *dmCoarsened = coarseDM;
14575421bac9SToby Isaac     PetscFunctionReturn(0);
14585421bac9SToby Isaac   }
14595421bac9SToby Isaac   ierr = DMForestTemplate(dm,comm,dmCoarsened);CHKERRQ(ierr);
14605421bac9SToby Isaac   ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr);
14615421bac9SToby Isaac   if (!coarsen) {
14625421bac9SToby Isaac     ierr = DMCreateLabel(dm,"coarsen");CHKERRQ(ierr);
14635421bac9SToby Isaac     ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr);
14645421bac9SToby Isaac     ierr = DMLabelSetDefaultValue(coarsen,DM_FOREST_COARSEN);CHKERRQ(ierr);
14655421bac9SToby Isaac   }
14665421bac9SToby Isaac   ierr = DMForestSetAdaptivityLabel(*dmCoarsened,"coarsen");CHKERRQ(ierr);
14675421bac9SToby Isaac   PetscFunctionReturn(0);
14685421bac9SToby Isaac }
14695421bac9SToby Isaac 
14705421bac9SToby Isaac #undef __FUNCT__
1471d222f98bSToby Isaac #define __FUNCT__ "DMInitialize_Forest"
1472d222f98bSToby Isaac static PetscErrorCode DMInitialize_Forest(DM dm)
1473d222f98bSToby Isaac {
1474d222f98bSToby Isaac   PetscErrorCode ierr;
1475d222f98bSToby Isaac 
1476d222f98bSToby Isaac   PetscFunctionBegin;
1477d222f98bSToby Isaac   ierr = PetscMemzero(dm->ops,sizeof(*(dm->ops)));CHKERRQ(ierr);
1478d222f98bSToby Isaac 
1479d222f98bSToby Isaac   dm->ops->clone          = DMClone_Forest;
1480d222f98bSToby Isaac   dm->ops->setfromoptions = DMSetFromOptions_Forest;
1481d222f98bSToby Isaac   dm->ops->destroy        = DMDestroy_Forest;
1482d8984e3bSMatthew G. Knepley   dm->ops->createsubdm    = DMCreateSubDM_Forest;
14835421bac9SToby Isaac   dm->ops->refine         = DMRefine_Forest;
14845421bac9SToby Isaac   dm->ops->coarsen        = DMCoarsen_Forest;
1485d222f98bSToby Isaac   PetscFunctionReturn(0);
1486d222f98bSToby Isaac }
1487d222f98bSToby Isaac 
1488*9be51f97SToby Isaac /*MC
1489*9be51f97SToby Isaac   DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh.  Forests usually have a base DM (see DMForestGetBaseDM()), from which it is refined.  The refinement and partitioning of forests is considered immutable after DMSetUp() is called.  To adapt a mesh, one should call DMForestTemplate() to create a new mesh that will default to being identical to it, specify how that mesh should differ, and then calling DMSetUp() on the new mesh.
1490*9be51f97SToby Isaac 
1491*9be51f97SToby Isaac   To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the previous mesh whose values indicate which cells should be refined (DM_FOREST_REFINE) or coarsened (DM_FOREST_COARSEN) and how (subtypes are free to allow additional values for things like anisotropic refinement).  The name of the label should be given to the *new* mesh with DMForestSetAdaptivityLabel().
1492*9be51f97SToby Isaac 
1493*9be51f97SToby Isaac   Level: advanced
1494*9be51f97SToby Isaac 
1495*9be51f97SToby Isaac .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel()
1496*9be51f97SToby Isaac M*/
1497*9be51f97SToby Isaac 
1498d222f98bSToby Isaac #undef __FUNCT__
1499db4d5e8cSToby Isaac #define __FUNCT__ "DMCreate_Forest"
1500db4d5e8cSToby Isaac PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1501db4d5e8cSToby Isaac {
1502db4d5e8cSToby Isaac   DM_Forest      *forest;
1503db4d5e8cSToby Isaac   PetscErrorCode ierr;
1504db4d5e8cSToby Isaac 
1505db4d5e8cSToby Isaac   PetscFunctionBegin;
1506db4d5e8cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1507db4d5e8cSToby Isaac   ierr                        = PetscNewLog(dm,&forest);CHKERRQ(ierr);
1508db4d5e8cSToby Isaac   dm->dim                     = 0;
1509db4d5e8cSToby Isaac   dm->data                    = forest;
1510db4d5e8cSToby Isaac   forest->refct               = 1;
1511db4d5e8cSToby Isaac   forest->data                = NULL;
151258762b62SToby Isaac   forest->setfromoptionscalled      = PETSC_FALSE;
1513db4d5e8cSToby Isaac   forest->topology            = NULL;
1514db4d5e8cSToby Isaac   forest->base                = NULL;
1515db4d5e8cSToby Isaac   forest->adjDim              = PETSC_DEFAULT;
1516db4d5e8cSToby Isaac   forest->overlap             = PETSC_DEFAULT;
1517db4d5e8cSToby Isaac   forest->minRefinement       = PETSC_DEFAULT;
1518db4d5e8cSToby Isaac   forest->maxRefinement       = PETSC_DEFAULT;
151956ba9f64SToby Isaac   forest->initRefinement      = PETSC_DEFAULT;
1520c7eeac06SToby Isaac   forest->cStart              = PETSC_DETERMINE;
1521c7eeac06SToby Isaac   forest->cEnd                = PETSC_DETERMINE;
1522db4d5e8cSToby Isaac   forest->cellSF              = 0;
1523ebdf65a2SToby Isaac   forest->adaptLabel          = NULL;
1524db4d5e8cSToby Isaac   forest->gradeFactor         = 2;
1525db4d5e8cSToby Isaac   forest->cellWeights         = NULL;
1526db4d5e8cSToby Isaac   forest->cellWeightsCopyMode = PETSC_USE_POINTER;
1527db4d5e8cSToby Isaac   forest->weightsFactor       = 1.;
1528db4d5e8cSToby Isaac   forest->weightCapacity      = 1.;
1529a73e2921SToby Isaac   ierr = DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);CHKERRQ(ierr);
1530d222f98bSToby Isaac   ierr = DMInitialize_Forest(dm);CHKERRQ(ierr);
1531db4d5e8cSToby Isaac   PetscFunctionReturn(0);
1532db4d5e8cSToby Isaac }
1533db4d5e8cSToby Isaac 
1534