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