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