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