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