xref: /petsc/include/petscdmforest.h (revision bcc5ffd92252408533a007b566200b59df89ee0a)
1 /*
2   DMForest, for parallel, hierarchically refined, distributed mesh problems.
3 */
4 #if !defined(__PETSCDMFOREST_H)
5 #define __PETSCDMFOREST_H
6 
7 #include <petscdm.h>
8 
9 /*J
10     DMForestTopology - String with the name of a PETSc DMForest base mesh topology
11 
12    Level: beginner
13 
14 .seealso: DMForestSetType(), DMForest
15 J*/
16 typedef const char* DMForestTopology;
17 
18 /* Just a name for the shape of the domain */
19 PETSC_EXTERN PetscErrorCode DMForestSetTopology(DM, DMForestTopology);
20 PETSC_EXTERN PetscErrorCode DMForestGetTopology(DM, DMForestTopology *);
21 
22 /* this is the coarsest possible forest: can be any DM which we can
23  * convert to a DMForest (right now: plex) */
24 PETSC_EXTERN PetscErrorCode DMForestSetBaseDM(DM, DM);
25 PETSC_EXTERN PetscErrorCode DMForestGetBaseDM(DM, DM *);
26 
27 /* these are forests from which we can adapt */
28 PETSC_EXTERN PetscErrorCode DMForestSetCoarseForest(DM, DM);
29 PETSC_EXTERN PetscErrorCode DMForestGetCoarseForest(DM, DM *);
30 PETSC_EXTERN PetscErrorCode DMForestSetFineForest(DM, DM);
31 PETSC_EXTERN PetscErrorCode DMForestGetFineForest(DM, DM *);
32 
33 /* what we consider adjacent, for the purposes of cell grading, overlap, etc. */
34 PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyDimension(DM, PetscInt);
35 PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyDimension(DM, PetscInt *);
36 PETSC_EXTERN PetscErrorCode DMForestSetAdjacencyCodimension(DM, PetscInt);
37 PETSC_EXTERN PetscErrorCode DMForestGetAdjacencyCodimension(DM, PetscInt *);
38 
39 PETSC_EXTERN PetscErrorCode DMForestSetPartitionOverlap(DM, PetscInt);
40 PETSC_EXTERN PetscErrorCode DMForestGetPartitionOverlap(DM, PetscInt *);
41 
42 PETSC_EXTERN PetscErrorCode DMForestSetMinimumRefinement(DM, PetscInt);
43 PETSC_EXTERN PetscErrorCode DMForestGetMinimumRefinement(DM, PetscInt *);
44 
45 PETSC_EXTERN PetscErrorCode DMForestSetMaximumRefinement(DM, PetscInt);
46 PETSC_EXTERN PetscErrorCode DMForestGetMaximumRefinement(DM, PetscInt *);
47 
48 PETSC_EXTERN PetscErrorCode DMForestSetInitialRefinement(DM, PetscInt);
49 PETSC_EXTERN PetscErrorCode DMForestGetInitialRefinement(DM, PetscInt *);
50 
51 PETSC_EXTERN PetscErrorCode DMForestGetCellChart(DM, PetscInt *, PetscInt *);
52 PETSC_EXTERN PetscErrorCode DMForestGetCellSF(DM, PetscSF *);
53 
54 /* reserve some adaptivity types */
55 enum {DM_FOREST_KEEP = 0,
56       DM_FOREST_REFINE,
57       DM_FOREST_COARSEN,
58       DM_FOREST_RESERVED_ADAPTIVITY_COUNT};
59 
60 /* flag each cell with an adaptivity count: should match the cell section */
61 PETSC_EXTERN PetscErrorCode DMForestSetCellAdaptivityMarkers(DM, PetscInt[], PetscCopyMode);
62 PETSC_EXTERN PetscErrorCode DMForestGetCellAdaptivityMarkers(DM, PetscInt *[]);
63 
64 /*J
65     DMForestAdaptivityStrategy - String with the name of a PETSc DMForest adaptivity strategy
66 
67    Level: intermediate
68 
69 .seealso: DMForestSetType(), DMForest
70 J*/
71 typedef const char* DMForestAdaptivityStrategy;
72 #define DMFORESTADAPTALL "all"
73 #define DMFORESTADAPTANY "any"
74 
75 /* how to combine: -flags         from multiple processes,
76  *                 -coarsen flags from multiple children
77  */
78 PETSC_EXTERN PetscErrorCode DMForestSetAdaptivityStrategy(DM, DMForestAdaptivityStrategy);
79 PETSC_EXTERN PetscErrorCode DMForestGetAdaptivityStrategy(DM, DMForestAdaptivityStrategy *);
80 
81 /* for a quadtree/octree mesh, this is the x:1 condition: 1 indicates a uniform mesh,
82  *                                                        2 indicates typical 2:1,
83  */
84 PETSC_EXTERN PetscErrorCode DMForestSetGradeFactor(DM, PetscInt);
85 PETSC_EXTERN PetscErrorCode DMForestGetGradeFactor(DM, PetscInt *);
86 
87 /* weights for repartitioning */
88 PETSC_EXTERN PetscErrorCode DMForestSetCellWeights(DM, PetscReal[], PetscCopyMode);
89 PETSC_EXTERN PetscErrorCode DMForestGetCellWeights(DM, PetscReal *[]);
90 
91 /* weight multiplier for refinement level: useful for sub-cycling time steps */
92 PETSC_EXTERN PetscErrorCode DMForestSetCellWeightFactor(DM, PetscReal);
93 PETSC_EXTERN PetscErrorCode DMForestGetCellWeightFactor(DM, PetscReal *);
94 
95 /* this process's capacity when redistributing the cells */
96 PETSC_EXTERN PetscErrorCode DMForestSetWeightCapacity(DM, PetscReal);
97 PETSC_EXTERN PetscErrorCode DMForestGetWeightCapacity(DM, PetscReal *);
98 
99 PETSC_EXTERN PetscErrorCode DMForestSetFromOptions(DM);
100 PETSC_EXTERN PetscErrorCode DMForestSetUp(DM);
101 
102 PETSC_EXTERN PetscErrorCode DMForestGetFineProjector(DM,Mat *);
103 PETSC_EXTERN PetscErrorCode DMForestGetCoarseRestrictor(DM,Mat *);
104 
105 /* miscellaneous */
106 PETSC_EXTERN PetscErrorCode DMForestTemplate(DM,DM);
107 
108 /* p4est */
109 PETSC_EXTERN PetscErrorCode DMP4estGetPartitionForCoarsening(DM,PetscBool *);
110 PETSC_EXTERN PetscErrorCode DMP4estSetPartitionForCoarsening(DM,PetscBool);
111 PETSC_EXTERN PetscErrorCode DMP8estGetPartitionForCoarsening(DM,PetscBool *);
112 PETSC_EXTERN PetscErrorCode DMP8estSetPartitionForCoarsening(DM,PetscBool);
113 
114 #endif
115