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