xref: /petsc/src/dm/impls/forest/forest.c (revision 2f6eced2a19e978d64f27de66fbc6c26cc5c7934)
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   ierr = PetscObjectComposeFunction((PetscObject)dm,"DMAdaptLabel_C",NULL);CHKERRQ(ierr);
216   PetscFunctionReturn(0);
217 }
218 
219 /*@C
220   DMForestSetTopology - Set the topology of a DMForest during the pre-setup phase.  The topology is a string (e.g.
221   "cube", "shell") and can be interpreted by subtypes of DMFOREST) to construct the base DM of a forest durint
222   DMSetUp().
223 
224   Logically collective on dm
225 
226   Input parameters:
227 + dm - the forest
228 - topology - the topology of the forest
229 
230   Level: intermediate
231 
232 .seealso(): DMForestGetTopology(), DMForestSetBaseDM()
233 @*/
234 PetscErrorCode DMForestSetTopology(DM dm, DMForestTopology topology)
235 {
236   DM_Forest      *forest = (DM_Forest*) dm->data;
237   PetscErrorCode ierr;
238 
239   PetscFunctionBegin;
240   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
241   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the topology after setup");
242   ierr = PetscFree(forest->topology);CHKERRQ(ierr);
243   ierr = PetscStrallocpy((const char*)topology,(char**) &forest->topology);CHKERRQ(ierr);
244   PetscFunctionReturn(0);
245 }
246 
247 /*@C
248   DMForestGetTopology - Get a string describing the topology of a DMForest.
249 
250   Not collective
251 
252   Input parameter:
253 . dm - the forest
254 
255   Output parameter:
256 . topology - the topology of the forest (e.g., 'cube', 'shell')
257 
258   Level: intermediate
259 
260 .seealso: DMForestSetTopology()
261 @*/
262 PetscErrorCode DMForestGetTopology(DM dm, DMForestTopology *topology)
263 {
264   DM_Forest *forest = (DM_Forest*) dm->data;
265 
266   PetscFunctionBegin;
267   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
268   PetscValidPointer(topology,2);
269   *topology = forest->topology;
270   PetscFunctionReturn(0);
271 }
272 
273 /*@
274   DMForestSetBaseDM - During the pre-setup phase, set the DM that defines the base mesh of a DMForest forest.  The
275   forest will be hierarchically refined from the base, and all refinements/coarsenings of the forest will share its
276   base.  In general, two forest must share a bse to be comparable, to do things like construct interpolators.
277 
278   Logically collective on dm
279 
280   Input Parameters:
281 + dm - the forest
282 - base - the base DM of the forest
283 
284   Level: intermediate
285 
286 .seealso(): DMForestGetBaseDM()
287 @*/
288 PetscErrorCode DMForestSetBaseDM(DM dm, DM base)
289 {
290   DM_Forest      *forest = (DM_Forest*) dm->data;
291   PetscInt       dim, dimEmbed;
292   PetscErrorCode ierr;
293 
294   PetscFunctionBegin;
295   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
296   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the base after setup");
297   ierr         = PetscObjectReference((PetscObject)base);CHKERRQ(ierr);
298   ierr         = DMDestroy(&forest->base);CHKERRQ(ierr);
299   forest->base = base;
300   if (base) {
301     PetscValidHeaderSpecific(base, DM_CLASSID, 2);
302     ierr = DMGetDimension(base,&dim);CHKERRQ(ierr);
303     ierr = DMSetDimension(dm,dim);CHKERRQ(ierr);
304     ierr = DMGetCoordinateDim(base,&dimEmbed);CHKERRQ(ierr);
305     ierr = DMSetCoordinateDim(dm,dimEmbed);CHKERRQ(ierr);
306   }
307   PetscFunctionReturn(0);
308 }
309 
310 /*@
311   DMForestGetBaseDM - Get the base DM of a DMForest forest.  The forest will be hierarchically refined from the base,
312   and all refinements/coarsenings of the forest will share its base.  In general, two forest must share a bse to be
313   comparable, to do things like construct interpolators.
314 
315   Not collective
316 
317   Input Parameter:
318 . dm - the forest
319 
320   Output Parameter:
321 . base - the base DM of the forest
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 Parameter:
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   PetscErrorCode ierr;
384 
385   PetscFunctionBegin;
386   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
387   PetscValidHeaderSpecific(dm, DM_CLASSID, 2);
388   forest   = (DM_Forest*) dm->data;
389   ierr     = DMForestGetAdaptivityForest(dm,&oldAdapt);CHKERRQ(ierr);
390   if (adapt != NULL && dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
391   adaptForest    = (DM_Forest*) (adapt ? adapt->data : NULL);
392   oldAdaptForest = (DM_Forest*) (oldAdapt ? oldAdapt->data : NULL);
393   if (adaptForest != oldAdaptForest) {
394     ierr = PetscSFDestroy(&forest->preCoarseToFine);CHKERRQ(ierr);
395     ierr = PetscSFDestroy(&forest->coarseToPreFine);CHKERRQ(ierr);
396     if (forest->clearadaptivityforest) {ierr = (forest->clearadaptivityforest)(dm);CHKERRQ(ierr);}
397   }
398   switch (forest->adaptPurpose) {
399   case DM_ADAPT_DETERMINE:
400     ierr          = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr);
401     ierr          = DMDestroy(&(forest->adapt));CHKERRQ(ierr);
402     forest->adapt = adapt;
403     break;
404   case DM_ADAPT_REFINE:
405     ierr = DMSetCoarseDM(dm,adapt);CHKERRQ(ierr);
406     break;
407   case DM_ADAPT_COARSEN:
408     ierr = DMSetFineDM(dm,adapt);CHKERRQ(ierr);
409     break;
410   default:
411     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
412   }
413   PetscFunctionReturn(0);
414 }
415 
416 /*@
417   DMForestGetAdaptivityForest - Get the forest from which the current forest is adapted.
418 
419   Not collective
420 
421   Input Parameter:
422 . dm - the forest
423 
424   Output Parameter:
425 . adapt - the forest from which dm is/was adapted
426 
427   Level: intermediate
428 
429 .seealso: DMForestSetAdaptivityForest(), DMForestSetAdaptivityPurpose()
430 @*/
431 PetscErrorCode DMForestGetAdaptivityForest(DM dm, DM *adapt)
432 {
433   DM_Forest      *forest;
434   PetscErrorCode ierr;
435 
436   PetscFunctionBegin;
437   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
438   forest = (DM_Forest*) dm->data;
439   switch (forest->adaptPurpose) {
440   case DM_ADAPT_DETERMINE:
441     *adapt = forest->adapt;
442     break;
443   case DM_ADAPT_REFINE:
444     ierr = DMGetCoarseDM(dm,adapt);CHKERRQ(ierr);
445     break;
446   case DM_ADAPT_COARSEN:
447     ierr = DMGetFineDM(dm,adapt);CHKERRQ(ierr);
448     break;
449   default:
450     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"invalid adaptivity purpose");
451   }
452   PetscFunctionReturn(0);
453 }
454 
455 /*@
456   DMForestSetAdaptivityPurpose - During the pre-setup phase, set whether the current DM is being adapted from its
457   source (set with DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening
458   (DM_ADAPT_COARSEN), or undefined (DM_ADAPT_DETERMINE).  This only matters for the purposes of reference counting:
459   during DMDestroy(), cyclic references can be found between DMs only if the cyclic reference is due to a fine/coarse
460   relationship (see DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does
461   not maintain a reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can
462   cause a memory leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
463 
464   Logically collective on dm
465 
466   Input Parameters:
467 + dm - the forest
468 - purpose - the adaptivity purpose (DM_ADAPT_DETERMINE/DM_ADAPT_REFINE/DM_ADAPT_COARSEN)
469 
470   Level: advanced
471 
472 .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest()
473 @*/
474 PetscErrorCode DMForestSetAdaptivityPurpose(DM dm, DMAdaptFlag purpose)
475 {
476   DM_Forest      *forest;
477   PetscErrorCode ierr;
478 
479   PetscFunctionBegin;
480   forest = (DM_Forest*) dm->data;
481   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adaptation forest after setup");
482   if (purpose != forest->adaptPurpose) {
483     DM adapt;
484 
485     ierr = DMForestGetAdaptivityForest(dm,&adapt);CHKERRQ(ierr);
486     ierr = PetscObjectReference((PetscObject)adapt);CHKERRQ(ierr);
487     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
488 
489     forest->adaptPurpose = purpose;
490 
491     ierr = DMForestSetAdaptivityForest(dm,adapt);CHKERRQ(ierr);
492     ierr = DMDestroy(&adapt);CHKERRQ(ierr);
493   }
494   PetscFunctionReturn(0);
495 }
496 
497 /*@
498   DMForestGetAdaptivityPurpose - Get whether the current DM is being adapted from its source (set with
499   DMForestSetAdaptivityForest()) for the purpose of refinement (DM_ADAPT_REFINE), coarsening (DM_ADAPT_COARSEN), or
500   undefined (DM_ADAPT_DETERMINE).  This only matters for the purposes of reference counting: during DMDestroy(), cyclic
501   references can be found between DMs only if the cyclic reference is due to a fine/coarse relationship (see
502   DMSetFineDM()/DMSetCoarseDM()).  If the purpose is not refinement or coarsening, and the user does not maintain a
503   reference to the post-adaptation forest (i.e., the one created by DMForestTemplate()), then this can cause a memory
504   leak.  This method is used by subtypes of DMForest when automatically constructing mesh hierarchies.
505 
506   Not collective
507 
508   Input Parameter:
509 . dm - the forest
510 
511   Output Parameter:
512 . purpose - the adaptivity purpose (DM_ADAPT_DETERMINE/DM_ADAPT_REFINE/DM_ADAPT_COARSEN)
513 
514   Level: advanced
515 
516 .seealso: DMForestTemplate(), DMForestSetAdaptivityForest(), DMForestGetAdaptivityForest()
517 @*/
518 PetscErrorCode DMForestGetAdaptivityPurpose(DM dm, DMAdaptFlag *purpose)
519 {
520   DM_Forest *forest;
521 
522   PetscFunctionBegin;
523   forest   = (DM_Forest*) dm->data;
524   *purpose = forest->adaptPurpose;
525   PetscFunctionReturn(0);
526 }
527 
528 /*@
529   DMForestSetAdjacencyDimension - During the pre-setup phase, set the dimension of interface points that determine
530   cell adjacency (for the purposes of partitioning and overlap).
531 
532   Logically collective on dm
533 
534   Input Parameters:
535 + dm - the forest
536 - adjDim - default 0 (i.e., vertices determine adjacency)
537 
538   Level: intermediate
539 
540 .seealso: DMForestGetAdjacencyDimension(), DMForestSetAdjacencyCodimension(), DMForestSetPartitionOverlap()
541 @*/
542 PetscErrorCode DMForestSetAdjacencyDimension(DM dm, PetscInt adjDim)
543 {
544   PetscInt       dim;
545   DM_Forest      *forest = (DM_Forest*) dm->data;
546   PetscErrorCode ierr;
547 
548   PetscFunctionBegin;
549   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
550   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the adjacency dimension after setup");
551   if (adjDim < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be < 0: %d", adjDim);
552   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
553   if (adjDim > dim) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"adjacency dim cannot be > %d: %d", dim, adjDim);
554   forest->adjDim = adjDim;
555   PetscFunctionReturn(0);
556 }
557 
558 /*@
559   DMForestSetAdjacencyCodimension - Like DMForestSetAdjacencyDimension(), but specified as a co-dimension (so that,
560   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
561 
562   Logically collective on dm
563 
564   Input Parameters:
565 + dm - the forest
566 - adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
567 
568   Level: intermediate
569 
570 .seealso: DMForestGetAdjacencyCodimension(), DMForestSetAdjacencyDimension()
571 @*/
572 PetscErrorCode DMForestSetAdjacencyCodimension(DM dm, PetscInt adjCodim)
573 {
574   PetscInt       dim;
575   PetscErrorCode ierr;
576 
577   PetscFunctionBegin;
578   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
579   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
580   ierr = DMForestSetAdjacencyDimension(dm,dim-adjCodim);CHKERRQ(ierr);
581   PetscFunctionReturn(0);
582 }
583 
584 /*@
585   DMForestGetAdjacencyDimension - Get the dimension of interface points that determine cell adjacency (for the
586   purposes of partitioning and overlap).
587 
588   Not collective
589 
590   Input Parameter:
591 . dm - the forest
592 
593   Output Parameter:
594 . adjDim - default 0 (i.e., vertices determine adjacency)
595 
596   Level: intermediate
597 
598 .seealso: DMForestSetAdjacencyDimension(), DMForestGetAdjacencyCodimension(), DMForestSetPartitionOverlap()
599 @*/
600 PetscErrorCode DMForestGetAdjacencyDimension(DM dm, PetscInt *adjDim)
601 {
602   DM_Forest *forest = (DM_Forest*) dm->data;
603 
604   PetscFunctionBegin;
605   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
606   PetscValidIntPointer(adjDim,2);
607   *adjDim = forest->adjDim;
608   PetscFunctionReturn(0);
609 }
610 
611 /*@
612   DMForestGetAdjacencyCodimension - Like DMForestGetAdjacencyDimension(), but specified as a co-dimension (so that,
613   e.g., adjacency based on facets can be specified by codimension 1 in all cases)
614 
615   Not collective
616 
617   Input Parameter:
618 . dm - the forest
619 
620   Output Parameter:
621 . adjCodim - default isthe dimension of the forest (see DMGetDimension()), since this is the codimension of vertices
622 
623   Level: intermediate
624 
625 .seealso: DMForestSetAdjacencyCodimension(), DMForestGetAdjacencyDimension()
626 @*/
627 PetscErrorCode DMForestGetAdjacencyCodimension(DM dm, PetscInt *adjCodim)
628 {
629   DM_Forest      *forest = (DM_Forest*) dm->data;
630   PetscInt       dim;
631   PetscErrorCode ierr;
632 
633   PetscFunctionBegin;
634   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
635   PetscValidIntPointer(adjCodim,2);
636   ierr      = DMGetDimension(dm,&dim);CHKERRQ(ierr);
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   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the overlap after setup");
663   if (overlap < 0) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"overlap cannot be < 0: %d", 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   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the minimum refinement after setup");
717   forest->minRefinement = minRefinement;
718   PetscFunctionReturn(0);
719 }
720 
721 /*@
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   if (dm->setupcalled) SETERRQ(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 Paramater:
784 . initefinement - 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   if (dm->setupcalled) SETERRQ(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   PetscErrorCode ierr;
875 
876   PetscFunctionBegin;
877   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
878   ierr = PetscFree(forest->adaptStrategy);CHKERRQ(ierr);
879   ierr = PetscStrallocpy((const char*) adaptStrategy,(char**)&forest->adaptStrategy);CHKERRQ(ierr);
880   PetscFunctionReturn(0);
881 }
882 
883 /*@C
884   DMForestSetAdaptivityStrategy - Get the strategy for combining adaptivity labels from multiple processes.  Subtypes
885   of DMForest may define their own strategies.  Two default strategies are DMFORESTADAPTALL, which indicates that all
886   processes must agree for a refinement/coarsening flag to be valid, and DMFORESTADAPTANY, which indicates that only
887   one process needs to specify refinement/coarsening.
888 
889   Not collective
890 
891   Input Parameter:
892 . dm - the forest
893 
894   Output Parameter:
895 . adaptStrategy - the adaptivity strategy (default DMFORESTADAPTALL)
896 
897   Level: advanced
898 
899 .seealso: DMForestSetAdaptivityStrategy()
900 @*/
901 PetscErrorCode DMForestGetAdaptivityStrategy(DM dm, DMForestAdaptivityStrategy *adaptStrategy)
902 {
903   DM_Forest *forest = (DM_Forest*) dm->data;
904 
905   PetscFunctionBegin;
906   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
907   PetscValidPointer(adaptStrategy,2);
908   *adaptStrategy = forest->adaptStrategy;
909   PetscFunctionReturn(0);
910 }
911 
912 /*@
913   DMForestGetAdaptivitySuccess - Return whether the requested adaptation (refinement, coarsening, repartitioning,
914   etc.) was successful.  PETSC_FALSE indicates that the post-adaptation forest is the same as the pre-adpatation
915   forest.  A requested adaptation may have been unsuccessful if, for example, the requested refinement would have
916   exceeded the maximum refinement level.
917 
918   Collective on dm
919 
920   Input Parameter:
921 
922 . dm - the post-adaptation forest
923 
924   Output Parameter:
925 
926 . success - PETSC_TRUE if the post-adaptation forest is different from the pre-adaptation forest.
927 
928   Level: intermediate
929 
930 .see
931 @*/
932 PetscErrorCode DMForestGetAdaptivitySuccess(DM dm, PetscBool *success)
933 {
934   DM_Forest      *forest;
935   PetscErrorCode ierr;
936 
937   PetscFunctionBegin;
938   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
939   if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"DMSetUp() has not been called yet.");
940   forest = (DM_Forest *) dm->data;
941   ierr = (forest->getadaptivitysuccess)(dm,success);CHKERRQ(ierr);
942   PetscFunctionReturn(0);
943 }
944 
945 /*@
946   DMForestSetComputeAdaptivitySF - During the pre-setup phase, set whether transfer PetscSFs should be computed
947   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().
948 
949   Logically collective on dm
950 
951   Input Parameters:
952 + dm - the post-adaptation forest
953 - computeSF - default PETSC_TRUE
954 
955   Level: advanced
956 
957 .seealso: DMForestGetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
958 @*/
959 PetscErrorCode DMForestSetComputeAdaptivitySF(DM dm, PetscBool computeSF)
960 {
961   DM_Forest *forest;
962 
963   PetscFunctionBegin;
964   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
965   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot compute adaptivity PetscSFs after setup is called");
966   forest                 = (DM_Forest*) dm->data;
967   forest->computeAdaptSF = computeSF;
968   PetscFunctionReturn(0);
969 }
970 
971 PetscErrorCode DMForestTransferVec(DM dmIn, Vec vecIn, DM dmOut, Vec vecOut, PetscBool useBCs, PetscReal time)
972 {
973   DM_Forest      *forest;
974   PetscErrorCode ierr;
975 
976   PetscFunctionBegin;
977   PetscValidHeaderSpecific(dmIn   ,DM_CLASSID  ,1);
978   PetscValidHeaderSpecific(vecIn  ,VEC_CLASSID ,2);
979   PetscValidHeaderSpecific(dmOut  ,DM_CLASSID  ,3);
980   PetscValidHeaderSpecific(vecOut ,VEC_CLASSID ,4);
981   forest = (DM_Forest *) dmIn->data;
982   if (!forest->transfervec) SETERRQ(PetscObjectComm((PetscObject)dmIn),PETSC_ERR_SUP,"DMForestTransferVec() not implemented");
983   ierr = (forest->transfervec)(dmIn,vecIn,dmOut,vecOut,useBCs,time);CHKERRQ(ierr);
984   PetscFunctionReturn(0);
985 }
986 
987 /*@
988   DMForestGetComputeAdaptivitySF - Get whether transfer PetscSFs should be computed relating the cells of the
989   pre-adaptation forest to the post-adaptiation forest.  After DMSetUp() is called, these transfer PetscSFs can be
990   accessed with DMForestGetAdaptivitySF().
991 
992   Not collective
993 
994   Input Parameter:
995 . dm - the post-adaptation forest
996 
997   Output Parameter:
998 . computeSF - default PETSC_TRUE
999 
1000   Level: advanced
1001 
1002 .seealso: DMForestSetComputeAdaptivitySF(), DMForestGetAdaptivitySF()
1003 @*/
1004 PetscErrorCode DMForestGetComputeAdaptivitySF(DM dm, PetscBool *computeSF)
1005 {
1006   DM_Forest *forest;
1007 
1008   PetscFunctionBegin;
1009   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1010   forest     = (DM_Forest*) dm->data;
1011   *computeSF = forest->computeAdaptSF;
1012   PetscFunctionReturn(0);
1013 }
1014 
1015 /*@
1016   DMForestGetAdaptivitySF - Get PetscSFs that relate the pre-adaptation forest to the post-adaptation forest.
1017   Adaptation can be any combination of refinement, coarsening, repartition, and change of overlap, so there may be
1018   some cells of the pre-adaptation that are parents of post-adaptation cells, and vice versa.  Therefore there are two
1019   PetscSFs: one that relates pre-adaptation coarse cells to post-adaptation fine cells, and one that relates
1020   pre-adaptation fine cells to post-adaptation coarse cells.
1021 
1022   Not collective
1023 
1024   Input Parameter:
1025   dm - the post-adaptation forest
1026 
1027   Output Parameter:
1028   preCoarseToFine - pre-adaptation coarse cells to post-adaptation fine cells: BCast goes from pre- to post-
1029   coarseToPreFine - post-adaptation coarse cells to pre-adaptation fine cells: BCast goes from post- to pre-
1030 
1031   Level: advanced
1032 
1033 .seealso: DMForestGetComputeAdaptivitySF(), DMForestSetComputeAdaptivitySF()
1034 @*/
1035 PetscErrorCode DMForestGetAdaptivitySF(DM dm, PetscSF *preCoarseToFine, PetscSF *coarseToPreFine)
1036 {
1037   DM_Forest      *forest;
1038   PetscErrorCode ierr;
1039 
1040   PetscFunctionBegin;
1041   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1042   ierr   = DMSetUp(dm);CHKERRQ(ierr);
1043   forest = (DM_Forest*) dm->data;
1044   if (preCoarseToFine) *preCoarseToFine = forest->preCoarseToFine;
1045   if (coarseToPreFine) *coarseToPreFine = forest->coarseToPreFine;
1046   PetscFunctionReturn(0);
1047 }
1048 
1049 /*@
1050   DMForestSetGradeFactor - During the pre-setup phase, set the desired amount of grading in the mesh, e.g. give 2 to
1051   indicate that the diameter of neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may
1052   only support one particular choice of grading factor.
1053 
1054   Logically collective on dm
1055 
1056   Input Parameters:
1057 + dm - the forest
1058 - grade - the grading factor
1059 
1060   Level: advanced
1061 
1062 .seealso: DMForestGetGradeFactor()
1063 @*/
1064 PetscErrorCode DMForestSetGradeFactor(DM dm, PetscInt grade)
1065 {
1066   DM_Forest *forest = (DM_Forest*) dm->data;
1067 
1068   PetscFunctionBegin;
1069   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1070   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the grade factor after setup");
1071   forest->gradeFactor = grade;
1072   PetscFunctionReturn(0);
1073 }
1074 
1075 /*@
1076   DMForestGetGradeFactor - Get the desired amount of grading in the mesh, e.g. give 2 to indicate that the diameter of
1077   neighboring cells should differ by at most a factor of 2.  Subtypes of DMForest may only support one particular
1078   choice of grading factor.
1079 
1080   Not collective
1081 
1082   Input Parameter:
1083 . dm - the forest
1084 
1085   Output Parameter:
1086 . grade - the grading factor
1087 
1088   Level: advanced
1089 
1090 .seealso: DMForestSetGradeFactor()
1091 @*/
1092 PetscErrorCode DMForestGetGradeFactor(DM dm, PetscInt *grade)
1093 {
1094   DM_Forest *forest = (DM_Forest*) dm->data;
1095 
1096   PetscFunctionBegin;
1097   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1098   PetscValidIntPointer(grade,2);
1099   *grade = forest->gradeFactor;
1100   PetscFunctionReturn(0);
1101 }
1102 
1103 /*@
1104   DMForestSetCellWeightFactor - During the pre-setup phase, set the factor by which the level of refinement changes
1105   the cell weight (see DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be
1106   (cellWeight) * (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on
1107   its level; a factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the
1108   computation associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
1109 
1110   Logically collective on dm
1111 
1112   Input Parameters:
1113 + dm - the forest
1114 - weightsFactors - default 1.
1115 
1116   Level: advanced
1117 
1118 .seealso: DMForestGetCellWeightFactor(), DMForestSetCellWeights()
1119 @*/
1120 PetscErrorCode DMForestSetCellWeightFactor(DM dm, PetscReal weightsFactor)
1121 {
1122   DM_Forest *forest = (DM_Forest*) dm->data;
1123 
1124   PetscFunctionBegin;
1125   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1126   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weights factor after setup");
1127   forest->weightsFactor = weightsFactor;
1128   PetscFunctionReturn(0);
1129 }
1130 
1131 /*@
1132   DMForestGetCellWeightFactor - Get the factor by which the level of refinement changes the cell weight (see
1133   DMForestSetCellWeights()) when calculating partitions.  The final weight of a cell will be (cellWeight) *
1134   (weightFactor^refinementLevel).  A factor of 1 indicates that the weight of a cell does not depend on its level; a
1135   factor of 2, for example, might be appropriate for sub-cycling time-stepping methods, when the computation
1136   associated with a cell is multiplied by a factor of 2 for each additional level of refinement.
1137 
1138   Not collective
1139 
1140   Input Parameter:
1141 . dm - the forest
1142 
1143   Output Parameter:
1144 . weightsFactors - default 1.
1145 
1146   Level: advanced
1147 
1148 .seealso: DMForestSetCellWeightFactor(), DMForestSetCellWeights()
1149 @*/
1150 PetscErrorCode DMForestGetCellWeightFactor(DM dm, PetscReal *weightsFactor)
1151 {
1152   DM_Forest *forest = (DM_Forest*) dm->data;
1153 
1154   PetscFunctionBegin;
1155   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1156   PetscValidRealPointer(weightsFactor,2);
1157   *weightsFactor = forest->weightsFactor;
1158   PetscFunctionReturn(0);
1159 }
1160 
1161 /*@
1162   DMForestGetCellChart - After the setup phase, get the local half-open interval of the chart of cells on this process
1163 
1164   Not collective
1165 
1166   Input Parameter:
1167 . dm - the forest
1168 
1169   Output Parameters:
1170 + cStart - the first cell on this process
1171 - cEnd - one after the final cell on this process
1172 
1173   Level: intermediate
1174 
1175 .seealso: DMForestGetCellSF()
1176 @*/
1177 PetscErrorCode DMForestGetCellChart(DM dm, PetscInt *cStart, PetscInt *cEnd)
1178 {
1179   DM_Forest      *forest = (DM_Forest*) dm->data;
1180   PetscErrorCode ierr;
1181 
1182   PetscFunctionBegin;
1183   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1184   PetscValidIntPointer(cStart,2);
1185   PetscValidIntPointer(cEnd,2);
1186   if (((forest->cStart == PETSC_DETERMINE) || (forest->cEnd == PETSC_DETERMINE)) && forest->createcellchart) {
1187     ierr = forest->createcellchart(dm,&forest->cStart,&forest->cEnd);CHKERRQ(ierr);
1188   }
1189   *cStart =  forest->cStart;
1190   *cEnd   =  forest->cEnd;
1191   PetscFunctionReturn(0);
1192 }
1193 
1194 /*@
1195   DMForestGetCellSF - After the setup phase, get the PetscSF for overlapping cells between processes
1196 
1197   Not collective
1198 
1199   Input Parameter:
1200 . dm - the forest
1201 
1202   Output Parameter:
1203 . cellSF - the PetscSF
1204 
1205   Level: intermediate
1206 
1207 .seealso: DMForestGetCellChart()
1208 @*/
1209 PetscErrorCode DMForestGetCellSF(DM dm, PetscSF *cellSF)
1210 {
1211   DM_Forest      *forest = (DM_Forest*) dm->data;
1212   PetscErrorCode ierr;
1213 
1214   PetscFunctionBegin;
1215   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1216   PetscValidPointer(cellSF,2);
1217   if ((!forest->cellSF) && forest->createcellsf) {
1218     ierr = forest->createcellsf(dm,&forest->cellSF);CHKERRQ(ierr);
1219   }
1220   *cellSF = forest->cellSF;
1221   PetscFunctionReturn(0);
1222 }
1223 
1224 /*@C
1225   DMForestSetAdaptivityLabel - During the pre-setup phase, set the label of the pre-adaptation forest (see
1226   DMForestGetAdaptivityForest()) that holds the adaptation flags (refinement, coarsening, or some combination).  The
1227   interpretation of the label values is up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP,
1228   DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have been reserved as choices that should be accepted by all subtypes.
1229 
1230   Logically collective on dm
1231 
1232   Input Parameters:
1233 - dm - the forest
1234 + adaptLabel - the label in the pre-adaptation forest
1235 
1236   Level: intermediate
1237 
1238 .seealso DMForestGetAdaptivityLabel()
1239 @*/
1240 PetscErrorCode DMForestSetAdaptivityLabel(DM dm, DMLabel adaptLabel)
1241 {
1242   DM_Forest      *forest = (DM_Forest*) dm->data;
1243   PetscErrorCode ierr;
1244 
1245   PetscFunctionBegin;
1246   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1247   adaptLabel->refct++;
1248   if (forest->adaptLabel) {ierr = DMLabelDestroy(&forest->adaptLabel);CHKERRQ(ierr);}
1249   forest->adaptLabel = adaptLabel;
1250   PetscFunctionReturn(0);
1251 }
1252 
1253 /*@C
1254   DMForestGetAdaptivityLabel - Get the label of the pre-adaptation forest (see DMForestGetAdaptivityForest()) that
1255   holds the adaptation flags (refinement, coarsening, or some combination).  The interpretation of the label values is
1256   up to the subtype of DMForest, but DM_ADAPT_DETERMINE, DM_ADAPT_KEEP, DM_ADAPT_REFINE, and DM_ADAPT_COARSEN have
1257   been reserved as choices that should be accepted by all subtypes.
1258 
1259   Not collective
1260 
1261   Input Parameter:
1262 . dm - the forest
1263 
1264   Output Parameter:
1265 . adaptLabel - the name of the label in the pre-adaptation forest
1266 
1267   Level: intermediate
1268 
1269 .seealso DMForestSetAdaptivityLabel()
1270 @*/
1271 PetscErrorCode DMForestGetAdaptivityLabel(DM dm, DMLabel *adaptLabel)
1272 {
1273   DM_Forest *forest = (DM_Forest*) dm->data;
1274 
1275   PetscFunctionBegin;
1276   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1277   *adaptLabel = forest->adaptLabel;
1278   PetscFunctionReturn(0);
1279 }
1280 
1281 /*@
1282   DMForestSetCellWeights - Set the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
1283   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
1284   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
1285   of 1.
1286 
1287   Logically collective on dm
1288 
1289   Input Parameters:
1290 + dm - the forest
1291 . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
1292 - copyMode - how weights should reference weights
1293 
1294   Level: advanced
1295 
1296 .seealso: DMForestGetCellWeights(), DMForestSetWeightCapacity()
1297 @*/
1298 PetscErrorCode DMForestSetCellWeights(DM dm, PetscReal weights[], PetscCopyMode copyMode)
1299 {
1300   DM_Forest      *forest = (DM_Forest*) dm->data;
1301   PetscInt       cStart, cEnd;
1302   PetscErrorCode ierr;
1303 
1304   PetscFunctionBegin;
1305   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1306   ierr = DMForestGetCellChart(dm,&cStart,&cEnd);CHKERRQ(ierr);
1307   if (cEnd < cStart) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"cell chart [%d,%d) is not valid",cStart,cEnd);
1308   if (copyMode == PETSC_COPY_VALUES) {
1309     if (forest->cellWeightsCopyMode != PETSC_OWN_POINTER || forest->cellWeights == weights) {
1310       ierr = PetscMalloc1(cEnd-cStart,&forest->cellWeights);CHKERRQ(ierr);
1311     }
1312     ierr                        = PetscMemcpy(forest->cellWeights,weights,(cEnd-cStart)*sizeof(*weights));CHKERRQ(ierr);
1313     forest->cellWeightsCopyMode = PETSC_OWN_POINTER;
1314     PetscFunctionReturn(0);
1315   }
1316   if (forest->cellWeightsCopyMode == PETSC_OWN_POINTER) {
1317     ierr = PetscFree(forest->cellWeights);CHKERRQ(ierr);
1318   }
1319   forest->cellWeights         = weights;
1320   forest->cellWeightsCopyMode = copyMode;
1321   PetscFunctionReturn(0);
1322 }
1323 
1324 /*@
1325   DMForestGetCellWeights - Get the weights assigned to each of the cells (see DMForestGetCellChart()) of the current
1326   process: weights are used to determine parallel partitioning.  Partitions will be created so that each process's
1327   ratio of weight to capacity (see DMForestSetWeightCapacity()) is roughly equal. If NULL, each cell receives a weight
1328   of 1.
1329 
1330   Not collective
1331 
1332   Input Parameter:
1333 . dm - the forest
1334 
1335   Output Parameter:
1336 . weights - the array of weights for all cells, or NULL to indicate each cell has weight 1.
1337 
1338   Level: advanced
1339 
1340 .seealso: DMForestSetCellWeights(), DMForestSetWeightCapacity()
1341 @*/
1342 PetscErrorCode DMForestGetCellWeights(DM dm, PetscReal **weights)
1343 {
1344   DM_Forest *forest = (DM_Forest*) dm->data;
1345 
1346   PetscFunctionBegin;
1347   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1348   PetscValidPointer(weights,2);
1349   *weights = forest->cellWeights;
1350   PetscFunctionReturn(0);
1351 }
1352 
1353 /*@
1354   DMForestSetWeightCapacity - During the pre-setup phase, set the capacity of the current process when repartitioning
1355   a pre-adaptation forest (see DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each
1356   process's cells to the process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that
1357   the current process should not have any cells after repartitioning.
1358 
1359   Logically Collective on dm
1360 
1361   Input parameters:
1362 + dm - the forest
1363 - capacity - this process's capacity
1364 
1365   Level: advanced
1366 
1367 .seealso DMForestGetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
1368 @*/
1369 PetscErrorCode DMForestSetWeightCapacity(DM dm, PetscReal capacity)
1370 {
1371   DM_Forest *forest = (DM_Forest*) dm->data;
1372 
1373   PetscFunctionBegin;
1374   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1375   if (dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Cannot change the weight capacity after setup");
1376   if (capacity < 0.) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Cannot have negative weight capacity; %f",capacity);
1377   forest->weightCapacity = capacity;
1378   PetscFunctionReturn(0);
1379 }
1380 
1381 /*@
1382   DMForestGetWeightCapacity - Set the capacity of the current process when repartitioning a pre-adaptation forest (see
1383   DMForestGetAdaptivityForest()).  After partitioning, the ratio of the weight of each process's cells to the
1384   process's capacity will be roughly equal for all processes.  A capacity of 0 indicates that the current process
1385   should not have any cells after repartitioning.
1386 
1387   Not collective
1388 
1389   Input parameter:
1390 . dm - the forest
1391 
1392   Output parameter:
1393 . capacity - this process's capacity
1394 
1395   Level: advanced
1396 
1397 .seealso DMForestSetWeightCapacity(), DMForestSetCellWeights(), DMForestSetCellWeightFactor()
1398 @*/
1399 PetscErrorCode DMForestGetWeightCapacity(DM dm, PetscReal *capacity)
1400 {
1401   DM_Forest *forest = (DM_Forest*) dm->data;
1402 
1403   PetscFunctionBegin;
1404   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
1405   PetscValidRealPointer(capacity,2);
1406   *capacity = forest->weightCapacity;
1407   PetscFunctionReturn(0);
1408 }
1409 
1410 PETSC_EXTERN PetscErrorCode DMSetFromOptions_Forest(PetscOptionItems *PetscOptionsObject,DM dm)
1411 {
1412   DM_Forest                  *forest = (DM_Forest*) dm->data;
1413   PetscBool                  flg, flg1, flg2, flg3, flg4;
1414   DMForestTopology           oldTopo;
1415   char                       stringBuffer[256];
1416   PetscViewer                viewer;
1417   PetscViewerFormat          format;
1418   PetscInt                   adjDim, adjCodim, overlap, minRefinement, initRefinement, maxRefinement, grade;
1419   PetscReal                  weightsFactor;
1420   DMForestAdaptivityStrategy adaptStrategy;
1421   PetscErrorCode             ierr;
1422 
1423   PetscFunctionBegin;
1424   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1425   forest->setfromoptionscalled = PETSC_TRUE;
1426   ierr                         = DMForestGetTopology(dm, &oldTopo);CHKERRQ(ierr);
1427   ierr                         = PetscOptionsHead(PetscOptionsObject,"DMForest Options");CHKERRQ(ierr);
1428   ierr                         = PetscOptionsString("-dm_forest_topology","the topology of the forest's base mesh","DMForestSetTopology",oldTopo,stringBuffer,256,&flg1);CHKERRQ(ierr);
1429   ierr                         = PetscOptionsViewer("-dm_forest_base_dm","load the base DM from a viewer specification","DMForestSetBaseDM",&viewer,&format,&flg2);CHKERRQ(ierr);
1430   ierr                         = PetscOptionsViewer("-dm_forest_coarse_forest","load the coarse forest from a viewer specification","DMForestSetCoarseForest",&viewer,&format,&flg3);CHKERRQ(ierr);
1431   ierr                         = PetscOptionsViewer("-dm_forest_fine_forest","load the fine forest from a viewer specification","DMForestSetFineForest",&viewer,&format,&flg4);CHKERRQ(ierr);
1432   if ((PetscInt) flg1 + (PetscInt) flg2 + (PetscInt) flg3 + (PetscInt) flg4 > 1) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Specify only one of -dm_forest_{topology,base_dm,coarse_forest,fine_forest}");
1433   if (flg1) {
1434     ierr = DMForestSetTopology(dm,(DMForestTopology)stringBuffer);CHKERRQ(ierr);
1435     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1436     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
1437   }
1438   if (flg2) {
1439     DM base;
1440 
1441     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&base);CHKERRQ(ierr);
1442     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1443     ierr = DMLoad(base,viewer);CHKERRQ(ierr);
1444     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1445     ierr = DMForestSetBaseDM(dm,base);CHKERRQ(ierr);
1446     ierr = DMDestroy(&base);CHKERRQ(ierr);
1447     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
1448     ierr = DMForestSetAdaptivityForest(dm,NULL);CHKERRQ(ierr);
1449   }
1450   if (flg3) {
1451     DM coarse;
1452 
1453     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&coarse);CHKERRQ(ierr);
1454     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1455     ierr = DMLoad(coarse,viewer);CHKERRQ(ierr);
1456     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1457     ierr = DMForestSetAdaptivityForest(dm,coarse);CHKERRQ(ierr);
1458     ierr = DMDestroy(&coarse);CHKERRQ(ierr);
1459     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
1460     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1461   }
1462   if (flg4) {
1463     DM fine;
1464 
1465     ierr = DMCreate(PetscObjectComm((PetscObject)dm),&fine);CHKERRQ(ierr);
1466     ierr = PetscViewerPushFormat(viewer,format);CHKERRQ(ierr);
1467     ierr = DMLoad(fine,viewer);CHKERRQ(ierr);
1468     ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1469     ierr = DMForestSetAdaptivityForest(dm,fine);CHKERRQ(ierr);
1470     ierr = DMDestroy(&fine);CHKERRQ(ierr);
1471     ierr = DMForestSetTopology(dm,NULL);CHKERRQ(ierr);
1472     ierr = DMForestSetBaseDM(dm,NULL);CHKERRQ(ierr);
1473   }
1474   ierr = DMForestGetAdjacencyDimension(dm,&adjDim);CHKERRQ(ierr);
1475   ierr = PetscOptionsInt("-dm_forest_adjacency_dimension","set the dimension of points that define adjacency in the forest","DMForestSetAdjacencyDimension",adjDim,&adjDim,&flg);CHKERRQ(ierr);
1476   if (flg) {
1477     ierr = DMForestSetAdjacencyDimension(dm,adjDim);CHKERRQ(ierr);
1478   } else {
1479     ierr = DMForestGetAdjacencyCodimension(dm,&adjCodim);CHKERRQ(ierr);
1480     ierr = PetscOptionsInt("-dm_forest_adjacency_codimension","set the codimension of points that define adjacency in the forest","DMForestSetAdjacencyCodimension",adjCodim,&adjCodim,&flg);CHKERRQ(ierr);
1481     if (flg) {
1482       ierr = DMForestSetAdjacencyCodimension(dm,adjCodim);CHKERRQ(ierr);
1483     }
1484   }
1485   ierr = DMForestGetPartitionOverlap(dm,&overlap);CHKERRQ(ierr);
1486   ierr = PetscOptionsInt("-dm_forest_partition_overlap","set the degree of partition overlap","DMForestSetPartitionOverlap",overlap,&overlap,&flg);CHKERRQ(ierr);
1487   if (flg) {
1488     ierr = DMForestSetPartitionOverlap(dm,overlap);CHKERRQ(ierr);
1489   }
1490 #if 0
1491   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);
1492   if (flg) {
1493     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1494     ierr = DMForestSetInitialRefinement(dm,minRefinement);CHKERRQ(ierr);
1495   }
1496   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);
1497   if (flg) {
1498     ierr = DMForestSetMinimumRefinement(dm,0);CHKERRQ(ierr);
1499     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
1500   }
1501 #endif
1502   ierr = DMForestGetMinimumRefinement(dm,&minRefinement);CHKERRQ(ierr);
1503   ierr = PetscOptionsInt("-dm_forest_minimum_refinement","set the minimum level of refinement in the forest","DMForestSetMinimumRefinement",minRefinement,&minRefinement,&flg);CHKERRQ(ierr);
1504   if (flg) {
1505     ierr = DMForestSetMinimumRefinement(dm,minRefinement);CHKERRQ(ierr);
1506   }
1507   ierr = DMForestGetInitialRefinement(dm,&initRefinement);CHKERRQ(ierr);
1508   ierr = PetscOptionsInt("-dm_forest_initial_refinement","set the initial level of refinement in the forest","DMForestSetInitialRefinement",initRefinement,&initRefinement,&flg);CHKERRQ(ierr);
1509   if (flg) {
1510     ierr = DMForestSetInitialRefinement(dm,initRefinement);CHKERRQ(ierr);
1511   }
1512   ierr = DMForestGetMaximumRefinement(dm,&maxRefinement);CHKERRQ(ierr);
1513   ierr = PetscOptionsInt("-dm_forest_maximum_refinement","set the maximum level of refinement in the forest","DMForestSetMaximumRefinement",maxRefinement,&maxRefinement,&flg);CHKERRQ(ierr);
1514   if (flg) {
1515     ierr = DMForestSetMaximumRefinement(dm,maxRefinement);CHKERRQ(ierr);
1516   }
1517   ierr = DMForestGetAdaptivityStrategy(dm,&adaptStrategy);CHKERRQ(ierr);
1518   ierr = PetscOptionsString("-dm_forest_adaptivity_strategy","the forest's adaptivity-flag resolution strategy","DMForestSetAdaptivityStrategy",adaptStrategy,stringBuffer,256,&flg);CHKERRQ(ierr);
1519   if (flg) {
1520     ierr = DMForestSetAdaptivityStrategy(dm,(DMForestAdaptivityStrategy)stringBuffer);CHKERRQ(ierr);
1521   }
1522   ierr = DMForestGetGradeFactor(dm,&grade);CHKERRQ(ierr);
1523   ierr = PetscOptionsInt("-dm_forest_grade_factor","grade factor between neighboring cells","DMForestSetGradeFactor",grade,&grade,&flg);CHKERRQ(ierr);
1524   if (flg) {
1525     ierr = DMForestSetGradeFactor(dm,grade);CHKERRQ(ierr);
1526   }
1527   ierr = DMForestGetCellWeightFactor(dm,&weightsFactor);CHKERRQ(ierr);
1528   ierr = PetscOptionsReal("-dm_forest_cell_weight_factor","multiplying weight factor for cell refinement","DMForestSetCellWeightFactor",weightsFactor,&weightsFactor,&flg);CHKERRQ(ierr);
1529   if (flg) {
1530     ierr = DMForestSetCellWeightFactor(dm,weightsFactor);CHKERRQ(ierr);
1531   }
1532   ierr = PetscOptionsTail();CHKERRQ(ierr);
1533   PetscFunctionReturn(0);
1534 }
1535 
1536 PetscErrorCode DMCreateSubDM_Forest(DM dm, PetscInt numFields, PetscInt fields[], IS *is, DM *subdm)
1537 {
1538   PetscErrorCode ierr;
1539 
1540   PetscFunctionBegin;
1541   if (subdm) {ierr = DMClone(dm, subdm);CHKERRQ(ierr);}
1542   ierr = DMCreateSubDM_Section_Private(dm, numFields, fields, is, subdm);CHKERRQ(ierr);
1543   PetscFunctionReturn(0);
1544 }
1545 
1546 PetscErrorCode DMRefine_Forest(DM dm, MPI_Comm comm, DM *dmRefined)
1547 {
1548   DMLabel        refine;
1549   DM             fineDM;
1550   PetscErrorCode ierr;
1551 
1552   PetscFunctionBegin;
1553   ierr = DMGetFineDM(dm,&fineDM);CHKERRQ(ierr);
1554   if (fineDM) {
1555     ierr       = PetscObjectReference((PetscObject)fineDM);CHKERRQ(ierr);
1556     *dmRefined = fineDM;
1557     PetscFunctionReturn(0);
1558   }
1559   ierr = DMForestTemplate(dm,comm,dmRefined);CHKERRQ(ierr);
1560   ierr = DMGetLabel(dm,"refine",&refine);CHKERRQ(ierr);
1561   if (!refine) {
1562     ierr = DMLabelCreate("refine",&refine);CHKERRQ(ierr);
1563     ierr = DMLabelSetDefaultValue(refine,DM_ADAPT_REFINE);CHKERRQ(ierr);
1564   }
1565   else {
1566     refine->refct++;
1567   }
1568   ierr = DMForestSetAdaptivityLabel(*dmRefined,refine);CHKERRQ(ierr);
1569   ierr = DMLabelDestroy(&refine);CHKERRQ(ierr);
1570   PetscFunctionReturn(0);
1571 }
1572 
1573 PetscErrorCode DMCoarsen_Forest(DM dm, MPI_Comm comm, DM *dmCoarsened)
1574 {
1575   DMLabel        coarsen;
1576   DM             coarseDM;
1577   PetscErrorCode ierr;
1578 
1579   PetscFunctionBegin;
1580   {
1581     PetscMPIInt mpiComparison;
1582     MPI_Comm    dmcomm = PetscObjectComm((PetscObject)dm);
1583 
1584     ierr = MPI_Comm_compare(comm, dmcomm, &mpiComparison);CHKERRQ(ierr);
1585     if (mpiComparison != MPI_IDENT && mpiComparison != MPI_CONGRUENT) SETERRQ(dmcomm,PETSC_ERR_SUP,"No support for different communicators yet");
1586   }
1587   ierr = DMGetCoarseDM(dm,&coarseDM);CHKERRQ(ierr);
1588   if (coarseDM) {
1589     ierr         = PetscObjectReference((PetscObject)coarseDM);CHKERRQ(ierr);
1590     *dmCoarsened = coarseDM;
1591     PetscFunctionReturn(0);
1592   }
1593   ierr = DMForestTemplate(dm,comm,dmCoarsened);CHKERRQ(ierr);
1594   ierr = DMForestSetAdaptivityPurpose(coarseDM,DM_ADAPT_COARSEN);CHKERRQ(ierr);
1595   ierr = DMGetLabel(dm,"coarsen",&coarsen);CHKERRQ(ierr);
1596   if (!coarsen) {
1597     ierr = DMLabelCreate("coarsen",&coarsen);CHKERRQ(ierr);
1598     ierr = DMLabelSetDefaultValue(coarsen,DM_ADAPT_COARSEN);CHKERRQ(ierr);
1599   } else {
1600     coarsen->refct++;
1601   }
1602   ierr = DMForestSetAdaptivityLabel(*dmCoarsened,coarsen);CHKERRQ(ierr);
1603   ierr = DMLabelDestroy(&coarsen);CHKERRQ(ierr);
1604   PetscFunctionReturn(0);
1605 }
1606 
1607 static PetscErrorCode DMAdaptLabel_Forest(DM dm, DMLabel label, DM *adaptedDM)
1608 {
1609   PetscBool      success;
1610   PetscErrorCode ierr;
1611 
1612   PetscFunctionBegin;
1613   ierr = DMForestTemplate(dm,PetscObjectComm((PetscObject)dm),adaptedDM);CHKERRQ(ierr);
1614   ierr = DMForestSetAdaptivityLabel(*adaptedDM,label);CHKERRQ(ierr);
1615   ierr = DMSetUp(*adaptedDM);CHKERRQ(ierr);
1616   ierr = DMForestGetAdaptivitySuccess(*adaptedDM,&success);CHKERRQ(ierr);
1617   if (!success) {
1618     ierr = DMDestroy(adaptedDM);CHKERRQ(ierr);
1619     *adaptedDM = NULL;
1620   }
1621   PetscFunctionReturn(0);
1622 }
1623 
1624 static PetscErrorCode DMInitialize_Forest(DM dm)
1625 {
1626   PetscErrorCode ierr;
1627 
1628   PetscFunctionBegin;
1629   ierr = PetscMemzero(dm->ops,sizeof(*(dm->ops)));CHKERRQ(ierr);
1630 
1631   dm->ops->clone          = DMClone_Forest;
1632   dm->ops->setfromoptions = DMSetFromOptions_Forest;
1633   dm->ops->destroy        = DMDestroy_Forest;
1634   dm->ops->createsubdm    = DMCreateSubDM_Forest;
1635   dm->ops->refine         = DMRefine_Forest;
1636   dm->ops->coarsen        = DMCoarsen_Forest;
1637   ierr                    = PetscObjectComposeFunction((PetscObject)dm,"DMAdaptLabel_C",DMAdaptLabel_Forest);CHKERRQ(ierr);
1638   PetscFunctionReturn(0);
1639 }
1640 
1641 /*MC
1642   DMFOREST = "forest" - A DM object that encapsulates a hierarchically refined mesh.  Forests usually have a base DM (see DMForestGetBaseDM()), from which it is refined.  The refinement and partitioning of forests is considered immutable after DMSetUp() is called.  To adapt a mesh, one should call DMForestTemplate() to create a new mesh that will default to being identical to it, specify how that mesh should differ, and then calling DMSetUp() on the new mesh.
1643 
1644   To specify that a mesh should be refined or coarsened from the previous mesh, a label should be defined on the previous mesh whose values indicate which cells should be refined (DM_ADAPT_REFINE) or coarsened (DM_ADAPT_COARSEN) and how (subtypes are free to allow additional values for things like anisotropic refinement).  The label should be given to the *new* mesh with DMForestSetAdaptivityLabel().
1645 
1646   Level: advanced
1647 
1648 .seealso: DMType, DMCreate(), DMSetType(), DMForestGetBaseDM(), DMForestSetBaseDM(), DMForestTemplate(), DMForestSetAdaptivityLabel()
1649 M*/
1650 
1651 PETSC_EXTERN PetscErrorCode DMCreate_Forest(DM dm)
1652 {
1653   DM_Forest      *forest;
1654   PetscErrorCode ierr;
1655 
1656   PetscFunctionBegin;
1657   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1658   ierr                         = PetscNewLog(dm,&forest);CHKERRQ(ierr);
1659   dm->dim                      = 0;
1660   dm->data                     = forest;
1661   forest->refct                = 1;
1662   forest->data                 = NULL;
1663   forest->setfromoptionscalled = PETSC_FALSE;
1664   forest->topology             = NULL;
1665   forest->adapt                = NULL;
1666   forest->base                 = NULL;
1667   forest->adaptPurpose         = DM_ADAPT_DETERMINE;
1668   forest->adjDim               = PETSC_DEFAULT;
1669   forest->overlap              = PETSC_DEFAULT;
1670   forest->minRefinement        = PETSC_DEFAULT;
1671   forest->maxRefinement        = PETSC_DEFAULT;
1672   forest->initRefinement       = PETSC_DEFAULT;
1673   forest->cStart               = PETSC_DETERMINE;
1674   forest->cEnd                 = PETSC_DETERMINE;
1675   forest->cellSF               = NULL;
1676   forest->adaptLabel           = NULL;
1677   forest->gradeFactor          = 2;
1678   forest->cellWeights          = NULL;
1679   forest->cellWeightsCopyMode  = PETSC_USE_POINTER;
1680   forest->weightsFactor        = 1.;
1681   forest->weightCapacity       = 1.;
1682   ierr                         = DMForestSetAdaptivityStrategy(dm,DMFORESTADAPTALL);CHKERRQ(ierr);
1683   ierr                         = DMInitialize_Forest(dm);CHKERRQ(ierr);
1684   PetscFunctionReturn(0);
1685 }
1686 
1687