xref: /petsc/src/dm/impls/forest/forest.c (revision bff67a9b4439570f10deaeebe845aeaede4dac2a)
1 #include <petsc/private/dmforestimpl.h> /*I "petscdmforest.h" I*/
2 #include <petsc/private/dmimpl.h>       /*I "petscdm.h" I*/
3 #include <petscsf.h>
4 
5 PetscBool DMForestPackageInitialized = PETSC_FALSE;
6 
7 typedef struct _DMForestTypeLink*DMForestTypeLink;
8 
9 struct _DMForestTypeLink
10 {
11   char             *name;
12   DMForestTypeLink next;
13 };
14 
15 DMForestTypeLink DMForestTypeList;
16 
17 #undef __FUNCT__
18 #define __FUNCT__ "DMForestPackageFinalize"
19 static PetscErrorCode DMForestPackageFinalize(void)
20 {
21   DMForestTypeLink oldLink, link = DMForestTypeList;
22   PetscErrorCode   ierr;
23 
24   PetscFunctionBegin;
25   while (link) {
26     oldLink = link;
27     ierr    = PetscFree(oldLink->name);CHKERRQ(ierr);
28     link    = oldLink->next;
29     ierr    = PetscFree(oldLink);CHKERRQ(ierr);
30   }
31   PetscFunctionReturn(0);
32 }
33 
34 #undef __FUNCT__
35 #define __FUNCT__ "DMForestPackageInitialize"
36 static PetscErrorCode DMForestPackageInitialize(void)
37 {
38   PetscErrorCode ierr;
39 
40   PetscFunctionBegin;
41   if (DMForestPackageInitialized) PetscFunctionReturn(0);
42   DMForestPackageInitialized = PETSC_TRUE;
43 
44   ierr = DMForestRegisterType(DMFOREST);CHKERRQ(ierr);
45   ierr = PetscRegisterFinalize(DMForestPackageFinalize);CHKERRQ(ierr);
46   PetscFunctionReturn(0);
47 }
48 
49 #undef __FUNCT__
50 #define __FUNCT__ "DMForestRegisterType"
51 /*@C
52   DMForestRegisterType - Registers a DMType as a subtype of DMFOREST (so that DMIsForest() will be correct)
53 
54   Not Collective
55 
56   Input parameter:
57 . name - the name of the type
58 
59   Level: advanced
60 
61 .seealso: DMFOREST, DMIsForest()
62 @*/
63 PetscErrorCode DMForestRegisterType(DMType name)
64 {
65   DMForestTypeLink link;
66   PetscErrorCode   ierr;
67 
68   PetscFunctionBegin;
69   ierr             = DMForestPackageInitialize();CHKERRQ(ierr);
70   ierr             = PetscNew(&link);CHKERRQ(ierr);
71   ierr             = PetscStrallocpy(name,&link->name);CHKERRQ(ierr);
72   link->next       = DMForestTypeList;
73   DMForestTypeList = link;
74   PetscFunctionReturn(0);
75 }
76 
77 #undef __FUNCT__
78 #define __FUNCT__ "DMIsForest"
79 /*@
80   DMIsForest - Check whether a DM uses the DMFOREST interface for hierarchically-refined meshes
81 
82   Not Collective
83 
84   Input parameter:
85 . dm - the DM object
86 
87   Output parameter:
88 . isForest - whether dm is a subtype of DMFOREST
89 
90   Level: intermediate
91 
92 .seealso: DMFOREST, DMForestRegisterType()
93 @*/
94 PetscErrorCode DMIsForest(DM dm, PetscBool *isForest)
95 {
96   DMForestTypeLink link = DMForestTypeList;
97   PetscErrorCode   ierr;
98 
99   PetscFunctionBegin;
100   while (link) {
101     PetscBool sameType;
102     ierr = PetscObjectTypeCompare((PetscObject)dm,link->name,&sameType);CHKERRQ(ierr);
103     if (sameType) {
104       *isForest = PETSC_TRUE;
105       PetscFunctionReturn(0);
106     }
107     link = link->next;
108   }
109   *isForest = PETSC_FALSE;
110   PetscFunctionReturn(0);
111 }
112 
113 #undef __FUNCT__
114 #define __FUNCT__ "DMForestTemplate"
115 /*@
116   DMForestTemplate - Create a new DM that will be adapted from a source DM.  The new DM reproduces the configuration
117   of the source, but is not yet setup, so that the user can then define only the ways that the new DM should differ
118   (by, e.g., refinement or repartitioning).  The source DM is also set as the adaptivity source DM of the new DM (see
119   DMForestSetAdaptivityForest()).
120 
121   Collective on dm
122 
123   Input Parameters:
124 + dm - the source DM object
125 - comm - the communicator for the new DM (this communicator is currently ignored, but is present so that DMForestTemplate() can be used within DMCoarsen())
126 
127   Output Parameter:
128 . tdm - the new DM object
129 
130   Level: intermediate
131 
132 .seealso: DMForestSetAdaptivityForest()
133 @*/
134 PetscErrorCode DMForestTemplate(DM dm, MPI_Comm comm, DM *tdm)
135 {
136   DM_Forest                  *forest = (DM_Forest*) dm->data;
137   DMType                     type;
138   DM                         base;
139   DMForestTopology           topology;
140   PetscInt                   dim, overlap, ref, factor;
141   DMForestAdaptivityStrategy strat;
142   PetscDS                    ds;
143   void                       *ctx;
144   PetscErrorCode             (*map)(DM, PetscInt, PetscInt, const PetscReal[], PetscReal[], void*);
145   void                       *mapCtx;
146   PetscErrorCode             ierr;
147 
148   PetscFunctionBegin;
149   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
150   ierr = DMCreate(PetscObjectComm((PetscObject)dm),tdm);CHKERRQ(ierr);
151   ierr = DMGetType(dm,&type);CHKERRQ(ierr);
152   ierr = DMSetType(*tdm,type);CHKERRQ(ierr);
153   ierr = DMForestGetBaseDM(dm,&base);CHKERRQ(ierr);
154   ierr = DMForestSetBaseDM(*tdm,base);CHKERRQ(ierr);
155   ierr = DMForestGetTopology(dm,&topology);CHKERRQ(ierr);
156   ierr = DMForestSetTopology(*tdm,topology);CHKERRQ(ierr);
157   ierr = DMForestGetAdjacencyDimension(dm,&dim);CHKERRQ(ierr);
158   ierr = DMForestSetAdjacencyDimension(*tdm,dim);CHKERRQ(ierr);
159   ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr);
160   ierr = DMForestSetPartitionOverlap(*tdm,overlap);CHKERRQ(ierr);
161   ierr = DMForestGetMinimumRefinement(dm,&ref);CHKERRQ(ierr);
162   ierr = DMForestSetMinimumRefinement(*tdm,ref);CHKERRQ(ierr);
163   ierr = DMForestGetMaximumRefinement(dm,&ref);CHKERRQ(ierr);
164   ierr = DMForestSetMaximumRefinement(*tdm,ref);CHKERRQ(ierr);
165   ierr = DMForestGetAdaptivityStrategy(dm,&strat);CHKERRQ(ierr);
166   ierr = DMForestSetAdaptivityStrategy(*tdm,strat);CHKERRQ(ierr);
167   ierr = DMForestGetGradeFactor(dm,&factor);CHKERRQ(ierr);
168   ierr = DMForestSetGradeFactor(*tdm,factor);CHKERRQ(ierr);
169   ierr = DMForestGetBaseCoordinateMapping(dm,&map,&mapCtx);CHKERRQ(ierr);
170   ierr = DMForestSetBaseCoordinateMapping(*tdm,map,mapCtx);CHKERRQ(ierr);
171   if (forest->ftemplate) {
172     ierr = (forest->ftemplate)(dm, *tdm);CHKERRQ(ierr);
173   }
174   ierr = DMForestSetAdaptivityForest(*tdm,dm);CHKERRQ(ierr);
175   ierr = DMGetDS(dm,&ds);CHKERRQ(ierr);
176   ierr = DMSetDS(*tdm,ds);CHKERRQ(ierr);
177   ierr = DMGetApplicationContext(dm,&ctx);CHKERRQ(ierr);
178   ierr = DMSetApplicationContext(*tdm,&ctx);CHKERRQ(ierr);
179   if (dm->maxCell) {
180     const PetscReal      *maxCell, *L;
181     const DMBoundaryType *bd;
182 
183     ierr = DMGetPeriodicity(dm,&maxCell,&L,&bd);CHKERRQ(ierr);
184     ierr = DMSetPeriodicity(*tdm,maxCell,L,bd);CHKERRQ(ierr);
185   }
186   ierr = DMCopyBoundary(dm,*tdm);CHKERRQ(ierr);
187   PetscFunctionReturn(0);
188 }
189 
190 static PetscErrorCode DMInitialize_Forest(DM dm);
191 
192 #undef __FUNCT__
193 #define __FUNCT__ "DMClone_Forest"
194 PETSC_EXTERN PetscErrorCode DMClone_Forest(DM dm, DM *newdm)
195 {
196   DM_Forest      *forest = (DM_Forest*) dm->data;
197   const char     *type;
198   PetscErrorCode ierr;
199 
200   PetscFunctionBegin;
201   forest->refct++;
202   (*newdm)->data = forest;
203   ierr           = PetscObjectGetType((PetscObject) dm, &type);CHKERRQ(ierr);
204   ierr           = PetscObjectChangeTypeName((PetscObject) *newdm, type);CHKERRQ(ierr);
205   ierr           = DMInitialize_Forest(*newdm);CHKERRQ(ierr);
206   PetscFunctionReturn(0);
207 }
208 
209 #undef __FUNCT__
210 #define __FUNCT__ "DMDestroy_Forest"
211 static PetscErrorCode DMDestroy_Forest(DM dm)
212 {
213   DM_Forest      *forest = (DM_Forest*) dm->data;
214   PetscErrorCode ierr;
215 
216   PetscFunctionBegin;
217   if (--forest->refct > 0) PetscFunctionReturn(0);
218   if (forest->destroy) {ierr = forest->destroy(dm);CHKERRQ(ierr);}
219   ierr = PetscSFDestroy(&forest->cellSF);CHKERRQ(ierr);
220   ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr);
221   ierr = PetscSFDestroy(&forest->coarseToPreFine);CHKERRQ(ierr);
222   ierr = PetscFree(forest->adaptLabel);CHKERRQ(ierr);
223   ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr);
224   ierr = DMDestroy(&forest->base);CHKERRQ(ierr);
225   ierr = DMDestroy(&forest->adapt);CHKERRQ(ierr);
226   ierr = PetscFree(forest->topology);CHKERRQ(ierr);
227   ierr = PetscFree(forest);CHKERRQ(ierr);
228   PetscFunctionReturn(0);
229 }
230 
231 #undef __FUNCT__
232 #define __FUNCT__ "DMForestSetTopology"
233 /*@C
234   DMForestSetTopology - Set the topology of a DMForest during the pre-setup phase.  The topology is a string (e.g.
235   "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest durint
236   DMSetUp().
237 
238   Logically collective on dm
239 
240   Input parameters:
241 + dm - the forest
242 - topology - the topology of the forest
243 
244   Level: intermediate
245 
246 .seealso(): DMForestGetTopology(), DMForestSetBaseDM()
247 @*/
248 PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology)
249 {
250   DM_Forest      *forest = (DM_Forest*) dm->data;
251   PetscErrorCode ierr;
252 
253   PetscFunctionBegin;
254   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
255   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the topology after setup");
256   ierr = PetscFree(forest->topology);CHKERRQ(ierr);
257   ierr = PetscStrallocpy((const char*)topology,(char**) &forest->topology);CHKERRQ(ierr);
258   PetscFunctionReturn(0);
259 }
260 
261 #undef __FUNCT__
262 #define __FUNCT__ "DMForestGetTopology"
263 /*@C
264   DMForestGetTopology - Get a string describing the topology of a DMForest.
265 
266   Not collective
267 
268   Input parameter:
269 . dm - the forest
270 
271   Output parameter:
272 . topology - the topology of the forest (e.g., 'cube', 'shell')
273 
274   Level: intermediate
275 
276 .seealso: DMForestSetTopology()
277 @*/
278 PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology)
279 {
280   DM_Forest *forest = (DM_Forest*) dm->data;
281 
282   PetscFunctionBegin;
283   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
284   PetscValidPointer(topology,2);
285   *topology = forest->topology;
286   PetscFunctionReturn(0);
287 }
288 
289 #undef __FUNCT__
290 #define __FUNCT__ "DMForestSetBaseDM"
291 /*@
292   DMForestSetBaseDM - During the pre-setup phase, set the DM that defines the base mesh of a DMForest forest.  The
293   forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its
294   base.  In general, two forest must share a bse to be comparable, to do things like construct interpolators.
295 
296   Logically collective on dm
297 
298   Input Parameters:
299 + dm - the forest
300 - base - the base DM of the forest
301 
302   Level: intermediate
303 
304 .seealso(): DMForestGetBaseDM()
305 @*/
306 PetscErrorCode DMForestSetBaseDM(DM dm, DM base)
307 {
308   DM_Forest      *forest = (DM_Forest*) dm->data;
309   PetscInt       dim, dimEmbed;
310   PetscErrorCode ierr;
311 
312   PetscFunctionBegin;
313   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
314   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the base after setup");
315   ierr         = PetscObjectReference((PetscObject)base);CHKERRQ(ierr);
316   ierr         = DMDestroy(&forest->base);CHKERRQ(ierr);
317   forest->base = base;
318   if (base) {
319     PetscValidHeaderSpecific(base, DM_CLASSID, 2);
320     ierr = DMGetDimension(base,&dim);CHKERRQ(ierr);
321     ierr = DMSetDimension(dm,dim);CHKERRQ(ierr);
322     ierr = DMGetCoordinateDim(base,&dimEmbed);CHKERRQ(ierr);
323     ierr = DMSetCoordinateDim(dm,dimEmbed);CHKERRQ(ierr);
324   }
325   PetscFunctionReturn(0);
326 }
327 
328 #undef __FUNCT__
329 #define __FUNCT__ "DMForestGetBaseDM"
330 /*@
331   DMForestGetBaseDM - Get the base DM of a DMForest forest.  The forest will be hierarchically refined from the base,
332   and all refinements/coarsenings of the forest will share its base.  In general, two forest must share a bse to be
333   comparable, to do things like construct interpolators.
334 
335   Not collective
336 
337   Input Parameter:
338 . dm - the forest
339 
340   Output Parameter:
341 . base - the base DM of the forest
342 
343   Level: intermediate
344 
345 .seealso(); DMForestSetBaseDM()
346 @*/
347 PetscErrorCode DMForestGetBaseDM(DM dm, DM *base)
348 {
349   DM_Forest *forest = (DM_Forest*) dm->data;
350 
351   PetscFunctionBegin;
352   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
353   PetscValidPointer(base, 2);
354   *base = forest->base;
355   PetscFunctionReturn(0);
356 }
357 
358 #undef __FUNCT__
359 #define __FUNCT__ "DMForestSetBaseCoordinateMapping"
360 PetscErrorCode DMForestSetBaseCoordinateMapping(DM dm, PetscErrorCode (*func)(DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx)
361 {
362   DM_Forest *forest = (DM_Forest*) dm->data;
363 
364   PetscFunctionBegin;
365   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
366   forest->mapcoordinates    = func;
367   forest->mapcoordinatesctx = ctx;
368   PetscFunctionReturn(0);
369 }
370 
371 #undef __FUNCT__
372 #define __FUNCT__ "DMForestGetBaseCoordinateMapping"
373 PetscErrorCode DMForestGetBaseCoordinateMapping(DM dm, PetscErrorCode (**func) (DM,PetscInt,PetscInt,const PetscReal [],PetscReal [],void*),void *ctx)
374 {
375   DM_Forest *forest = (DM_Forest*) dm->data;
376 
377   PetscFunctionBegin;
378   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
379   if (func) *func = forest->mapcoordinates;
380   if (ctx) *((void**) ctx) = forest->mapcoordinatesctx;
381   PetscFunctionReturn(0);
382 }
383 
384 #undef __FUNCT__
385 #define __FUNCT__ "DMForestSetAdaptivityForest"
386 /*@
387   DMForestSetAdaptivityForest - During the pre-setup phase, set the forest from which the current forest will be
388   adapted (e.g., the current forest will be refined/coarsened/repartitioned from it) im DMSetUp().  Usually not needed
389   by users directly: DMForestTemplate() constructs a new forest to be adapted from an old forest and calls this
390   routine.
391 
392   Note that this can be called after setup with adapt = NULL, which will clear all internal data related to the
393   adaptivity forest from dm.  This way, repeatedly adapting does not leave stale DM objects in memory.
394 
395   Logically collective on dm
396 
397   Input Parameter:
398 + dm - the new forest, which will be constructed from adapt
399 - adapt - the old forest
400 
401   Level: intermediate
402 
403 .seealso: DMForestGetAdaptivityForest(), DMForestSetAdaptivityPurpose()
404 @*/
405 PetscErrorCode DMForestSetAdaptivityForest(DM dm,DM adapt)
406 {
407   DM_Forest      *forest, *adaptForest, *oldAdaptForest;
408   DM             oldAdapt;
409   PetscErrorCode ierr;
410 
411   PetscFunctionBegin;
412   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
413   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
414   forest   = (DM_Forest*) dm->data;
415   ierr     = DMForestGetAdaptivityForest(dm,&oldAdapt);CHKERRQ(ierr);
416   if (adapt != NULL && dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
417   adaptForest    = (DM_Forest*) adapt ? adapt->data : NULL;
418   oldAdaptForest = (DM_Forest*) oldAdapt ? oldAdapt->data : NULL;
419   if (adaptForest != oldAdaptForest) {
420     ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr);
421     ierr = PetscSFDestroy(&forest->coarseToPreFine);CHKERRQ(ierr);
422     if (forest->clearadaptivityforest) {ierr = (forest->clearadaptivityforest)(dm);CHKERRQ(ierr);}
423   }
424   switch (forest->adaptPurpose) {
425   case DM_FOREST_KEEP:
426     ierr          = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr);
427     ierr          = DMDestroy(&(forest->adapt));CHKERRQ(ierr);
428     forest->adapt = adapt;
429     break;
430   case DM_FOREST_REFINE:
431     ierr = DMSetCoarseDM(dm,adapt);CHKERRQ(ierr);
432     break;
433   case DM_FOREST_COARSEN:
434     ierr = DMSetFineDM(dm,adapt);CHKERRQ(ierr);
435     break;
436   default:
437     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
438   }
439   PetscFunctionReturn(0);
440 }
441 
442 #undef __FUNCT__
443 #define __FUNCT__ "DMForestGetAdaptivityForest"
444 /*@
445   DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted.
446 
447   Not collective
448 
449   Input Parameter:
450 . dm - the forest
451 
452   Output Parameter:
453 . adapt - the forest from which dm is/was adapted
454 
455   Level: intermediate
456 
457 .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose()
458 @*/
459 PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt)
460 {
461   DM_Forest      *forest;
462   PetscErrorCode ierr;
463 
464   PetscFunctionBegin;
465   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
466   forest = (DM_Forest*) dm->data;
467   switch (forest->adaptPurpose) {
468   case DM_FOREST_KEEP:
469     *adapt = forest->adapt;
470     break;
471   case DM_FOREST_REFINE:
472     ierr = DMGetCoarseDM(dm,adapt);CHKERRQ(ierr);
473     break;
474   case DM_FOREST_COARSEN:
475     ierr = DMGetFineDM(dm,adapt);CHKERRQ(ierr);
476     break;
477   default:
478     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
479   }
480   PetscFunctionReturn(0);
481 }
482 
483 #undef __FUNCT__
484 #define __FUNCT__ "DMForestSetAdaptivityPurpose"
485 /*@
486   DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its
487   source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_FOREST_REFINE), coarsening
488   (DM_FOREST_COARSEN), or undefined (DM_FOREST_NONE).  This only matters for the purposes of reference counting:
489   during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse
490   relationship (see DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does
491   not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can
492   cause a memory leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
493 
494   Logically collective on dm
495 
496   Input Parameters:
497 + dm - the forest
498 - purpose - the adaptivity purpose (DM_FOREST_NONE/DM_FOREST_REFINE/DM_FOREST_COARSEN)
499 
500   Level: advanced
501 
502 .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest()
503 @*/
504 PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMForestAdaptivityPurpose purpose)
505 {
506   DM_Forest      *forest;
507   PetscErrorCode ierr;
508 
509   PetscFunctionBegin;
510   forest = (DM_Forest*) dm->data;
511   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
512   if (purpose != forest->adaptPurpose) {
513     DM adapt;
514 
515     ierr = DMForestGetAdaptivityForest(dm,&adapt);CHKERRQ(ierr);
516     ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr);
517     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
518 
519     forest->adaptPurpose = purpose;
520 
521     ierr = DMForestSetAdaptivityForest(dm,adapt);CHKERRQ(ierr);
522     ierr = DMDestroy(&adapt);CHKERRQ(ierr);
523   }
524   PetscFunctionReturn(0);
525 }
526 
527 #undef __FUNCT__
528 #define __FUNCT__ "DMForestGetAdaptivityPurpose"
529 /*@
530   DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with
531   DMForestSetAdaptivityForest()) for the purpose of refinement (DM_FOREST_REFINE), coarsening (DM_FOREST_COARSEN), or
532   undefined (DM_FOREST_NONE).  This only matters for the purposes of reference counting: during DMDestroy(), cyclic
533   references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see
534   DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does not maintain a
535   reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory
536   leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
537 
538   Not collective
539 
540   Input Parameter:
541 . dm - the forest
542 
543   Output Parameter:
544 . purpose - the adaptivity purpose (DM_FOREST_NONE/DM_FOREST_REFINE/DM_FOREST_COARSEN)
545 
546   Level: advanced
547 
548 .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest()
549 @*/
550 PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMForestAdaptivityPurpose *purpose)
551 {
552   DM_Forest *forest;
553 
554   PetscFunctionBegin;
555   forest   = (DM_Forest*) dm->data;
556   *purpose = forest->adaptPurpose;
557   PetscFunctionReturn(0);
558 }
559 
560 #undef __FUNCT__
561 #define __FUNCT__ "DMForestSetAdjacencyDimension"
562 /*@
563   DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine
564   cell adjacency (for the purposes of partitioning and overlap).
565 
566   Logically collective on dm
567 
568   Input Parameters:
569 + dm - the forest
570 - adjDim - default 0 (i.e., vertices determine adjacency)
571 
572   Level: intermediate
573 
574 .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap()
575 @*/
576 PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim)
577 {
578   PetscInt       dim;
579   DM_Forest      *forest = (DM_Forest*) dm->data;
580   PetscErrorCode ierr;
581 
582   PetscFunctionBegin;
583   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
584   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup");
585   if (adjDim < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %d", adjDim);
586   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
587   if (adjDim > dim) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %d: %d", dim, adjDim);
588   forest->adjDim = adjDim;
589   PetscFunctionReturn(0);
590 }
591 
592 #undef __FUNCT__
593 #define __FUNCT__ "DMForestSetAdjacencyCodimension"
594 /*@
595   DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that,
596   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
597 
598   Logically collective on dm
599 
600   Input Parameters:
601 + dm - the forest
602 - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
603 
604   Level: intermediate
605 
606 .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension()
607 @*/
608 PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim)
609 {
610   PetscInt       dim;
611   PetscErrorCode ierr;
612 
613   PetscFunctionBegin;
614   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
615   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
616   ierr = DMForestSetAdjacencyDimension(dm,dim-adjCodim);CHKERRQ(ierr);
617   PetscFunctionReturn(0);
618 }
619 
620 #undef __FUNCT__
621 #define __FUNCT__ "DMForestGetAdjacencyDimension"
622 /*@
623   DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the
624   purposes of partitioning and overlap).
625 
626   Not collective
627 
628   Input Parameter:
629 . dm - the forest
630 
631   Output Parameter:
632 . adjDim - default 0 (i.e., vertices determine adjacency)
633 
634   Level: intermediate
635 
636 .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap()
637 @*/
638 PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim)
639 {
640   DM_Forest *forest = (DM_Forest*) dm->data;
641 
642   PetscFunctionBegin;
643   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
644   PetscValidIntPointer(adjDim,2);
645   *adjDim = forest->adjDim;
646   PetscFunctionReturn(0);
647 }
648 
649 #undef __FUNCT__
650 #define __FUNCT__ "DMForestGetAdjacencyCodimension"
651 /*@
652   DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that,
653   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
654 
655   Not collective
656 
657   Input Parameter:
658 . dm - the forest
659 
660   Output Parameter:
661 . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
662 
663   Level: intermediate
664 
665 .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension()
666 @*/
667 PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim)
668 {
669   DM_Forest      *forest = (DM_Forest*) dm->data;
670   PetscInt       dim;
671   PetscErrorCode ierr;
672 
673   PetscFunctionBegin;
674   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
675   PetscValidIntPointer(adjCodim,2);
676   ierr      = DMGetDimension(dm,&dim);CHKERRQ(ierr);
677   *adjCodim = dim - forest->adjDim;
678   PetscFunctionReturn(0);
679 }
680 
681 #undef __FUNCT__
682 #define __FUNCT__ "DMForestSetPartitionOverlap"
683 /*@
684   DMForestSetPartitionOverlap - During the pre-setup phase, set the amount of cell-overlap present in parallel
685   partitions of a forest, with values > 0 indicating subdomains that are expanded by that many iterations of adding
686   adjacent cells
687 
688   Logically collective on dm
689 
690   Input Parameters:
691 + dm - the forest
692 - overlap - default 0
693 
694   Level: intermediate
695 
696 .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
697 @*/
698 PetscErrorCode DMForestSetPartitionOverlap(DM dm, PetscInt overlap)
699 {
700   DM_Forest *forest = (DM_Forest*) dm->data;
701 
702   PetscFunctionBegin;
703   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
704   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup");
705   if (overlap < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %d", overlap);
706   forest->overlap = overlap;
707   PetscFunctionReturn(0);
708 }
709 
710 #undef __FUNCT__
711 #define __FUNCT__ "DMForestGetPartitionOverlap"
712 /*@
713   DMForestGetPartitionOverlap - Get the amount of cell-overlap present in parallel partitions of a forest, with values
714   > 0 indicating subdomains that are expanded by that many iterations of adding adjacent cells
715 
716   Not collective
717 
718   Input Parameter:
719 . dm - the forest
720 
721   Output Parameter:
722 . overlap - default 0
723 
724   Level: intermediate
725 
726 .seealso: DMForestGetPartitionOverlap(), DMForestSetAdjacencyDimension(), DMForestSetAdjacencyCodimension()
727 @*/
728 PetscErrorCode DMForestGetPartitionOverlap(DM dm, PetscInt *overlap)
729 {
730   DM_Forest *forest = (DM_Forest*) dm->data;
731 
732   PetscFunctionBegin;
733   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
734   PetscValidIntPointer(overlap,2);
735   *overlap = forest->overlap;
736   PetscFunctionReturn(0);
737 }
738 
739 #undef __FUNCT__
740 #define __FUNCT__ "DMForestSetMinimumRefinement"
741 /*@
742   DMForestSetMinimumRefinement - During the pre-setup phase, set the minimum level of refinement (relative to the base
743   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest
744   (see DMForestGetAdaptivityForest()) this limits the amount of coarsening.
745 
746   Logically collective on dm
747 
748   Input Parameters:
749 + dm - the forest
750 - minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
751 
752   Level: intermediate
753 
754 .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
755 @*/
756 PetscErrorCode DMForestSetMinimumRefinement(DM dm, PetscInt minRefinement)
757 {
758   DM_Forest *forest = (DM_Forest*) dm->data;
759 
760   PetscFunctionBegin;
761   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
762   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup");
763   forest->minRefinement = minRefinement;
764   PetscFunctionReturn(0);
765 }
766 
767 #undef __FUNCT__
768 #define __FUNCT__ "DMForestGetMinimumRefinement"
769 /*@
770   DMForestGetMinimumRefinement - Get the minimum level of refinement (relative to the base DM, see
771   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by coarsening a previous forest (see
772   DMForestGetAdaptivityForest()), this limits the amount of coarsening.
773 
774   Not collective
775 
776   Input Parameter:
777 . dm - the forest
778 
779   Output Parameter:
780 . minRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
781 
782   Level: intermediate
783 
784 .seealso: DMForestSetMinimumRefinement(), DMForestGetMaximumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
785 @*/
786 PetscErrorCode DMForestGetMinimumRefinement(DM dm, PetscInt *minRefinement)
787 {
788   DM_Forest *forest = (DM_Forest*) dm->data;
789 
790   PetscFunctionBegin;
791   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
792   PetscValidIntPointer(minRefinement,2);
793   *minRefinement = forest->minRefinement;
794   PetscFunctionReturn(0);
795 }
796 
797 #undef __FUNCT__
798 #define __FUNCT__ "DMForestSetInitialRefinement"
799 /*@
800   DMForestSetInitialRefinement - During the pre-setup phase, set the initial level of refinement (relative to the base
801   DM, see DMForestGetBaseDM()) allowed in the forest.
802 
803   Logically collective on dm
804 
805   Input Parameters:
806 + dm - the forest
807 - initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
808 
809   Level: intermediate
810 
811 .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
812 @*/
813 PetscErrorCode DMForestSetInitialRefinement(DM dm, PetscInt initRefinement)
814 {
815   DM_Forest *forest = (DM_Forest*) dm->data;
816 
817   PetscFunctionBegin;
818   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
819   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the initial refinement after setup");
820   forest->initRefinement = initRefinement;
821   PetscFunctionReturn(0);
822 }
823 
824 #undef __FUNCT__
825 #define __FUNCT__ "DMForestGetInitialRefinement"
826 /*@
827   DMForestGetInitialRefinement - Get the initial level of refinement (relative to the base DM, see
828   DMForestGetBaseDM()) allowed in the forest.
829 
830   Not collective
831 
832   Input Parameter:
833 . dm - the forest
834 
835   Output Paramater:
836 . initefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
837 
838   Level: intermediate
839 
840 .seealso: DMForestSetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestGetBaseDM()
841 @*/
842 PetscErrorCode DMForestGetInitialRefinement(DM dm, PetscInt *initRefinement)
843 {
844   DM_Forest *forest = (DM_Forest*) dm->data;
845 
846   PetscFunctionBegin;
847   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
848   PetscValidIntPointer(initRefinement,2);
849   *initRefinement = forest->initRefinement;
850   PetscFunctionReturn(0);
851 }
852 
853 #undef __FUNCT__
854 #define __FUNCT__ "DMForestSetMaximumRefinement"
855 /*@
856   DMForestSetMaximumRefinement - During the pre-setup phase, set the maximum level of refinement (relative to the base
857   DM, see DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest
858   (see DMForestGetAdaptivityForest()), this limits the amount of refinement.
859 
860   Logically collective on dm
861 
862   Input Parameters:
863 + dm - the forest
864 - maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
865 
866   Level: intermediate
867 
868 .seealso: DMForestGetMinimumRefinement(), DMForestSetMaximumRefinement(), DMForestSetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityDM()
869 @*/
870 PetscErrorCode DMForestSetMaximumRefinement(DM dm, PetscInt maxRefinement)
871 {
872   DM_Forest *forest = (DM_Forest*) dm->data;
873 
874   PetscFunctionBegin;
875   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
876   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the maximum refinement after setup");
877   forest->maxRefinement = maxRefinement;
878   PetscFunctionReturn(0);
879 }
880 
881 #undef __FUNCT__
882 #define __FUNCT__ "DMForestGetMaximumRefinement"
883 /*@
884   DMForestGetMaximumRefinement - Get the maximum level of refinement (relative to the base DM, see
885   DMForestGetBaseDM()) allowed in the forest.  If the forest is being created by refining a previous forest (see
886   DMForestGetAdaptivityForest()), this limits the amount of refinement.
887 
888   Not collective
889 
890   Input Parameter:
891 . dm - the forest
892 
893   Output Parameter:
894 . maxRefinement - default PETSC_DEFAULT (interpreted by the subtype of DMForest)
895 
896   Level: intermediate
897 
898 .seealso: DMForestSetMaximumRefinement(), DMForestGetMinimumRefinement(), DMForestGetInitialRefinement(), DMForestGetBaseDM(), DMForestGetAdaptivityForest()
899 @*/
900 PetscErrorCode DMForestGetMaximumRefinement(DM dm, PetscInt *maxRefinement)
901 {
902   DM_Forest *forest = (DM_Forest*) dm->data;
903 
904   PetscFunctionBegin;
905   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
906   PetscValidIntPointer(maxRefinement,2);
907   *maxRefinement = forest->maxRefinement;
908   PetscFunctionReturn(0);
909 }
910 
911 #undef __FUNCT__
912 #define __FUNCT__ "DMForestSetAdaptivityStrategy"
913 /*@C
914   DMForestSetAdaptivityStrategy - During the pre-setup phase, set the strategy for combining adaptivity labels from multiple processes.
915   Subtypes of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all processes must agree
916   for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only one process needs to
917   specify refinement/coarsening.
918 
919   Logically collective on dm
920 
921   Input Parameters:
922 + dm - the forest
923 - adaptStrategy - default DMFORESTADAPTALL
924 
925   Level: advanced
926 
927 .seealso: DMForestGetAdaptivityStrategy()
928 @*/
929 PetscErrorCode DMForestSetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy adaptStrategy)
930 {
931   DM_Forest      *forest = (DM_Forest*) dm->data;
932   PetscErrorCode ierr;
933 
934   PetscFunctionBegin;
935   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
936   ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr);
937   ierr = PetscStrallocpy((const char*) adaptStrategy,(char**)&forest->adaptStrategy);CHKERRQ(ierr);
938   PetscFunctionReturn(0);
939 }
940 
941 #undef __FUNCT__
942 #define __FUNCT__ "DMForestGetAdaptivityStrategy"
943 /*@C
944   DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes.  Subtypes
945   of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all
946   processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only
947   one process needs to specify refinement/coarsening.
948 
949   Not collective
950 
951   Input Parameter:
952 . dm - the forest
953 
954   Output Parameter:
955 . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL)
956 
957   Level: advanced
958 
959 .seealso: DMForestSetAdaptivityStrategy()
960 @*/
961 PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy)
962 {
963   DM_Forest *forest = (DM_Forest*) dm->data;
964 
965   PetscFunctionBegin;
966   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
967   PetscValidPointer(adaptStrategy,2);
968   *adaptStrategy = forest->adaptStrategy;
969   PetscFunctionReturn(0);
970 }
971 
972 #undef __FUNCT__
973 #define __FUNCT__ "DMForestSetComputeAdaptivitySF"
974 /*@
975   DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed
976   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().
977 
978   Logically collective on dm
979 
980   Input Parameters:
981 + dm - the post-adaptation forest
982 - computeSF - default PETSC_TRUE
983 
984   Level: advanced
985 
986 .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
987 @*/
988 PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF)
989 {
990   DM_Forest *forest;
991 
992   PetscFunctionBegin;
993   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
994   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called");
995   forest                 = (DM_Forest*) dm->data;
996   forest->computeAdaptSF = computeSF;
997   PetscFunctionReturn(0);
998 }
999 
1000 #undef __FUNCT__
1001 #define __FUNCT__ "DMForestTransferVec"
1002 PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time)
1003 {
1004   DM_Forest      *forest;
1005   PetscErrorCode ierr;
1006 
1007   PetscFunctionBegin;
1008   PetscValidHeaderSpecific(dmIn   ,DM_CLASSID  ,1);
1009   PetscValidHeaderSpecific(vecIn  ,VEC_CLASSID ,2);
1010   PetscValidHeaderSpecific(dmOut  ,DM_CLASSID  ,3);
1011   PetscValidHeaderSpecific(vecOut ,VEC_CLASSID ,4);
1012   forest = (DM_Forest *) dmIn->data;
1013   if (!forest->transfervec) SETERRQ(PetscObjectComm((PetscObject)dmIn),PETSC_ERR_SUP,"DMForestTransferVec() not implemented");
1014   ierr = (forest->transfervec)(dmIn,vecIn,dmOut,vecOut,useBCs,time);CHKERRQ(ierr);
1015   PetscFunctionReturn(0);
1016 }
1017 
1018 #undef __FUNCT__
1019 #define __FUNCT__ "DMForestGetComputeAdaptivitySF"
1020 /*@
1021   DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the
1022   pre-adaptation forest to the post-adaptiation forest.  After DMSetUp() is called, these transfer PetscSFs can be
1023   accessed with DMForestGetAdaptivitySF().
1024 
1025   Not collective
1026 
1027   Input Parameter:
1028 . dm - the post-adaptation forest
1029 
1030   Output Parameter:
1031 . computeSF - default PETSC_TRUE
1032 
1033   Level: advanced
1034 
1035 .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
1036 @*/
1037 PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF)
1038 {
1039   DM_Forest *forest;
1040 
1041   PetscFunctionBegin;
1042   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1043   forest     = (DM_Forest*) dm->data;
1044   *computeSF = forest->computeAdaptSF;
1045   PetscFunctionReturn(0);
1046 }
1047 
1048 #undef __FUNCT__
1049 #define __FUNCT__ "DMForestGetAdaptivitySF"
1050 /*@
1051   DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest.
1052   Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be
1053   some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa.  Therefore there are two
1054   PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates
1055   pre-adaptation fine cells to post-adaptation coarse cells.
1056 
1057   Not collective
1058 
1059   Input Parameter:
1060   dm - the post-adaptation forest
1061 
1062   Output Parameter:
1063   preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post-
1064   coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre-
1065 
1066   Level: advanced
1067 
1068 .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF()
1069 @*/
1070 PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine)
1071 {
1072   DM_Forest      *forest;
1073   PetscErrorCode ierr;
1074 
1075   PetscFunctionBegin;
1076   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1077   ierr   = DMSetUp(dm);CHKERRQ(ierr);
1078   forest = (DM_Forest*) dm->data;
1079   if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine;
1080   if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine;
1081   PetscFunctionReturn(0);
1082 }
1083 
1084 #undef __FUNCT__
1085 #define __FUNCT__ "DMForestSetGradeFactor"
1086 /*@
1087   DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to
1088   indicate that the diameter of neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may
1089   only support one particular choice of grading factor.
1090 
1091   Logically collective on dm
1092 
1093   Input Parameters:
1094 + dm - the forest
1095 - grade - the grading factor
1096 
1097   Level: advanced
1098 
1099 .seealso: DMForestGetGradeFactor()
1100 @*/
1101 PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade)
1102 {
1103   DM_Forest *forest = (DM_Forest*) dm->data;
1104 
1105   PetscFunctionBegin;
1106   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1107   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup");
1108   forest->gradeFactor = grade;
1109   PetscFunctionReturn(0);
1110 }
1111 
1112 #undef __FUNCT__
1113 #define __FUNCT__ "DMForestGetGradeFactor"
1114 /*@
1115   DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of
1116   neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may only support one particular
1117   choice of grading factor.
1118 
1119   Not collective
1120 
1121   Input Parameter:
1122 . dm - the forest
1123 
1124   Output Parameter:
1125 . grade - the grading factor
1126 
1127   Level: advanced
1128 
1129 .seealso: DMForestSetGradeFactor()
1130 @*/
1131 PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade)
1132 {
1133   DM_Forest *forest = (DM_Forest*) dm->data;
1134 
1135   PetscFunctionBegin;
1136   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1137   PetscValidIntPointer(grade,2);
1138   *grade = forest->gradeFactor;
1139   PetscFunctionReturn(0);
1140 }
1141 
1142 #undef __FUNCT__
1143 #define __FUNCT__ "DMForestSetCellWeightFactor"
1144 /*@
1145   DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes
1146   the cell weight (see DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be
1147   (cellWeight) * (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on
1148   its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the
1149   computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
1150 
1151   Logically collective on dm
1152 
1153   Input Parameters:
1154 + dm - the forest
1155 - weightsFactors - default 1.
1156 
1157   Level: advanced
1158 
1159 .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights()
1160 @*/
1161 PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor)
1162 {
1163   DM_Forest *forest = (DM_Forest*) dm->data;
1164 
1165   PetscFunctionBegin;
1166   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1167   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup");
1168   forest->weightsFactor = weightsFactor;
1169   PetscFunctionReturn(0);
1170 }
1171 
1172 #undef __FUNCT__
1173 #define __FUNCT__ "DMForestGetCellWeightFactor"
1174 /*@
1175   DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see
1176   DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be (cellWeight) *
1177   (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on its level; a
1178   factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation
1179   associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
1180 
1181   Not collective
1182 
1183   Input Parameter:
1184 . dm - the forest
1185 
1186   Output Parameter:
1187 . weightsFactors - default 1.
1188 
1189   Level: advanced
1190 
1191 .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights()
1192 @*/
1193 PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor)
1194 {
1195   DM_Forest *forest = (DM_Forest*) dm->data;
1196 
1197   PetscFunctionBegin;
1198   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1199   PetscValidRealPointer(weightsFactor,2);
1200   *weightsFactor = forest->weightsFactor;
1201   PetscFunctionReturn(0);
1202 }
1203 
1204 #undef __FUNCT__
1205 #define __FUNCT__ "DMForestGetCellChart"
1206 /*@
1207   DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process
1208 
1209   Not collective
1210 
1211   Input Parameter:
1212 . dm - the forest
1213 
1214   Output Parameters:
1215 + cStart - the first cell on this process
1216 - cEnd - one after the final cell on this process
1217 
1218   Level: intermediate
1219 
1220 .seealso: DMForestGetCellSF()
1221 @*/
1222 PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd)
1223 {
1224   DM_Forest      *forest = (DM_Forest*) dm->data;
1225   PetscErrorCode ierr;
1226 
1227   PetscFunctionBegin;
1228   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1229   PetscValidIntPointer(cStart,2);
1230   PetscValidIntPointer(cEnd,2);
1231   if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) {
1232     ierr = forest->createcellchart(dm,&forest->cStart,&forest->cEnd);CHKERRQ(ierr);
1233   }
1234   *cStart =  forest->cStart;
1235   *cEnd   =  forest->cEnd;
1236   PetscFunctionReturn(0);
1237 }
1238 
1239 #undef __FUNCT__
1240 #define __FUNCT__ "DMForestGetCellSF"
1241 /*@
1242   DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes
1243 
1244   Not collective
1245 
1246   Input Parameter:
1247 . dm - the forest
1248 
1249   Output Parameter:
1250 . cellSF - the PetscSF
1251 
1252   Level: intermediate
1253 
1254 .seealso: DMForestGetCellChart()
1255 @*/
1256 PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF)
1257 {
1258   DM_Forest      *forest = (DM_Forest*) dm->data;
1259   PetscErrorCode ierr;
1260 
1261   PetscFunctionBegin;
1262   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1263   PetscValidPointer(cellSF,2);
1264   if ((!forest->cellSF) && forest->createcellsf) {
1265     ierr = forest->createcellsf(dm,&forest->cellSF);CHKERRQ(ierr);
1266   }
1267   *cellSF = forest->cellSF;
1268   PetscFunctionReturn(0);
1269 }
1270 
1271 #undef __FUNCT__
1272 #define __FUNCT__ "DMForestSetAdaptivityLabel"
1273 /*@C
1274   DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see
1275   DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination).  The
1276   interpretation of the label values is up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and
1277   DM_FOREST_COARSEN have been reserved as choices that should be accepted by all subtypes.
1278 
1279   Logically collective on dm
1280 
1281   Input Parameters:
1282 - dm - the forest
1283 + adaptLabel - the name of the label in the pre-adaptation forest
1284 
1285   Level: intermediate
1286 
1287 .seealso DMForestGetAdaptivityLabel()
1288 @*/
1289 PetscErrorCode DMForestSetAdaptivityLabel(DM dm, const char * adaptLabel)
1290 {
1291   DM_Forest      *forest = (DM_Forest*) dm->data;
1292   PetscErrorCode ierr;
1293 
1294   PetscFunctionBegin;
1295   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1296   ierr = PetscFree(forest->adaptLabel);CHKERRQ(ierr);
1297   ierr = PetscStrallocpy(adaptLabel,&forest->adaptLabel);CHKERRQ(ierr);
1298   PetscFunctionReturn(0);
1299 }
1300 
1301 #undef __FUNCT__
1302 #define __FUNCT__ "DMForestGetAdaptivityLabel"
1303 /*@C
1304   DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that
1305   holds the adaptation flags (refinement, coarsening, or some combination).  The interpretation of the label values is
1306   up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and DM_FOREST_COARSEN have been reserved as
1307   choices that should be accepted by all subtypes.
1308 
1309   Not collective
1310 
1311   Input Parameter:
1312 . dm - the forest
1313 
1314   Output Parameter:
1315 . adaptLabel - the name of the label in the pre-adaptation forest
1316 
1317   Level: intermediate
1318 
1319 .seealso DMForestSetAdaptivityLabel()
1320 @*/
1321 PetscErrorCode DMForestGetAdaptivityLabel(DM dm, const char ** adaptLabel)
1322 {
1323   DM_Forest *forest = (DM_Forest*) dm->data;
1324 
1325   PetscFunctionBegin;
1326   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1327   *adaptLabel = forest->adaptLabel;
1328   PetscFunctionReturn(0);
1329 }
1330 
1331 #undef __FUNCT__
1332 #define __FUNCT__ "DMForestSetCellWeights"
1333 /*@
1334   DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
1335   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
1336   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
1337   of 1.
1338 
1339   Logically collective on dm
1340 
1341   Input Parameters:
1342 + dm - the forest
1343 . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
1344 - copyMode - how weights should reference weights
1345 
1346   Level: advanced
1347 
1348 .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity()
1349 @*/
1350 PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode)
1351 {
1352   DM_Forest      *forest = (DM_Forest*) dm->data;
1353   PetscInt       cStart, cEnd;
1354   PetscErrorCode ierr;
1355 
1356   PetscFunctionBegin;
1357   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1358   ierr = DMForestGetCellChart(dm,&cStart,&cEnd);CHKERRQ(ierr);
1359   if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd);
1360   if (copyMode == PETSC_COPY_VALUES) {
1361     if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) {
1362       ierr = PetscMalloc1(cEnd-cStart,&forest->cellWeights);CHKERRQ(ierr);
1363     }
1364     ierr                        = PetscMemcpy(forest->cellWeights,weights,(cEnd-cStart)*sizeof(*weights));CHKERRQ(ierr);
1365     forest->cellWeightsCopyMode = PETSC_OWN_POINTER;
1366     PetscFunctionReturn(0);
1367   }
1368   if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) {
1369     ierr = PetscFree(forest->cellWeights);CHKERRQ(ierr);
1370   }
1371   forest->cellWeights         = weights;
1372   forest->cellWeightsCopyMode = copyMode;
1373   PetscFunctionReturn(0);
1374 }
1375 
1376 #undef __FUNCT__
1377 #define __FUNCT__ "DMForestGetCellWeights"
1378 /*@
1379   DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
1380   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
1381   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
1382   of 1.
1383 
1384   Not collective
1385 
1386   Input Parameter:
1387 . dm - the forest
1388 
1389   Output Parameter:
1390 . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
1391 
1392   Level: advanced
1393 
1394 .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity()
1395 @*/
1396 PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights)
1397 {
1398   DM_Forest *forest = (DM_Forest*) dm->data;
1399 
1400   PetscFunctionBegin;
1401   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1402   PetscValidPointer(weights,2);
1403   *weights = forest->cellWeights;
1404   PetscFunctionReturn(0);
1405 }
1406 
1407 #undef __FUNCT__
1408 #define __FUNCT__ "DMForestSetWeightCapacity"
1409 /*@
1410   DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning
1411   a pre-adaptation forest (see DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each
1412   process's cells to the process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that
1413   the current process should not have any cells after repartitioning.
1414 
1415   Logically Collective on dm
1416 
1417   Input parameters:
1418 + dm - the forest
1419 - capacity - this process's capacity
1420 
1421   Level: advanced
1422 
1423 .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
1424 @*/
1425 PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity)
1426 {
1427   DM_Forest *forest = (DM_Forest*) dm->data;
1428 
1429   PetscFunctionBegin;
1430   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1431   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup");
1432   if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity);
1433   forest->weightCapacity = capacity;
1434   PetscFunctionReturn(0);
1435 }
1436 
1437 #undef __FUNCT__
1438 #define __FUNCT__ "DMForestGetWeightCapacity"
1439 /*@
1440   DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see
1441   DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each process's cells to the
1442   process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that the current process
1443   should not have any cells after repartitioning.
1444 
1445   Not collective
1446 
1447   Input parameter:
1448 . dm - the forest
1449 
1450   Output parameter:
1451 . capacity - this process's capacity
1452 
1453   Level: advanced
1454 
1455 .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
1456 @*/
1457 PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity)
1458 {
1459   DM_Forest *forest = (DM_Forest*) dm->data;
1460 
1461   PetscFunctionBegin;
1462   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1463   PetscValidRealPointer(capacity,2);
1464   *capacity = forest->weightCapacity;
1465   PetscFunctionReturn(0);
1466 }
1467 
1468 #undef __FUNCT__
1469 #define __FUNCT__ "DMSetFromOptions_Forest"
1470 PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm)
1471 {
1472   DM_Forest                  *forest = (DM_Forest*) dm->data;
1473   PetscBool                  flg, flg1, flg2, flg3, flg4;
1474   DMForestTopology           oldTopo;
1475   char                       stringBuffer[256];
1476   PetscViewer                viewer;
1477   PetscViewerFormat          format;
1478   PetscInt                   adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade;
1479   PetscReal                  weightsFactor;
1480   DMForestAdaptivityStrategy adaptStrategy;
1481   PetscErrorCode             ierr;
1482 
1483   PetscFunctionBegin;
1484   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1485   forest->setfromoptionscalled = PETSC_TRUE;
1486   ierr                         = DMForestGetTopology(dm, &oldTopo);CHKERRQ(ierr);
1487   ierr                         = PetscOptionsHead(PetscOptionsObject,"DMForest Options");CHKERRQ(ierr);
1488   ierr                         = PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);CHKERRQ(ierr);
1489   ierr                         = PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);CHKERRQ(ierr);
1490   ierr                         = PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);CHKERRQ(ierr);
1491   ierr                         = PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);CHKERRQ(ierr);
1492   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}");
1493   if (flg1) {
1494     ierr = DMForestSetTopology(dm,(DMForestTopology)stringBuffer);CHKERRQ(ierr);
1495     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1496     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
1497   }
1498   if (flg2) {
1499     DM base;
1500 
1501     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&base);CHKERRQ(ierr);
1502     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1503     ierr = DMLoad(base,viewer);CHKERRQ(ierr);
1504     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1505     ierr = DMForestSetBaseDM(dm,base);CHKERRQ(ierr);
1506     ierr = DMDestroy(&base);CHKERRQ(ierr);
1507     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
1508     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
1509   }
1510   if (flg3) {
1511     DM coarse;
1512 
1513     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&coarse);CHKERRQ(ierr);
1514     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1515     ierr = DMLoad(coarse,viewer);CHKERRQ(ierr);
1516     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1517     ierr = DMForestSetAdaptivityForest(dm,coarse);CHKERRQ(ierr);
1518     ierr = DMDestroy(&coarse);CHKERRQ(ierr);
1519     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
1520     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1521   }
1522   if (flg4) {
1523     DM fine;
1524 
1525     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&fine);CHKERRQ(ierr);
1526     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1527     ierr = DMLoad(fine,viewer);CHKERRQ(ierr);
1528     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1529     ierr = DMForestSetAdaptivityForest(dm,fine);CHKERRQ(ierr);
1530     ierr = DMDestroy(&fine);CHKERRQ(ierr);
1531     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
1532     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1533   }
1534   ierr = DMForestGetAdjacencyDimension(dm,&adjDim);CHKERRQ(ierr);
1535   ierr = PetscOptionsInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg);CHKERRQ(ierr);
1536   if (flg) {
1537     ierr = DMForestSetAdjacencyDimension(dm,adjDim);CHKERRQ(ierr);
1538   } else {
1539     ierr = DMForestGetAdjacencyCodimension(dm,&adjCodim);CHKERRQ(ierr);
1540     ierr = PetscOptionsInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg);CHKERRQ(ierr);
1541     if (flg) {
1542       ierr = DMForestSetAdjacencyCodimension(dm,adjCodim);CHKERRQ(ierr);
1543     }
1544   }
1545   ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr);
1546   ierr = PetscOptionsInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg);CHKERRQ(ierr);
1547   if (flg) {
1548     ierr = DMForestSetPartitionOverlap(dm,overlap);CHKERRQ(ierr);
1549   }
1550 #if 0
1551   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);
1552   if (flg) {
1553     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1554     ierr = DMForestSetInitialRefinement(dm,minRefinement);CHKERRQ(ierr);
1555   }
1556   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);
1557   if (flg) {
1558     ierr = DMForestSetMinimumRefinement(dm,0);CHKERRQ(ierr);
1559     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
1560   }
1561 #endif
1562   ierr = DMForestGetMinimumRefinement(dm,&minRefinement);CHKERRQ(ierr);
1563   ierr = PetscOptionsInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg);CHKERRQ(ierr);
1564   if (flg) {
1565     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1566   }
1567   ierr = DMForestGetInitialRefinement(dm,&initRefinement);CHKERRQ(ierr);
1568   ierr = PetscOptionsInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg);CHKERRQ(ierr);
1569   if (flg) {
1570     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
1571   }
1572   ierr = DMForestGetMaximumRefinement(dm,&maxRefinement);CHKERRQ(ierr);
1573   ierr = PetscOptionsInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg);CHKERRQ(ierr);
1574   if (flg) {
1575     ierr = DMForestSetMaximumRefinement(dm,maxRefinement);CHKERRQ(ierr);
1576   }
1577   ierr = DMForestGetAdaptivityStrategy(dm,&adaptStrategy);CHKERRQ(ierr);
1578   ierr = PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,256,&flg);CHKERRQ(ierr);
1579   if (flg) {
1580     ierr = DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);CHKERRQ(ierr);
1581   }
1582   ierr = DMForestGetGradeFactor(dm,&grade);CHKERRQ(ierr);
1583   ierr = PetscOptionsInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg);CHKERRQ(ierr);
1584   if (flg) {
1585     ierr = DMForestSetGradeFactor(dm,grade);CHKERRQ(ierr);
1586   }
1587   ierr = DMForestGetCellWeightFactor(dm,&weightsFactor);CHKERRQ(ierr);
1588   ierr = PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);CHKERRQ(ierr);
1589   if (flg) {
1590     ierr = DMForestSetCellWeightFactor(dm,weightsFactor);CHKERRQ(ierr);
1591   }
1592   ierr = PetscOptionsTail();CHKERRQ(ierr);
1593   PetscFunctionReturn(0);
1594 }
1595 
1596 #undef __FUNCT__
1597 #define __FUNCT__ "DMCreateSubDM_Forest"
1598 PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
1599 {
1600   PetscErrorCode ierr;
1601 
1602   PetscFunctionBegin;
1603   if (subdm) {ierr = DMClone(dm, subdm);CHKERRQ(ierr);}
1604   ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr);
1605   PetscFunctionReturn(0);
1606 }
1607 
1608 #undef __FUNCT__
1609 #define __FUNCT__ "DMRefine_Forest"
1610 PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
1611 {
1612   DMLabel        refine;
1613   DM             fineDM;
1614   PetscErrorCode ierr;
1615 
1616   PetscFunctionBegin;
1617   ierr = DMGetFineDM(dm,&fineDM);CHKERRQ(ierr);
1618   if (fineDM) {
1619     ierr       = PetscObjectReference((PetscObject)fineDM);CHKERRQ(ierr);
1620     *dmRefined = fineDM;
1621     PetscFunctionReturn(0);
1622   }
1623   ierr = DMForestTemplate(dm,comm,dmRefined);CHKERRQ(ierr);
1624   ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr);
1625   if (!refine) {
1626     ierr = DMCreateLabel(dm,"refine");CHKERRQ(ierr);
1627     ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr);
1628     ierr = DMLabelSetDefaultValue(refine,DM_FOREST_REFINE);CHKERRQ(ierr);
1629   }
1630   ierr = DMForestSetAdaptivityLabel(*dmRefined,"refine");CHKERRQ(ierr);
1631   PetscFunctionReturn(0);
1632 }
1633 
1634 #undef __FUNCT__
1635 #define __FUNCT__ "DMCoarsen_Forest"
1636 PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
1637 {
1638   DMLabel        coarsen;
1639   DM             coarseDM;
1640   PetscErrorCode ierr;
1641 
1642   PetscFunctionBegin;
1643   {
1644     PetscMPIInt mpiComparison;
1645     MPI_Comm    dmcomm = PetscObjectComm((PetscObject)dm);
1646 
1647     ierr = MPI_Comm_compare(comm, dmcomm, &mpiComparison);CHKERRQ(ierr);
1648     if (mpiComparison != MPI_IDENT && mpiComparison != MPI_CONGRUENT) SETERRQ(dmcomm,PETSC_ERR_SUP,"No support for different communicators yet");
1649   }
1650   ierr = DMGetCoarseDM(dm,&coarseDM);CHKERRQ(ierr);
1651   if (coarseDM) {
1652     ierr         = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr);
1653     *dmCoarsened = coarseDM;
1654     PetscFunctionReturn(0);
1655   }
1656   ierr = DMForestTemplate(dm,comm,dmCoarsened);CHKERRQ(ierr);
1657   ierr = DMForestSetAdaptivityPurpose(coarseDM,DM_FOREST_COARSEN);CHKERRQ(ierr);
1658   ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr);
1659   if (!coarsen) {
1660     ierr = DMCreateLabel(dm,"coarsen");CHKERRQ(ierr);
1661     ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr);
1662     ierr = DMLabelSetDefaultValue(coarsen,DM_FOREST_COARSEN);CHKERRQ(ierr);
1663   }
1664   ierr = DMForestSetAdaptivityLabel(*dmCoarsened,"coarsen");CHKERRQ(ierr);
1665   PetscFunctionReturn(0);
1666 }
1667 
1668 #undef __FUNCT__
1669 #define __FUNCT__ "DMInitialize_Forest"
1670 static PetscErrorCode DMInitialize_Forest(DM dm)
1671 {
1672   PetscErrorCode ierr;
1673 
1674   PetscFunctionBegin;
1675   ierr = PetscMemzero(dm->ops,sizeof(*(dm->ops)));CHKERRQ(ierr);
1676 
1677   dm->ops->clone          = DMClone_Forest;
1678   dm->ops->setfromoptions = DMSetFromOptions_Forest;
1679   dm->ops->destroy        = DMDestroy_Forest;
1680   dm->ops->createsubdm    = DMCreateSubDM_Forest;
1681   dm->ops->refine         = DMRefine_Forest;
1682   dm->ops->coarsen        = DMCoarsen_Forest;
1683   PetscFunctionReturn(0);
1684 }
1685 
1686 /*MC
1687   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.
1688 
1689   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().
1690 
1691   Level: advanced
1692 
1693 .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel()
1694 M*/
1695 
1696 #undef __FUNCT__
1697 #define __FUNCT__ "DMCreate_Forest"
1698 PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1699 {
1700   DM_Forest      *forest;
1701   PetscErrorCode ierr;
1702 
1703   PetscFunctionBegin;
1704   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1705   ierr                         = PetscNewLog(dm,&forest);CHKERRQ(ierr);
1706   dm->dim                      = 0;
1707   dm->data                     = forest;
1708   forest->refct                = 1;
1709   forest->data                 = NULL;
1710   forest->setfromoptionscalled = PETSC_FALSE;
1711   forest->topology             = NULL;
1712   forest->base                 = NULL;
1713   forest->adjDim               = PETSC_DEFAULT;
1714   forest->overlap              = PETSC_DEFAULT;
1715   forest->minRefinement        = PETSC_DEFAULT;
1716   forest->maxRefinement        = PETSC_DEFAULT;
1717   forest->initRefinement       = PETSC_DEFAULT;
1718   forest->cStart               = PETSC_DETERMINE;
1719   forest->cEnd                 = PETSC_DETERMINE;
1720   forest->cellSF               = 0;
1721   forest->adaptLabel           = NULL;
1722   forest->gradeFactor          = 2;
1723   forest->cellWeights          = NULL;
1724   forest->cellWeightsCopyMode  = PETSC_USE_POINTER;
1725   forest->weightsFactor        = 1.;
1726   forest->weightCapacity       = 1.;
1727   ierr                         = DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);CHKERRQ(ierr);
1728   ierr                         = DMInitialize_Forest(dm);CHKERRQ(ierr);
1729   PetscFunctionReturn(0);
1730 }
1731 
1732