xref: /petsc/src/dm/impls/forest/forest.c (revision ccc738f73cdca4da6a86aae33987f346ef28fee4)
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__ "DMForestGetComputeAdaptivitySF"
989 /*@
990   DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the
991   pre-adaptation forest to the post-adaptiation forest.  After DMSetUp() is called, these transfer PetscSFs can be
992   accessed with DMForestGetAdaptivitySF().
993 
994   Not collective
995 
996   Input Parameter:
997 . dm - the post-adaptation forest
998 
999   Output Parameter:
1000 . computeSF - default PETSC_TRUE
1001 
1002   Level: advanced
1003 
1004 .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
1005 @*/
1006 PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF)
1007 {
1008   DM_Forest *forest;
1009 
1010   PetscFunctionBegin;
1011   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1012   forest     = (DM_Forest*) dm->data;
1013   *computeSF = forest->computeAdaptSF;
1014   PetscFunctionReturn(0);
1015 }
1016 
1017 #undef __FUNCT__
1018 #define __FUNCT__ "DMForestGetAdaptivitySF"
1019 /*@
1020   DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest.
1021   Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be
1022   some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa.  Therefore there are two
1023   PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates
1024   pre-adaptation fine cells to post-adaptation coarse cells.
1025 
1026   Not collective
1027 
1028   Input Parameter:
1029   dm - the post-adaptation forest
1030 
1031   Output Parameter:
1032   preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post-
1033   coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre-
1034 
1035   Level: advanced
1036 
1037 .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF()
1038 @*/
1039 PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine)
1040 {
1041   DM_Forest      *forest;
1042   PetscErrorCode ierr;
1043 
1044   PetscFunctionBegin;
1045   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1046   ierr   = DMSetUp(dm);CHKERRQ(ierr);
1047   forest = (DM_Forest*) dm->data;
1048   if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine;
1049   if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine;
1050   PetscFunctionReturn(0);
1051 }
1052 
1053 #undef __FUNCT__
1054 #define __FUNCT__ "DMForestSetGradeFactor"
1055 /*@
1056   DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to
1057   indicate that the diameter of neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may
1058   only support one particular choice of grading factor.
1059 
1060   Logically collective on dm
1061 
1062   Input Parameters:
1063 + dm - the forest
1064 - grade - the grading factor
1065 
1066   Level: advanced
1067 
1068 .seealso: DMForestGetGradeFactor()
1069 @*/
1070 PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade)
1071 {
1072   DM_Forest *forest = (DM_Forest*) dm->data;
1073 
1074   PetscFunctionBegin;
1075   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1076   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup");
1077   forest->gradeFactor = grade;
1078   PetscFunctionReturn(0);
1079 }
1080 
1081 #undef __FUNCT__
1082 #define __FUNCT__ "DMForestGetGradeFactor"
1083 /*@
1084   DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of
1085   neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may only support one particular
1086   choice of grading factor.
1087 
1088   Not collective
1089 
1090   Input Parameter:
1091 . dm - the forest
1092 
1093   Output Parameter:
1094 . grade - the grading factor
1095 
1096   Level: advanced
1097 
1098 .seealso: DMForestSetGradeFactor()
1099 @*/
1100 PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade)
1101 {
1102   DM_Forest *forest = (DM_Forest*) dm->data;
1103 
1104   PetscFunctionBegin;
1105   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1106   PetscValidIntPointer(grade,2);
1107   *grade = forest->gradeFactor;
1108   PetscFunctionReturn(0);
1109 }
1110 
1111 #undef __FUNCT__
1112 #define __FUNCT__ "DMForestSetCellWeightFactor"
1113 /*@
1114   DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes
1115   the cell weight (see DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be
1116   (cellWeight) * (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on
1117   its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the
1118   computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
1119 
1120   Logically collective on dm
1121 
1122   Input Parameters:
1123 + dm - the forest
1124 - weightsFactors - default 1.
1125 
1126   Level: advanced
1127 
1128 .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights()
1129 @*/
1130 PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor)
1131 {
1132   DM_Forest *forest = (DM_Forest*) dm->data;
1133 
1134   PetscFunctionBegin;
1135   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1136   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup");
1137   forest->weightsFactor = weightsFactor;
1138   PetscFunctionReturn(0);
1139 }
1140 
1141 #undef __FUNCT__
1142 #define __FUNCT__ "DMForestGetCellWeightFactor"
1143 /*@
1144   DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see
1145   DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be (cellWeight) *
1146   (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on its level; a
1147   factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation
1148   associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
1149 
1150   Not collective
1151 
1152   Input Parameter:
1153 . dm - the forest
1154 
1155   Output Parameter:
1156 . weightsFactors - default 1.
1157 
1158   Level: advanced
1159 
1160 .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights()
1161 @*/
1162 PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor)
1163 {
1164   DM_Forest *forest = (DM_Forest*) dm->data;
1165 
1166   PetscFunctionBegin;
1167   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1168   PetscValidRealPointer(weightsFactor,2);
1169   *weightsFactor = forest->weightsFactor;
1170   PetscFunctionReturn(0);
1171 }
1172 
1173 #undef __FUNCT__
1174 #define __FUNCT__ "DMForestGetCellChart"
1175 /*@
1176   DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process
1177 
1178   Not collective
1179 
1180   Input Parameter:
1181 . dm - the forest
1182 
1183   Output Parameters:
1184 + cStart - the first cell on this process
1185 - cEnd - one after the final cell on this process
1186 
1187   Level: intermediate
1188 
1189 .seealso: DMForestGetCellSF()
1190 @*/
1191 PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd)
1192 {
1193   DM_Forest      *forest = (DM_Forest*) dm->data;
1194   PetscErrorCode ierr;
1195 
1196   PetscFunctionBegin;
1197   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1198   PetscValidIntPointer(cStart,2);
1199   PetscValidIntPointer(cEnd,2);
1200   if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) {
1201     ierr = forest->createcellchart(dm,&forest->cStart,&forest->cEnd);CHKERRQ(ierr);
1202   }
1203   *cStart =  forest->cStart;
1204   *cEnd   =  forest->cEnd;
1205   PetscFunctionReturn(0);
1206 }
1207 
1208 #undef __FUNCT__
1209 #define __FUNCT__ "DMForestGetCellSF"
1210 /*@
1211   DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes
1212 
1213   Not collective
1214 
1215   Input Parameter:
1216 . dm - the forest
1217 
1218   Output Parameter:
1219 . cellSF - the PetscSF
1220 
1221   Level: intermediate
1222 
1223 .seealso: DMForestGetCellChart()
1224 @*/
1225 PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF)
1226 {
1227   DM_Forest      *forest = (DM_Forest*) dm->data;
1228   PetscErrorCode ierr;
1229 
1230   PetscFunctionBegin;
1231   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1232   PetscValidPointer(cellSF,2);
1233   if ((!forest->cellSF) && forest->createcellsf) {
1234     ierr = forest->createcellsf(dm,&forest->cellSF);CHKERRQ(ierr);
1235   }
1236   *cellSF = forest->cellSF;
1237   PetscFunctionReturn(0);
1238 }
1239 
1240 #undef __FUNCT__
1241 #define __FUNCT__ "DMForestSetAdaptivityLabel"
1242 /*@C
1243   DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see
1244   DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination).  The
1245   interpretation of the label values is up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and
1246   DM_FOREST_COARSEN have been reserved as choices that should be accepted by all subtypes.
1247 
1248   Logically collective on dm
1249 
1250   Input Parameters:
1251 - dm - the forest
1252 + adaptLabel - the name of the label in the pre-adaptation forest
1253 
1254   Level: intermediate
1255 
1256 .seealso DMForestGetAdaptivityLabel()
1257 @*/
1258 PetscErrorCode DMForestSetAdaptivityLabel(DM dm, const char * adaptLabel)
1259 {
1260   DM_Forest      *forest = (DM_Forest*) dm->data;
1261   PetscErrorCode ierr;
1262 
1263   PetscFunctionBegin;
1264   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1265   ierr = PetscFree(forest->adaptLabel);CHKERRQ(ierr);
1266   ierr = PetscStrallocpy(adaptLabel,&forest->adaptLabel);CHKERRQ(ierr);
1267   PetscFunctionReturn(0);
1268 }
1269 
1270 #undef __FUNCT__
1271 #define __FUNCT__ "DMForestGetAdaptivityLabel"
1272 /*@C
1273   DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that
1274   holds the adaptation flags (refinement, coarsening, or some combination).  The interpretation of the label values is
1275   up to the subtype of DMForest, but DM_FOREST_KEEP, DM_FOREST_REFINE, and DM_FOREST_COARSEN have been reserved as
1276   choices that should be accepted by all subtypes.
1277 
1278   Not collective
1279 
1280   Input Parameter:
1281 . dm - the forest
1282 
1283   Output Parameter:
1284 . adaptLabel - the name of the label in the pre-adaptation forest
1285 
1286   Level: intermediate
1287 
1288 .seealso DMForestSetAdaptivityLabel()
1289 @*/
1290 PetscErrorCode DMForestGetAdaptivityLabel(DM dm, const char ** adaptLabel)
1291 {
1292   DM_Forest *forest = (DM_Forest*) dm->data;
1293 
1294   PetscFunctionBegin;
1295   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1296   *adaptLabel = forest->adaptLabel;
1297   PetscFunctionReturn(0);
1298 }
1299 
1300 #undef __FUNCT__
1301 #define __FUNCT__ "DMForestSetCellWeights"
1302 /*@
1303   DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
1304   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
1305   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
1306   of 1.
1307 
1308   Logically collective on dm
1309 
1310   Input Parameters:
1311 + dm - the forest
1312 . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
1313 - copyMode - how weights should reference weights
1314 
1315   Level: advanced
1316 
1317 .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity()
1318 @*/
1319 PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode)
1320 {
1321   DM_Forest      *forest = (DM_Forest*) dm->data;
1322   PetscInt       cStart, cEnd;
1323   PetscErrorCode ierr;
1324 
1325   PetscFunctionBegin;
1326   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1327   ierr = DMForestGetCellChart(dm,&cStart,&cEnd);CHKERRQ(ierr);
1328   if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd);
1329   if (copyMode == PETSC_COPY_VALUES) {
1330     if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) {
1331       ierr = PetscMalloc1(cEnd-cStart,&forest->cellWeights);CHKERRQ(ierr);
1332     }
1333     ierr                        = PetscMemcpy(forest->cellWeights,weights,(cEnd-cStart)*sizeof(*weights));CHKERRQ(ierr);
1334     forest->cellWeightsCopyMode = PETSC_OWN_POINTER;
1335     PetscFunctionReturn(0);
1336   }
1337   if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) {
1338     ierr = PetscFree(forest->cellWeights);CHKERRQ(ierr);
1339   }
1340   forest->cellWeights         = weights;
1341   forest->cellWeightsCopyMode = copyMode;
1342   PetscFunctionReturn(0);
1343 }
1344 
1345 #undef __FUNCT__
1346 #define __FUNCT__ "DMForestGetCellWeights"
1347 /*@
1348   DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
1349   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
1350   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
1351   of 1.
1352 
1353   Not collective
1354 
1355   Input Parameter:
1356 . dm - the forest
1357 
1358   Output Parameter:
1359 . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
1360 
1361   Level: advanced
1362 
1363 .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity()
1364 @*/
1365 PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights)
1366 {
1367   DM_Forest *forest = (DM_Forest*) dm->data;
1368 
1369   PetscFunctionBegin;
1370   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1371   PetscValidPointer(weights,2);
1372   *weights = forest->cellWeights;
1373   PetscFunctionReturn(0);
1374 }
1375 
1376 #undef __FUNCT__
1377 #define __FUNCT__ "DMForestSetWeightCapacity"
1378 /*@
1379   DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning
1380   a pre-adaptation forest (see DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each
1381   process's cells to the process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that
1382   the current process should not have any cells after repartitioning.
1383 
1384   Logically Collective on dm
1385 
1386   Input parameters:
1387 + dm - the forest
1388 - capacity - this process's capacity
1389 
1390   Level: advanced
1391 
1392 .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
1393 @*/
1394 PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity)
1395 {
1396   DM_Forest *forest = (DM_Forest*) dm->data;
1397 
1398   PetscFunctionBegin;
1399   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1400   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup");
1401   if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity);
1402   forest->weightCapacity = capacity;
1403   PetscFunctionReturn(0);
1404 }
1405 
1406 #undef __FUNCT__
1407 #define __FUNCT__ "DMForestGetWeightCapacity"
1408 /*@
1409   DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see
1410   DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each process's cells to the
1411   process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that the current process
1412   should not have any cells after repartitioning.
1413 
1414   Not collective
1415 
1416   Input parameter:
1417 . dm - the forest
1418 
1419   Output parameter:
1420 . capacity - this process's capacity
1421 
1422   Level: advanced
1423 
1424 .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
1425 @*/
1426 PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity)
1427 {
1428   DM_Forest *forest = (DM_Forest*) dm->data;
1429 
1430   PetscFunctionBegin;
1431   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1432   PetscValidRealPointer(capacity,2);
1433   *capacity = forest->weightCapacity;
1434   PetscFunctionReturn(0);
1435 }
1436 
1437 #undef __FUNCT__
1438 #define __FUNCT__ "DMSetFromOptions_Forest"
1439 PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm)
1440 {
1441   DM_Forest                  *forest = (DM_Forest*) dm->data;
1442   PetscBool                  flg, flg1, flg2, flg3, flg4;
1443   DMForestTopology           oldTopo;
1444   char                       stringBuffer[256];
1445   PetscViewer                viewer;
1446   PetscViewerFormat          format;
1447   PetscInt                   adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade;
1448   PetscReal                  weightsFactor;
1449   DMForestAdaptivityStrategy adaptStrategy;
1450   PetscErrorCode             ierr;
1451 
1452   PetscFunctionBegin;
1453   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1454   forest->setfromoptionscalled = PETSC_TRUE;
1455   ierr                         = DMForestGetTopology(dm, &oldTopo);CHKERRQ(ierr);
1456   ierr                         = PetscOptionsHead(PetscOptionsObject,"DMForest Options");CHKERRQ(ierr);
1457   ierr                         = PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);CHKERRQ(ierr);
1458   ierr                         = PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);CHKERRQ(ierr);
1459   ierr                         = PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);CHKERRQ(ierr);
1460   ierr                         = PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);CHKERRQ(ierr);
1461   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}");
1462   if (flg1) {
1463     ierr = DMForestSetTopology(dm,(DMForestTopology)stringBuffer);CHKERRQ(ierr);
1464     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1465     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
1466   }
1467   if (flg2) {
1468     DM base;
1469 
1470     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&base);CHKERRQ(ierr);
1471     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1472     ierr = DMLoad(base,viewer);CHKERRQ(ierr);
1473     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1474     ierr = DMForestSetBaseDM(dm,base);CHKERRQ(ierr);
1475     ierr = DMDestroy(&base);CHKERRQ(ierr);
1476     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
1477     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
1478   }
1479   if (flg3) {
1480     DM coarse;
1481 
1482     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&coarse);CHKERRQ(ierr);
1483     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1484     ierr = DMLoad(coarse,viewer);CHKERRQ(ierr);
1485     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1486     ierr = DMForestSetAdaptivityForest(dm,coarse);CHKERRQ(ierr);
1487     ierr = DMDestroy(&coarse);CHKERRQ(ierr);
1488     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
1489     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1490   }
1491   if (flg4) {
1492     DM fine;
1493 
1494     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&fine);CHKERRQ(ierr);
1495     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1496     ierr = DMLoad(fine,viewer);CHKERRQ(ierr);
1497     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1498     ierr = DMForestSetAdaptivityForest(dm,fine);CHKERRQ(ierr);
1499     ierr = DMDestroy(&fine);CHKERRQ(ierr);
1500     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
1501     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1502   }
1503   ierr = DMForestGetAdjacencyDimension(dm,&adjDim);CHKERRQ(ierr);
1504   ierr = PetscOptionsInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg);CHKERRQ(ierr);
1505   if (flg) {
1506     ierr = DMForestSetAdjacencyDimension(dm,adjDim);CHKERRQ(ierr);
1507   } else {
1508     ierr = DMForestGetAdjacencyCodimension(dm,&adjCodim);CHKERRQ(ierr);
1509     ierr = PetscOptionsInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg);CHKERRQ(ierr);
1510     if (flg) {
1511       ierr = DMForestSetAdjacencyCodimension(dm,adjCodim);CHKERRQ(ierr);
1512     }
1513   }
1514   ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr);
1515   ierr = PetscOptionsInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg);CHKERRQ(ierr);
1516   if (flg) {
1517     ierr = DMForestSetPartitionOverlap(dm,overlap);CHKERRQ(ierr);
1518   }
1519 #if 0
1520   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);
1521   if (flg) {
1522     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1523     ierr = DMForestSetInitialRefinement(dm,minRefinement);CHKERRQ(ierr);
1524   }
1525   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);
1526   if (flg) {
1527     ierr = DMForestSetMinimumRefinement(dm,0);CHKERRQ(ierr);
1528     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
1529   }
1530 #endif
1531   ierr = DMForestGetMinimumRefinement(dm,&minRefinement);CHKERRQ(ierr);
1532   ierr = PetscOptionsInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg);CHKERRQ(ierr);
1533   if (flg) {
1534     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1535   }
1536   ierr = DMForestGetInitialRefinement(dm,&initRefinement);CHKERRQ(ierr);
1537   ierr = PetscOptionsInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg);CHKERRQ(ierr);
1538   if (flg) {
1539     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
1540   }
1541   ierr = DMForestGetMaximumRefinement(dm,&maxRefinement);CHKERRQ(ierr);
1542   ierr = PetscOptionsInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg);CHKERRQ(ierr);
1543   if (flg) {
1544     ierr = DMForestSetMaximumRefinement(dm,maxRefinement);CHKERRQ(ierr);
1545   }
1546   ierr = DMForestGetAdaptivityStrategy(dm,&adaptStrategy);CHKERRQ(ierr);
1547   ierr = PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,256,&flg);CHKERRQ(ierr);
1548   if (flg) {
1549     ierr = DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);CHKERRQ(ierr);
1550   }
1551   ierr = DMForestGetGradeFactor(dm,&grade);CHKERRQ(ierr);
1552   ierr = PetscOptionsInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg);CHKERRQ(ierr);
1553   if (flg) {
1554     ierr = DMForestSetGradeFactor(dm,grade);CHKERRQ(ierr);
1555   }
1556   ierr = DMForestGetCellWeightFactor(dm,&weightsFactor);CHKERRQ(ierr);
1557   ierr = PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);CHKERRQ(ierr);
1558   if (flg) {
1559     ierr = DMForestSetCellWeightFactor(dm,weightsFactor);CHKERRQ(ierr);
1560   }
1561   ierr = PetscOptionsTail();CHKERRQ(ierr);
1562   PetscFunctionReturn(0);
1563 }
1564 
1565 #undef __FUNCT__
1566 #define __FUNCT__ "DMCreateSubDM_Forest"
1567 PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
1568 {
1569   PetscErrorCode ierr;
1570 
1571   PetscFunctionBegin;
1572   if (subdm) {ierr = DMClone(dm, subdm);CHKERRQ(ierr);}
1573   ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr);
1574   PetscFunctionReturn(0);
1575 }
1576 
1577 #undef __FUNCT__
1578 #define __FUNCT__ "DMRefine_Forest"
1579 PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
1580 {
1581   DMLabel        refine;
1582   DM             fineDM;
1583   PetscErrorCode ierr;
1584 
1585   PetscFunctionBegin;
1586   ierr = DMGetFineDM(dm,&fineDM);CHKERRQ(ierr);
1587   if (fineDM) {
1588     ierr       = PetscObjectReference((PetscObject)fineDM);CHKERRQ(ierr);
1589     *dmRefined = fineDM;
1590     PetscFunctionReturn(0);
1591   }
1592   ierr = DMForestTemplate(dm,comm,dmRefined);CHKERRQ(ierr);
1593   ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr);
1594   if (!refine) {
1595     ierr = DMCreateLabel(dm,"refine");CHKERRQ(ierr);
1596     ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr);
1597     ierr = DMLabelSetDefaultValue(refine,DM_FOREST_REFINE);CHKERRQ(ierr);
1598   }
1599   ierr = DMForestSetAdaptivityLabel(*dmRefined,"refine");CHKERRQ(ierr);
1600   PetscFunctionReturn(0);
1601 }
1602 
1603 #undef __FUNCT__
1604 #define __FUNCT__ "DMCoarsen_Forest"
1605 PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
1606 {
1607   DMLabel        coarsen;
1608   DM             coarseDM;
1609   PetscErrorCode ierr;
1610 
1611   PetscFunctionBegin;
1612   {
1613     PetscMPIInt mpiComparison;
1614     MPI_Comm    dmcomm = PetscObjectComm((PetscObject)dm);
1615 
1616     ierr = MPI_Comm_compare(comm, dmcomm, &mpiComparison);CHKERRQ(ierr);
1617     if (mpiComparison != MPI_IDENT && mpiComparison != MPI_CONGRUENT) SETERRQ(dmcomm,PETSC_ERR_SUP,"No support for different communicators yet");
1618   }
1619   ierr = DMGetCoarseDM(dm,&coarseDM);CHKERRQ(ierr);
1620   if (coarseDM) {
1621     ierr         = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr);
1622     *dmCoarsened = coarseDM;
1623     PetscFunctionReturn(0);
1624   }
1625   ierr = DMForestTemplate(dm,comm,dmCoarsened);CHKERRQ(ierr);
1626   ierr = DMForestSetAdaptivityPurpose(coarseDM,DM_FOREST_COARSEN);CHKERRQ(ierr);
1627   ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr);
1628   if (!coarsen) {
1629     ierr = DMCreateLabel(dm,"coarsen");CHKERRQ(ierr);
1630     ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr);
1631     ierr = DMLabelSetDefaultValue(coarsen,DM_FOREST_COARSEN);CHKERRQ(ierr);
1632   }
1633   ierr = DMForestSetAdaptivityLabel(*dmCoarsened,"coarsen");CHKERRQ(ierr);
1634   PetscFunctionReturn(0);
1635 }
1636 
1637 #undef __FUNCT__
1638 #define __FUNCT__ "DMInitialize_Forest"
1639 static PetscErrorCode DMInitialize_Forest(DM dm)
1640 {
1641   PetscErrorCode ierr;
1642 
1643   PetscFunctionBegin;
1644   ierr = PetscMemzero(dm->ops,sizeof(*(dm->ops)));CHKERRQ(ierr);
1645 
1646   dm->ops->clone          = DMClone_Forest;
1647   dm->ops->setfromoptions = DMSetFromOptions_Forest;
1648   dm->ops->destroy        = DMDestroy_Forest;
1649   dm->ops->createsubdm    = DMCreateSubDM_Forest;
1650   dm->ops->refine         = DMRefine_Forest;
1651   dm->ops->coarsen        = DMCoarsen_Forest;
1652   PetscFunctionReturn(0);
1653 }
1654 
1655 /*MC
1656   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.
1657 
1658   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().
1659 
1660   Level: advanced
1661 
1662 .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel()
1663 M*/
1664 
1665 #undef __FUNCT__
1666 #define __FUNCT__ "DMCreate_Forest"
1667 PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1668 {
1669   DM_Forest      *forest;
1670   PetscErrorCode ierr;
1671 
1672   PetscFunctionBegin;
1673   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1674   ierr                         = PetscNewLog(dm,&forest);CHKERRQ(ierr);
1675   dm->dim                      = 0;
1676   dm->data                     = forest;
1677   forest->refct                = 1;
1678   forest->data                 = NULL;
1679   forest->setfromoptionscalled = PETSC_FALSE;
1680   forest->topology             = NULL;
1681   forest->base                 = NULL;
1682   forest->adjDim               = PETSC_DEFAULT;
1683   forest->overlap              = PETSC_DEFAULT;
1684   forest->minRefinement        = PETSC_DEFAULT;
1685   forest->maxRefinement        = PETSC_DEFAULT;
1686   forest->initRefinement       = PETSC_DEFAULT;
1687   forest->cStart               = PETSC_DETERMINE;
1688   forest->cEnd                 = PETSC_DETERMINE;
1689   forest->cellSF               = 0;
1690   forest->adaptLabel           = NULL;
1691   forest->gradeFactor          = 2;
1692   forest->cellWeights          = NULL;
1693   forest->cellWeightsCopyMode  = PETSC_USE_POINTER;
1694   forest->weightsFactor        = 1.;
1695   forest->weightCapacity       = 1.;
1696   ierr                         = DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);CHKERRQ(ierr);
1697   ierr                         = DMInitialize_Forest(dm);CHKERRQ(ierr);
1698   PetscFunctionReturn(0);
1699 }
1700 
1701