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