Actual source code: petscdmtypes.h

  1: #ifndef PETSCDMTYPES_H
  2: #define PETSCDMTYPES_H

  4: /* SUBMANSEC = DM */

  6: /*S
  7:      DM - Abstract PETSc object that manages an abstract grid-like object and its interactions with the algebraic solvers

  9:    Level: intermediate

 11: .seealso: [](chapter_dmbase), `DMType`, `DMGetType()`, `DMCompositeCreate()`, `DMDACreate()`, `DMSetType()`, `DMType`, `DMDA`, `DMPLEX`
 12: S*/
 13: typedef struct _p_DM *DM;

 15: /*E
 16:   DMBoundaryType - Describes the choice for the filling of ghost cells on physical domain boundaries.

 18:   Values:
 19: + `DM_BOUNDARY_NONE` - no ghost nodes
 20: . `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
 21: . `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
 22:                          the ghost point value; not yet implemented for 3d
 23: . `DM_BOUNDARY_PERIODIC` - ghost vertices/cells filled by the opposite edge of the domain
 24: - `DM_BOUNDARY_TWIST` - like periodic, only glued backwards like a Mobius strip

 26:   Level: beginner

 28:   Notes:
 29:   This is information for the boundary of the __PHYSICAL__ domain. It has nothing to do with boundaries between
 30:   processes. That width is always determined by the stencil width; see `DMDASetStencilWidth()`.

 32:   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 .

 34:   Developer Notes:
 35:     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
 36:   as the 0th grid point where the physical boundary serves as the mirror?

 38:   References:
 39: . * -  https://scicomp.stackexchange.com/questions/5355/writing-the-poisson-equation-finite-difference-matrix-with-neumann-boundary-cond

 41: .seealso: `DM`, `DMDA`, `DMDASetBoundaryType()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDACreate()`
 42: E*/
 43: typedef enum {
 44:   DM_BOUNDARY_NONE,
 45:   DM_BOUNDARY_GHOSTED,
 46:   DM_BOUNDARY_MIRROR,
 47:   DM_BOUNDARY_PERIODIC,
 48:   DM_BOUNDARY_TWIST
 49: } DMBoundaryType;

 51: /*E
 52:   DMBoundaryConditionType - indicates what type of boundary condition is to be imposed

 54:   Values:
 55: + `DM_BC_ESSENTIAL`       - A Dirichlet condition using a function of the coordinates
 56: . `DM_BC_ESSENTIAL_FIELD` - A Dirichlet condition using a function of the coordinates and auxiliary field data
 57: . `DM_BC_ESSENTIAL_BD_FIELD` - A Dirichlet condition using a function of the coordinates, facet normal, and auxiliary field data
 58: . `DM_BC_NATURAL`         - A Neumann condition using a function of the coordinates
 59: . `DM_BC_NATURAL_FIELD`   - A Neumann condition using a function of the coordinates and auxiliary field data
 60: - `DM_BC_NATURAL_RIEMANN` - A flux condition which determines the state in ghost cells

 62:   Level: beginner

 64:   Note:
 65:   The user can check whether a boundary condition is essential using (type & `DM_BC_ESSENTIAL`), and similarly for
 66:   natural conditions (type & `DM_BC_NATURAL`)

 68: .seealso: `DM`, `DMAddBoundary()`, `DSAddBoundary()`, `DSGetBoundary()`
 69: E*/
 70: typedef enum {
 71:   DM_BC_ESSENTIAL          = 1,
 72:   DM_BC_ESSENTIAL_FIELD    = 5,
 73:   DM_BC_NATURAL            = 2,
 74:   DM_BC_NATURAL_FIELD      = 6,
 75:   DM_BC_ESSENTIAL_BD_FIELD = 9,
 76:   DM_BC_NATURAL_RIEMANN    = 10
 77: } DMBoundaryConditionType;

 79: /*E
 80:   DMPointLocationType - Describes the method to handle point location failure

 82:   Values:
 83: +  `DM_POINTLOCATION_NONE` - return a negative cell number
 84: .  `DM_POINTLOCATION_NEAREST` - the (approximate) nearest point in the mesh is used
 85: -  `DM_POINTLOCATION_REMOVE` - returns values only for points which were located

 87:   Level: intermediate

 89: .seealso: `DM`, `DMLocatePoints()`
 90: E*/
 91: typedef enum {
 92:   DM_POINTLOCATION_NONE,
 93:   DM_POINTLOCATION_NEAREST,
 94:   DM_POINTLOCATION_REMOVE
 95: } DMPointLocationType;

 97: /*E
 98:   DMBlockingType - Describes how to choose variable block sizes

100:   Values:
101: +  `DM_BLOCKING_TOPOLOGICAL_POINT` - select all fields at a topological point (cell center, at a face, etc)
102: -  `DM_BLOCKING_FIELD_NODE` - using a separate block for each field at a topological point

104:   Level: intermediate

106:   Note:
107:   When using `PCVPBJACOBI`, one can choose to block by topological point (all fields at a cell center, at a face, etc.)
108:   or by field nodes (using number of components per field to identify "nodes"). Field nodes lead to smaller blocks, but
109:   may converge more slowly. For example, a cubic Lagrange hexahedron will have one node at vertices, two at edges, four
110:   at faces, and eight at cell centers. If using point blocking, the `PCVPBJACOBI` preconditioner will work with block
111:   sizes up to 8 Lagrange nodes. For 5-component CFD, this produces matrices up to 40x40, which increases memory
112:   footprint and may harm performance. With field node blocking, the maximum block size will correspond to one Lagrange node,
113:   or 5x5 blocks for the CFD example.

115: .seealso: `PCVPBJACOBI`, `MatSetVariableBlockSizes()`, `DMSetBlockingType()`
116: E*/
117: typedef enum {
118:   DM_BLOCKING_TOPOLOGICAL_POINT,
119:   DM_BLOCKING_FIELD_NODE
120: } DMBlockingType;

122: /*E
123:   DMAdaptationStrategy - Describes the strategy used for adaptive solves

125:   Values:
126: +  `DM_ADAPTATION_INITIAL` - refine a mesh based on an initial guess
127: .  `DM_ADAPTATION_SEQUENTIAL` - refine the mesh based on a sequence of solves, much like grid sequencing
128: -  `DM_ADAPTATION_MULTILEVEL` - use the sequence of constructed meshes in a multilevel solve, much like the Systematic Upscaling of Brandt

130:   Level: beginner

132: .seealso: `DM`, `DMAdaptor`, `DMAdaptationCriterion`, `DMAdaptorSolve()`
133: E*/
134: typedef enum {
135:   DM_ADAPTATION_INITIAL,
136:   DM_ADAPTATION_SEQUENTIAL,
137:   DM_ADAPTATION_MULTILEVEL
138: } DMAdaptationStrategy;

140: /*E
141:   DMAdaptationCriterion - Describes the test used to decide whether to coarsen or refine parts of the mesh

143:   Values:
144: + `DM_ADAPTATION_REFINE` - uniformly refine a mesh, much like grid sequencing
145: . `DM_ADAPTATION_LABEL` - adapt the mesh based upon a label of the cells filled with `DMAdaptFlag` markers.
146: . `DM_ADAPTATION_METRIC` - try to mesh the manifold described by the input metric tensor uniformly. PETSc can also construct such a metric based
147:                            upon an input primal or a gradient field.
148: - `DM_ADAPTATION_NONE` - do no adaptation

150:   Level: beginner

152: .seealso: `DM`, `DMAdaptor`, `DMAdaptationStrategy`, `DMAdaptorSolve()`
153: E*/
154: typedef enum {
155:   DM_ADAPTATION_NONE,
156:   DM_ADAPTATION_REFINE,
157:   DM_ADAPTATION_LABEL,
158:   DM_ADAPTATION_METRIC
159: } DMAdaptationCriterion;

161: /*E
162:   DMAdaptFlag - Marker in the label prescribing what adaptation to perform

164:   Values:
165: +  `DM_ADAPT_DETERMINE` - undocumented
166: .  `DM_ADAPT_KEEP` - undocumented
167: .  `DM_ADAPT_REFINE` - undocumented
168: .  `DM_ADAPT_COARSEN` - undocumented
169: -  `DM_ADAPT_COARSEN_LAST` - undocumented

171:   Level: beginner

173: .seealso: `DM`, `DMAdaptor`, `DMAdaptationStrategy`, `DMAdaptationCriterion`, `DMAdaptorSolve()`, `DMAdaptLabel()`
174: E*/
175: typedef enum {
176:   DM_ADAPT_DETERMINE = PETSC_DETERMINE,
177:   DM_ADAPT_KEEP      = 0,
178:   DM_ADAPT_REFINE,
179:   DM_ADAPT_COARSEN,
180:   DM_ADAPT_COARSEN_LAST,
181:   DM_ADAPT_RESERVED_COUNT
182: } DMAdaptFlag;

184: /*E
185:   DMDirection - Indicates a coordinate direction

187:    Values:
188: +  `DM_X` - the x coordinate direction
189: .  `DM_Y` - the y coordinate direction
190: -  `DM_Z` - the z coordinate direction

192:   Level: beginner

194: .seealso: `DM`, `DMDA`, `DMDAGetRay()`, `DMDAGetProcessorSubset()`, `DMPlexShearGeometry()`
195: E*/
196: typedef enum {
197:   DM_X,
198:   DM_Y,
199:   DM_Z
200: } DMDirection;

202: /*E
203:   DMEnclosureType - The type of enclosure relation between one `DM` and another

205:    Values:
206: +  `DM_ENC_SUBMESH` - the `DM` is the boundary of another `DM`
207: .  `DM_ENC_SUPERMESH`  - the `DM` has the boundary of another `DM` (the reverse situation to `DM_ENC_SUBMESH`)
208: .  `DM_ENC_EQUALITY` - unknown what this means
209: .  `DM_ENC_NONE` - no relationship can be determined
210: -  `DM_ENC_UNKNOWN`  - the relationship is unknown

212:   Level: beginner

214: .seealso: `DM`, `DMGetEnclosureRelation()`
215: E*/
216: typedef enum {
217:   DM_ENC_EQUALITY,
218:   DM_ENC_SUPERMESH,
219:   DM_ENC_SUBMESH,
220:   DM_ENC_NONE,
221:   DM_ENC_UNKNOWN
222: } DMEnclosureType;

224: /*E
225:   DMPolytopeType - This describes the polytope represented by each cell.

227:   Level: beginner

229:   While most operations only need the topology information in the `DMPLEX`, we must sometimes have the
230:   user specify a polytope. For instance, when interpolating from a cell-vertex mesh, the type of
231:   polytope can be ambiguous. Also, `DMPLEX` allows different symmetries of a prism cell with the same
232:   constituent points. Normally these types are autoamtically inferred and the user does not specify
233:   them.

235: .seealso: `DM`, `DMPlexComputeCellTypes()`
236: E*/
237: typedef enum {
238:   DM_POLYTOPE_POINT,
239:   DM_POLYTOPE_SEGMENT,
240:   DM_POLYTOPE_POINT_PRISM_TENSOR,
241:   DM_POLYTOPE_TRIANGLE,
242:   DM_POLYTOPE_QUADRILATERAL,
243:   DM_POLYTOPE_SEG_PRISM_TENSOR,
244:   DM_POLYTOPE_TETRAHEDRON,
245:   DM_POLYTOPE_HEXAHEDRON,
246:   DM_POLYTOPE_TRI_PRISM,
247:   DM_POLYTOPE_TRI_PRISM_TENSOR,
248:   DM_POLYTOPE_QUAD_PRISM_TENSOR,
249:   DM_POLYTOPE_PYRAMID,
250:   DM_POLYTOPE_FV_GHOST,
251:   DM_POLYTOPE_INTERIOR_GHOST,
252:   DM_POLYTOPE_UNKNOWN,
253:   DM_NUM_POLYTOPES
254: } DMPolytopeType;
255: PETSC_EXTERN const char *const DMPolytopeTypes[];

257: /*E
258:   PetscUnit - The seven fundamental SI units

260:   Level: beginner

262: .seealso: `DMPlexGetScale()`, `DMPlexSetScale()`
263: E*/
264: typedef enum {
265:   PETSC_UNIT_LENGTH,
266:   PETSC_UNIT_MASS,
267:   PETSC_UNIT_TIME,
268:   PETSC_UNIT_CURRENT,
269:   PETSC_UNIT_TEMPERATURE,
270:   PETSC_UNIT_AMOUNT,
271:   PETSC_UNIT_LUMINOSITY,
272:   NUM_PETSC_UNITS
273: } PetscUnit;

275: /*S
276:     DMField - PETSc object for defining a field on a mesh topology

278:     Level: intermediate
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
286: S*/
287: typedef struct _p_UniversalLabel *DMUniversalLabel;

289: typedef struct _n_DMGeneratorFunctionList *DMGeneratorFunctionList;

291: #endif