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