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