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