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