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:    Note:
 11:    `DM` is an orphan initialism or orphan acronym, the letters have no meaning and never did.

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

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

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

 28:   Level: beginner

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

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

 36:   See <https://scicomp.stackexchange.com/questions/5355/writing-the-poisson-equation-finite-difference-matrix-with-neumann-boundary-cond>

 38:   Developer Note:
 39:   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
 40:   as the 0th grid point where the physical boundary serves as the mirror?

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

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

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

 63:   Level: beginner

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

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

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

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

 88:   Level: intermediate

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

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

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

105:   Level: intermediate

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

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

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

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

131:   Level: beginner

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

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

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

151:   Level: beginner

153: .seealso: [](ch_dmbase), `DM`, `DMAdaptor`, `DMAdaptationStrategy`, `DMAdaptorSolve()`
154: E*/
155: typedef enum {
156:   DM_ADAPTATION_NONE,
157:   DM_ADAPTATION_REFINE,
158:   DM_ADAPTATION_LABEL,
159:   DM_ADAPTATION_METRIC
160: } DMAdaptationCriterion;
161: PETSC_EXTERN const char *const DMAdaptationCriteria[];

163: /*E
164:   DMAdaptFlag - Marker in the label prescribing what adaptation to perform

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

173:   Level: beginner

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

186: /*E
187:   DMDirection - Indicates a coordinate direction

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

194:   Level: beginner

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

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

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

214:   Level: beginner

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

226: /*E
227:   DMPolytopeType - This describes the polytope represented by each cell.

229:   Level: beginner

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

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

261: /*E
262:   PetscUnit - The seven fundamental SI units

264:   Level: beginner

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

279: /*E
280:    DMReorderDefaultFlag - Flag indicating whether the `DM` should be reordered by default

282:    Values:
283: +  `DM_REORDER_DEFAULT_NOTSET` - Flag not set.
284: .  `DM_REORDER_DEFAULT_FALSE`  - Do not reorder by default.
285: -  `DM_REORDER_DEFAULT_TRUE`   - Reorder by default.

287:    Level: intermediate

289:    Developer Note:
290:    Could be replaced with `PETSC_BOOL3`

292: .seealso: `DMPlexReorderSetDefault()`, `DMPlexReorderGetDefault()`, `DMPlexGetOrdering()`, `DMPlexPermute()`
293: E*/
294: typedef enum {
295:   DM_REORDER_DEFAULT_NOTSET = -1,
296:   DM_REORDER_DEFAULT_FALSE  = 0,
297:   DM_REORDER_DEFAULT_TRUE
298: } DMReorderDefaultFlag;

300: /*S
301:     DMField - PETSc object for defining a field on a mesh topology

303:     Level: intermediate

305: .seealso: [](ch_dmbase), `DM`, `DMUniversalLabel`, `DMLabelCreate()`
306: S*/
307: typedef struct _p_DMField *DMField;

309: /*S
310:     DMUniversalLabel - A label that encodes a set of `DMLabel`s, bijectively

312:     Level: developer

314: .seealso: [](ch_dmbase), `DM`, `DMLabel`, `DMUniversalLabelCreate()`
315: S*/
316: typedef struct _p_UniversalLabel *DMUniversalLabel;

318: typedef struct _PETSc_DMCEED *DMCeed;

320: typedef struct _n_DMGeneratorFunctionList *DMGeneratorFunctionList;