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