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