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