xref: /petsc/include/petscdmtypes.h (revision df4cd43f92eaa320656440c40edb1046daee8f75)
1 #ifndef PETSCDMTYPES_H
2 #define PETSCDMTYPES_H
3 
4 /* SUBMANSEC = DM */
5 
6 /*S
7      DM - Abstract PETSc object that manages an abstract grid object and its interactions with the algebraic solvers
8 
9    Level: intermediate
10 
11 .seealso: `DMType`, `DMGetType()`, `DMCompositeCreate()`, `DMDACreate()`, `DMSetType()`, `DMType`, `DMDA`, `DMPLEX`
12 S*/
13 typedef struct _p_DM *DM;
14 
15 /*E
16   DMBoundaryType - Describes the choice for fill of ghost cells on physical domain boundaries.
17 
18   Level: beginner
19 
20   A boundary may be of type `DM_BOUNDARY_NONE` (no ghost nodes), `DM_BOUNDARY_GHOSTED` (ghost vertices/cells
21   exist but aren't filled; you can put values into them and then apply a stencil that uses those ghost locations),
22   `DM_BOUNDARY_MIRROR` (the ghost value is the same as the value 1 grid point in; that is, the 0th grid point in the real mesh acts like a mirror to define the ghost point value;
23   not yet implemented for 3d), `DM_BOUNDARY_PERIODIC` (ghost vertices/cells filled by the opposite
24   edge of the domain), or `DM_BOUNDARY_TWIST` (like periodic, only glued backwards like a Mobius strip).
25 
26   Notes:
27   This is information for the boundary of the __PHYSICAL__ domain. It has nothing to do with boundaries between
28   processes. That width is always determined by the stencil width; see `DMDASetStencilWidth()`.
29 
30   If the physical grid points have values 0 1 2 3 with `DM_BOUNDARY_MIRROR` then the local vector with ghost points has the values 1 0 1 2 3 2 .
31 
32   Developer Notes:
33     Should` DM_BOUNDARY_MIRROR` have the same meaning with DMDA_Q0, that is a staggered grid? In that case should the ghost point have the same value
34   as the 0th grid point where the physical boundary serves as the mirror?
35 
36   References:
37 . * -  https://scicomp.stackexchange.com/questions/5355/writing-the-poisson-equation-finite-difference-matrix-with-neumann-boundary-cond
38 
39 .seealso: `DMDASetBoundaryType()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDACreate()`
40 E*/
41 typedef enum {
42   DM_BOUNDARY_NONE,
43   DM_BOUNDARY_GHOSTED,
44   DM_BOUNDARY_MIRROR,
45   DM_BOUNDARY_PERIODIC,
46   DM_BOUNDARY_TWIST
47 } DMBoundaryType;
48 /*E
49   DMBoundaryConditionType - indicates what type of boundary condition is to be imposed
50 
51   Note: This flag indicates the type of function which will define the condition:
52 $ DM_BC_ESSENTIAL       - A Dirichlet condition using a function of the coordinates
53 $ DM_BC_ESSENTIAL_FIELD - A Dirichlet condition using a function of the coordinates and auxiliary field data
54 $ DM_BC_ESSENTIAL_BD_FIELD - A Dirichlet condition using a function of the coordinates, facet normal, and auxiliary field data
55 $ DM_BC_NATURAL         - A Neumann condition using a function of the coordinates
56 $ DM_BC_NATURAL_FIELD   - A Neumann condition using a function of the coordinates and auxiliary field data
57 $ DM_BC_NATURAL_RIEMANN - A flux condition which determines the state in ghost cells
58 The user can check whether a boundary condition is essential using (type & DM_BC_ESSENTIAL), and similarly for
59 natural conditions (type & DM_BC_NATURAL)
60 
61   Level: beginner
62 
63 .seealso: `DMAddBoundary()`, `DSAddBoundary()`, `DSGetBoundary()`
64 E*/
65 typedef enum {
66   DM_BC_ESSENTIAL          = 1,
67   DM_BC_ESSENTIAL_FIELD    = 5,
68   DM_BC_NATURAL            = 2,
69   DM_BC_NATURAL_FIELD      = 6,
70   DM_BC_ESSENTIAL_BD_FIELD = 9,
71   DM_BC_NATURAL_RIEMANN    = 10
72 } DMBoundaryConditionType;
73 
74 /*E
75   DMPointLocationType - Describes the method to handle point location failure
76 
77   Level: beginner
78 
79   If a search using `DM_POINTLOCATION_NONE` fails, the failure is signaled with a negative cell number. On the
80   other hand, if `DM_POINTLOCATION_NEAREST` is used, on failure, the (approximate) nearest point in the mesh is
81   used, replacing the given point in the input vector. `DM_POINTLOCATION_REMOVE` returns values only for points
82   which were located.
83 
84 .seealso: `DMLocatePoints()`
85 E*/
86 typedef enum {
87   DM_POINTLOCATION_NONE,
88   DM_POINTLOCATION_NEAREST,
89   DM_POINTLOCATION_REMOVE
90 } DMPointLocationType;
91 
92 /*E
93   DMBlockingType - Describes how to choose variable block sizes
94 
95   Level: intermediate
96 
97   When using `PCVPBJACOBI`, one can choose to block by topological point (all fields at a cell center, at a face, etc.)
98   or by field nodes (using number of components per field to identify "nodes"). Field nodes lead to smaller blocks, but
99   may converge more slowly. For example, a cubic Lagrange hexahedron will have one node at vertices, two at edges, four
100   at faces, and eight at cell centers. If using point blocking, the `PCVPBJACOBI` preconditioner will work with block
101   sizes up to 8 Lagrange nodes. For 5-component CFD, this produces matrices up to 40x40, which increases memory
102   footprint and may harm performance. With field node blocking, the max block size will correspond to one Lagrange node,
103   or 5x5 blocks for the CFD example.
104 
105 .seealso: `PCVPBJACOBI`, `MatSetVariableBlockSizes()`, `DMSetBlockingType()`
106 E*/
107 typedef enum {
108   DM_BLOCKING_TOPOLOGICAL_POINT,
109   DM_BLOCKING_FIELD_NODE
110 } DMBlockingType;
111 
112 /*E
113   DMAdaptationStrategy - Describes the strategy used for adaptive solves
114 
115   Level: beginner
116 
117   DM_ADAPTATION_INITIAL will refine a mesh based on an initial guess. DM_ADAPTATION_SEQUENTIAL will refine the
118   mesh based on a sequence of solves, much like grid sequencing. DM_ADAPTATION_MULTILEVEL will use the sequence
119   of constructed meshes in a multilevel solve, much like the Systematic Upscaling of Brandt.
120 
121 .seealso: `DMAdaptorSolve()`
122 E*/
123 typedef enum {
124   DM_ADAPTATION_INITIAL,
125   DM_ADAPTATION_SEQUENTIAL,
126   DM_ADAPTATION_MULTILEVEL
127 } DMAdaptationStrategy;
128 
129 /*E
130   DMAdaptationCriterion - Describes the test used to decide whether to coarsen or refine parts of the mesh
131 
132   Level: beginner
133 
134   `DM_ADAPTATION_REFINE` will uniformly refine a mesh, much like grid sequencing. `DM_ADAPTATION_LABEL` will adapt
135   the mesh based upon a label of the cells filled with `DMAdaptFlag` markers. `DM_ADAPTATION_METRIC` will try to
136   mesh the manifold described by the input metric tensor uniformly. PETSc can also construct such a metric based
137   upon an input primal or a gradient field.
138 
139 .seealso: `DMAdaptorSolve()`
140 E*/
141 typedef enum {
142   DM_ADAPTATION_NONE,
143   DM_ADAPTATION_REFINE,
144   DM_ADAPTATION_LABEL,
145   DM_ADAPTATION_METRIC
146 } DMAdaptationCriterion;
147 
148 /*E
149   DMAdaptFlag - Marker in the label prescribing adaptation
150 
151   Level: beginner
152 
153 .seealso: `DMAdaptLabel()`
154 E*/
155 typedef enum {
156   DM_ADAPT_DETERMINE = PETSC_DETERMINE,
157   DM_ADAPT_KEEP      = 0,
158   DM_ADAPT_REFINE,
159   DM_ADAPT_COARSEN,
160   DM_ADAPT_COARSEN_LAST,
161   DM_ADAPT_RESERVED_COUNT
162 } DMAdaptFlag;
163 
164 /*E
165   DMDirection - Indicates a coordinate direction
166 
167   Level: beginner
168 
169 .seealso: `DMDAGetRay()`, `DMDAGetProcessorSubset()`, `DMPlexShearGeometry()`
170 E*/
171 typedef enum {
172   DM_X,
173   DM_Y,
174   DM_Z
175 } DMDirection;
176 
177 /*E
178   DMEnclosureType - The type of enclosure relation between one `DM` and another
179 
180   Level: beginner
181 
182   Notes:
183   For example, one `DM` dmA may be the boundary of another dmB, in which case it would be labeled `DM_ENC_SUBMESH`.
184 
185   If the situation is reversed, and dmA has boundary dmB, it would be labeled `DM_ENC_SUPERMESH`.
186 
187   Likewise, if dmA was a subregion of dmB, it would be labeled `DM_ENC_SUBMESH`.
188 
189   If no relation can be determined, `DM_ENC_NONE` is used.
190 
191   If a relation is not yet known, then `DM_ENC_UNKNOWN` is used.
192 
193 .seealso: `DMGetEnclosureRelation()`
194 E*/
195 typedef enum {
196   DM_ENC_EQUALITY,
197   DM_ENC_SUPERMESH,
198   DM_ENC_SUBMESH,
199   DM_ENC_NONE,
200   DM_ENC_UNKNOWN
201 } DMEnclosureType;
202 
203 /*E
204   DMPolytopeType - This describes the polytope represented by each cell.
205 
206   Level: beginner
207 
208   While most operations only need the topology information in the `DMPLEX`, we must sometimes have the
209   user specify a polytope. For instance, when interpolating from a cell-vertex mesh, the type of
210   polytope can be ambiguous. Also, `DMPLEX` allows different symmetries of a prism cell with the same
211   constituent points. Normally these types are autoamtically inferred and the user does not specify
212   them.
213 
214 .seealso: `DMPlexComputeCellTypes()`
215 E*/
216 typedef enum {
217   DM_POLYTOPE_POINT,
218   DM_POLYTOPE_SEGMENT,
219   DM_POLYTOPE_POINT_PRISM_TENSOR,
220   DM_POLYTOPE_TRIANGLE,
221   DM_POLYTOPE_QUADRILATERAL,
222   DM_POLYTOPE_SEG_PRISM_TENSOR,
223   DM_POLYTOPE_TETRAHEDRON,
224   DM_POLYTOPE_HEXAHEDRON,
225   DM_POLYTOPE_TRI_PRISM,
226   DM_POLYTOPE_TRI_PRISM_TENSOR,
227   DM_POLYTOPE_QUAD_PRISM_TENSOR,
228   DM_POLYTOPE_PYRAMID,
229   DM_POLYTOPE_FV_GHOST,
230   DM_POLYTOPE_INTERIOR_GHOST,
231   DM_POLYTOPE_UNKNOWN,
232   DM_NUM_POLYTOPES
233 } DMPolytopeType;
234 PETSC_EXTERN const char *const DMPolytopeTypes[];
235 
236 /*E
237   PetscUnit - The seven fundamental SI units
238 
239   Level: beginner
240 
241 .seealso: `DMPlexGetScale()`, `DMPlexSetScale()`
242 E*/
243 typedef enum {
244   PETSC_UNIT_LENGTH,
245   PETSC_UNIT_MASS,
246   PETSC_UNIT_TIME,
247   PETSC_UNIT_CURRENT,
248   PETSC_UNIT_TEMPERATURE,
249   PETSC_UNIT_AMOUNT,
250   PETSC_UNIT_LUMINOSITY,
251   NUM_PETSC_UNITS
252 } PetscUnit;
253 
254 /*S
255     DMField - PETSc object for defining a field on a mesh topology
256 
257     Level: intermediate
258 S*/
259 typedef struct _p_DMField *DMField;
260 
261 /*S
262     DMUniversalLabel - A label that encodes a set of `DMLabel`s, bijectively
263 
264     Level: developer
265 S*/
266 typedef struct _p_UniversalLabel *DMUniversalLabel;
267 
268 typedef struct _n_DMGeneratorFunctionList *DMGeneratorFunctionList;
269 
270 #endif
271