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