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