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