xref: /petsc/include/petscdmtypes.h (revision 16a05f60a523f53ab316acaac9f77b7425611adc)
16524c165SJacob Faibussowitsch #ifndef PETSCDMTYPES_H
226bd1501SBarry Smith #define PETSCDMTYPES_H
31e25c274SJed Brown 
4ac09b921SBarry Smith /* SUBMANSEC = DM */
5ac09b921SBarry Smith 
61e25c274SJed Brown /*S
7*16a05f60SBarry Smith      DM - Abstract PETSc object that manages an abstract grid-like object and its interactions with the algebraic solvers
81e25c274SJed Brown 
91e25c274SJed Brown    Level: intermediate
101e25c274SJed Brown 
11*16a05f60SBarry Smith .seealso: [](chapter_dmbase), `DMType`, `DMGetType()`, `DMCompositeCreate()`, `DMDACreate()`, `DMSetType()`, `DMType`, `DMDA`, `DMPLEX`
121e25c274SJed Brown S*/
131e25c274SJed Brown typedef struct _p_DM *DM;
141e25c274SJed Brown 
15bff4a2f0SMatthew G. Knepley /*E
16*16a05f60SBarry Smith   DMBoundaryType - Describes the choice for the filling of ghost cells on physical domain boundaries.
17*16a05f60SBarry Smith 
18*16a05f60SBarry Smith   Values:
19*16a05f60SBarry Smith + `DM_BOUNDARY_NONE` - no ghost nodes
20*16a05f60SBarry Smith . `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*16a05f60SBarry Smith . `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*16a05f60SBarry Smith                          the ghost point value; not yet implemented for 3d
23*16a05f60SBarry Smith . `DM_BOUNDARY_PERIODIC` - ghost vertices/cells filled by the opposite edge of the domain
24*16a05f60SBarry Smith - `DM_BOUNDARY_TWIST` - like periodic, only glued backwards like a Mobius strip
25bff4a2f0SMatthew G. Knepley 
26bff4a2f0SMatthew G. Knepley   Level: beginner
27bff4a2f0SMatthew G. Knepley 
28dbb368e6SPatrick Sanan   Notes:
29dbb368e6SPatrick Sanan   This is information for the boundary of the __PHYSICAL__ domain. It has nothing to do with boundaries between
3087497f52SBarry Smith   processes. That width is always determined by the stencil width; see `DMDASetStencilWidth()`.
31bff4a2f0SMatthew G. Knepley 
3287497f52SBarry Smith   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 .
33288e7d53SBarry Smith 
3495452b02SPatrick Sanan   Developer Notes:
3587497f52SBarry Smith     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
36288e7d53SBarry Smith   as the 0th grid point where the physical boundary serves as the mirror?
37288e7d53SBarry Smith 
38dbb368e6SPatrick Sanan   References:
39606c0280SSatish Balay . * -  https://scicomp.stackexchange.com/questions/5355/writing-the-poisson-equation-finite-difference-matrix-with-neumann-boundary-cond
40288e7d53SBarry Smith 
41*16a05f60SBarry Smith .seealso: `DM`, `DMDA`, `DMDASetBoundaryType()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDACreate()`
42bff4a2f0SMatthew G. Knepley E*/
439371c9d4SSatish Balay typedef enum {
449371c9d4SSatish Balay   DM_BOUNDARY_NONE,
459371c9d4SSatish Balay   DM_BOUNDARY_GHOSTED,
469371c9d4SSatish Balay   DM_BOUNDARY_MIRROR,
479371c9d4SSatish Balay   DM_BOUNDARY_PERIODIC,
489371c9d4SSatish Balay   DM_BOUNDARY_TWIST
499371c9d4SSatish Balay } DMBoundaryType;
50*16a05f60SBarry Smith 
5162a38674SMatthew G. Knepley /*E
529dc85fa5SMatthew G. Knepley   DMBoundaryConditionType - indicates what type of boundary condition is to be imposed
539dc85fa5SMatthew G. Knepley 
54*16a05f60SBarry Smith   Values:
55*16a05f60SBarry Smith + `DM_BC_ESSENTIAL`       - A Dirichlet condition using a function of the coordinates
56*16a05f60SBarry Smith . `DM_BC_ESSENTIAL_FIELD` - A Dirichlet condition using a function of the coordinates and auxiliary field data
57*16a05f60SBarry Smith . `DM_BC_ESSENTIAL_BD_FIELD` - A Dirichlet condition using a function of the coordinates, facet normal, and auxiliary field data
58*16a05f60SBarry Smith . `DM_BC_NATURAL`         - A Neumann condition using a function of the coordinates
59*16a05f60SBarry Smith . `DM_BC_NATURAL_FIELD`   - A Neumann condition using a function of the coordinates and auxiliary field data
60*16a05f60SBarry Smith - `DM_BC_NATURAL_RIEMANN` - A flux condition which determines the state in ghost cells
619dc85fa5SMatthew G. Knepley 
629dc85fa5SMatthew G. Knepley   Level: beginner
639dc85fa5SMatthew G. Knepley 
64*16a05f60SBarry Smith   Note:
65*16a05f60SBarry Smith   The user can check whether a boundary condition is essential using (type & `DM_BC_ESSENTIAL`), and similarly for
66*16a05f60SBarry Smith   natural conditions (type & `DM_BC_NATURAL`)
67*16a05f60SBarry Smith 
68*16a05f60SBarry Smith .seealso: `DM`, `DMAddBoundary()`, `DSAddBoundary()`, `DSGetBoundary()`
699dc85fa5SMatthew G. Knepley E*/
709371c9d4SSatish Balay typedef enum {
719371c9d4SSatish Balay   DM_BC_ESSENTIAL          = 1,
729371c9d4SSatish Balay   DM_BC_ESSENTIAL_FIELD    = 5,
739371c9d4SSatish Balay   DM_BC_NATURAL            = 2,
749371c9d4SSatish Balay   DM_BC_NATURAL_FIELD      = 6,
759371c9d4SSatish Balay   DM_BC_ESSENTIAL_BD_FIELD = 9,
769371c9d4SSatish Balay   DM_BC_NATURAL_RIEMANN    = 10
779371c9d4SSatish Balay } DMBoundaryConditionType;
789dc85fa5SMatthew G. Knepley 
799dc85fa5SMatthew G. Knepley /*E
8062a38674SMatthew G. Knepley   DMPointLocationType - Describes the method to handle point location failure
8162a38674SMatthew G. Knepley 
82*16a05f60SBarry Smith   Values:
83*16a05f60SBarry Smith +  `DM_POINTLOCATION_NONE` - return a negative cell number
84*16a05f60SBarry Smith .  `DM_POINTLOCATION_NEAREST` - the (approximate) nearest point in the mesh is used
85*16a05f60SBarry Smith -  `DM_POINTLOCATION_REMOVE` - returns values only for points which were located
8662a38674SMatthew G. Knepley 
87*16a05f60SBarry Smith   Level: intermediate
8862a38674SMatthew G. Knepley 
89*16a05f60SBarry Smith .seealso: `DM`, `DMLocatePoints()`
9062a38674SMatthew G. Knepley E*/
919371c9d4SSatish Balay typedef enum {
929371c9d4SSatish Balay   DM_POINTLOCATION_NONE,
939371c9d4SSatish Balay   DM_POINTLOCATION_NEAREST,
949371c9d4SSatish Balay   DM_POINTLOCATION_REMOVE
959371c9d4SSatish Balay } DMPointLocationType;
9662a38674SMatthew G. Knepley 
975675c177SMatthew G. Knepley /*E
98863027abSJed Brown   DMBlockingType - Describes how to choose variable block sizes
99863027abSJed Brown 
100*16a05f60SBarry Smith   Values:
101*16a05f60SBarry Smith +  `DM_BLOCKING_TOPOLOGICAL_POINT` - select all fields at a topological point (cell center, at a face, etc)
102*16a05f60SBarry Smith -  `DM_BLOCKING_FIELD_NODE` - using a separate block for each field at a topological point
103*16a05f60SBarry Smith 
104863027abSJed Brown   Level: intermediate
105863027abSJed Brown 
106*16a05f60SBarry Smith   Note:
107863027abSJed Brown   When using `PCVPBJACOBI`, one can choose to block by topological point (all fields at a cell center, at a face, etc.)
108863027abSJed Brown   or by field nodes (using number of components per field to identify "nodes"). Field nodes lead to smaller blocks, but
109863027abSJed Brown   may converge more slowly. For example, a cubic Lagrange hexahedron will have one node at vertices, two at edges, four
110863027abSJed Brown   at faces, and eight at cell centers. If using point blocking, the `PCVPBJACOBI` preconditioner will work with block
111863027abSJed Brown   sizes up to 8 Lagrange nodes. For 5-component CFD, this produces matrices up to 40x40, which increases memory
112*16a05f60SBarry Smith   footprint and may harm performance. With field node blocking, the maximum block size will correspond to one Lagrange node,
113863027abSJed Brown   or 5x5 blocks for the CFD example.
114863027abSJed Brown 
115863027abSJed Brown .seealso: `PCVPBJACOBI`, `MatSetVariableBlockSizes()`, `DMSetBlockingType()`
116863027abSJed Brown E*/
117863027abSJed Brown typedef enum {
1180e762ea3SJed Brown   DM_BLOCKING_TOPOLOGICAL_POINT,
119863027abSJed Brown   DM_BLOCKING_FIELD_NODE
120863027abSJed Brown } DMBlockingType;
121863027abSJed Brown 
122863027abSJed Brown /*E
123174e7490SMatthew G. Knepley   DMAdaptationStrategy - Describes the strategy used for adaptive solves
1245675c177SMatthew G. Knepley 
125*16a05f60SBarry Smith   Values:
126*16a05f60SBarry Smith +  `DM_ADAPTATION_INITIAL` - refine a mesh based on an initial guess
127*16a05f60SBarry Smith .  `DM_ADAPTATION_SEQUENTIAL` - refine the mesh based on a sequence of solves, much like grid sequencing
128*16a05f60SBarry Smith -  `DM_ADAPTATION_MULTILEVEL` - use the sequence of constructed meshes in a multilevel solve, much like the Systematic Upscaling of Brandt
129*16a05f60SBarry Smith 
1305675c177SMatthew G. Knepley   Level: beginner
1315675c177SMatthew G. Knepley 
132*16a05f60SBarry Smith .seealso: `DM`, `DMAdaptor`, `DMAdaptationCriterion`, `DMAdaptorSolve()`
1335675c177SMatthew G. Knepley E*/
1349371c9d4SSatish Balay typedef enum {
1359371c9d4SSatish Balay   DM_ADAPTATION_INITIAL,
1369371c9d4SSatish Balay   DM_ADAPTATION_SEQUENTIAL,
1379371c9d4SSatish Balay   DM_ADAPTATION_MULTILEVEL
1389371c9d4SSatish Balay } DMAdaptationStrategy;
139174e7490SMatthew G. Knepley 
140174e7490SMatthew G. Knepley /*E
141174e7490SMatthew G. Knepley   DMAdaptationCriterion - Describes the test used to decide whether to coarsen or refine parts of the mesh
142174e7490SMatthew G. Knepley 
143*16a05f60SBarry Smith   Values:
144*16a05f60SBarry Smith + `DM_ADAPTATION_REFINE` - uniformly refine a mesh, much like grid sequencing
145*16a05f60SBarry Smith . `DM_ADAPTATION_LABEL` - adapt the mesh based upon a label of the cells filled with `DMAdaptFlag` markers.
146*16a05f60SBarry Smith . `DM_ADAPTATION_METRIC` - try to mesh the manifold described by the input metric tensor uniformly. PETSc can also construct such a metric based
147*16a05f60SBarry Smith                            upon an input primal or a gradient field.
148*16a05f60SBarry Smith - `DM_ADAPTATION_NONE` - do no adaptation
149*16a05f60SBarry Smith 
150174e7490SMatthew G. Knepley   Level: beginner
151174e7490SMatthew G. Knepley 
152*16a05f60SBarry Smith .seealso: `DM`, `DMAdaptor`, `DMAdaptationStrategy`, `DMAdaptorSolve()`
153174e7490SMatthew G. Knepley E*/
1549371c9d4SSatish Balay typedef enum {
1559371c9d4SSatish Balay   DM_ADAPTATION_NONE,
1569371c9d4SSatish Balay   DM_ADAPTATION_REFINE,
1579371c9d4SSatish Balay   DM_ADAPTATION_LABEL,
1589371c9d4SSatish Balay   DM_ADAPTATION_METRIC
1599371c9d4SSatish Balay } DMAdaptationCriterion;
1605675c177SMatthew G. Knepley 
1619dc85fa5SMatthew G. Knepley /*E
162*16a05f60SBarry Smith   DMAdaptFlag - Marker in the label prescribing what adaptation to perform
163*16a05f60SBarry Smith 
164*16a05f60SBarry Smith   Values:
165*16a05f60SBarry Smith +  `DM_ADAPT_DETERMINE` - undocumented
166*16a05f60SBarry Smith .  `DM_ADAPT_KEEP` - undocumented
167*16a05f60SBarry Smith .  `DM_ADAPT_REFINE` - undocumented
168*16a05f60SBarry Smith .  `DM_ADAPT_COARSEN` - undocumented
169*16a05f60SBarry Smith -  `DM_ADAPT_COARSEN_LAST` - undocumented
1709dc85fa5SMatthew G. Knepley 
1719dc85fa5SMatthew G. Knepley   Level: beginner
1729dc85fa5SMatthew G. Knepley 
173*16a05f60SBarry Smith .seealso: `DM`, `DMAdaptor`, `DMAdaptationStrategy`, `DMAdaptationCriterion`, `DMAdaptorSolve()`, `DMAdaptLabel()`
1749dc85fa5SMatthew G. Knepley E*/
1759371c9d4SSatish Balay typedef enum {
1769371c9d4SSatish Balay   DM_ADAPT_DETERMINE = PETSC_DETERMINE,
1779371c9d4SSatish Balay   DM_ADAPT_KEEP      = 0,
1789371c9d4SSatish Balay   DM_ADAPT_REFINE,
1799371c9d4SSatish Balay   DM_ADAPT_COARSEN,
1809371c9d4SSatish Balay   DM_ADAPT_COARSEN_LAST,
1819371c9d4SSatish Balay   DM_ADAPT_RESERVED_COUNT
1829371c9d4SSatish Balay } DMAdaptFlag;
1833ee9839eSMatthew G. Knepley 
1843ee9839eSMatthew G. Knepley /*E
1853ee9839eSMatthew G. Knepley   DMDirection - Indicates a coordinate direction
1863ee9839eSMatthew G. Knepley 
187*16a05f60SBarry Smith    Values:
188*16a05f60SBarry Smith +  `DM_X` - the x coordinate direction
189*16a05f60SBarry Smith .  `DM_Y` - the y coordinate direction
190*16a05f60SBarry Smith -  `DM_Z` - the z coordinate direction
191*16a05f60SBarry Smith 
1923ee9839eSMatthew G. Knepley   Level: beginner
1933ee9839eSMatthew G. Knepley 
194*16a05f60SBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetRay()`, `DMDAGetProcessorSubset()`, `DMPlexShearGeometry()`
1953ee9839eSMatthew G. Knepley E*/
1969371c9d4SSatish Balay typedef enum {
1979371c9d4SSatish Balay   DM_X,
1989371c9d4SSatish Balay   DM_Y,
1999371c9d4SSatish Balay   DM_Z
2009371c9d4SSatish Balay } DMDirection;
2019dc85fa5SMatthew G. Knepley 
202a6e0b375SMatthew G. Knepley /*E
20387497f52SBarry Smith   DMEnclosureType - The type of enclosure relation between one `DM` and another
204a6e0b375SMatthew G. Knepley 
205*16a05f60SBarry Smith    Values:
206*16a05f60SBarry Smith +  `DM_ENC_SUBMESH` - the `DM` is the boundary of another `DM`
207*16a05f60SBarry Smith .  `DM_ENC_SUPERMESH`  - the `DM` has the boundary of another `DM` (the reverse situation to `DM_ENC_SUBMESH`)
208*16a05f60SBarry Smith .  `DM_ENC_EQUALITY` - unkown what this means
209*16a05f60SBarry Smith .  `DM_ENC_NONE` - no relatiship can be determined
210*16a05f60SBarry Smith -  `DM_ENC_UNKNOWN`  - the relationship is unknown
211*16a05f60SBarry Smith 
212a6e0b375SMatthew G. Knepley   Level: beginner
213a6e0b375SMatthew G. Knepley 
214*16a05f60SBarry Smith .seealso: `DM`, `DMGetEnclosureRelation()`
215a6e0b375SMatthew G. Knepley E*/
2169371c9d4SSatish Balay typedef enum {
2179371c9d4SSatish Balay   DM_ENC_EQUALITY,
2189371c9d4SSatish Balay   DM_ENC_SUPERMESH,
2199371c9d4SSatish Balay   DM_ENC_SUBMESH,
2209371c9d4SSatish Balay   DM_ENC_NONE,
2219371c9d4SSatish Balay   DM_ENC_UNKNOWN
2229371c9d4SSatish Balay } DMEnclosureType;
223ba2698f1SMatthew G. Knepley 
224ba2698f1SMatthew G. Knepley /*E
225ba2698f1SMatthew G. Knepley   DMPolytopeType - This describes the polytope represented by each cell.
226ba2698f1SMatthew G. Knepley 
227ba2698f1SMatthew G. Knepley   Level: beginner
228ba2698f1SMatthew G. Knepley 
22987497f52SBarry Smith   While most operations only need the topology information in the `DMPLEX`, we must sometimes have the
230ba2698f1SMatthew G. Knepley   user specify a polytope. For instance, when interpolating from a cell-vertex mesh, the type of
23187497f52SBarry Smith   polytope can be ambiguous. Also, `DMPLEX` allows different symmetries of a prism cell with the same
232ba2698f1SMatthew G. Knepley   constituent points. Normally these types are autoamtically inferred and the user does not specify
233ba2698f1SMatthew G. Knepley   them.
234ba2698f1SMatthew G. Knepley 
235*16a05f60SBarry Smith .seealso: `DM`, `DMPlexComputeCellTypes()`
236ba2698f1SMatthew G. Knepley E*/
2379371c9d4SSatish Balay typedef enum {
2389371c9d4SSatish Balay   DM_POLYTOPE_POINT,
2399371c9d4SSatish Balay   DM_POLYTOPE_SEGMENT,
2409371c9d4SSatish Balay   DM_POLYTOPE_POINT_PRISM_TENSOR,
2419371c9d4SSatish Balay   DM_POLYTOPE_TRIANGLE,
2429371c9d4SSatish Balay   DM_POLYTOPE_QUADRILATERAL,
2439371c9d4SSatish Balay   DM_POLYTOPE_SEG_PRISM_TENSOR,
2449371c9d4SSatish Balay   DM_POLYTOPE_TETRAHEDRON,
2459371c9d4SSatish Balay   DM_POLYTOPE_HEXAHEDRON,
2469371c9d4SSatish Balay   DM_POLYTOPE_TRI_PRISM,
2479371c9d4SSatish Balay   DM_POLYTOPE_TRI_PRISM_TENSOR,
2489371c9d4SSatish Balay   DM_POLYTOPE_QUAD_PRISM_TENSOR,
2499371c9d4SSatish Balay   DM_POLYTOPE_PYRAMID,
2509371c9d4SSatish Balay   DM_POLYTOPE_FV_GHOST,
2519371c9d4SSatish Balay   DM_POLYTOPE_INTERIOR_GHOST,
2529371c9d4SSatish Balay   DM_POLYTOPE_UNKNOWN,
2539371c9d4SSatish Balay   DM_NUM_POLYTOPES
2549371c9d4SSatish Balay } DMPolytopeType;
255ba2698f1SMatthew G. Knepley PETSC_EXTERN const char *const DMPolytopeTypes[];
256a6e0b375SMatthew G. Knepley 
2579dc85fa5SMatthew G. Knepley /*E
2589dc85fa5SMatthew G. Knepley   PetscUnit - The seven fundamental SI units
2599dc85fa5SMatthew G. Knepley 
2609dc85fa5SMatthew G. Knepley   Level: beginner
2619dc85fa5SMatthew G. Knepley 
262db781477SPatrick Sanan .seealso: `DMPlexGetScale()`, `DMPlexSetScale()`
2639dc85fa5SMatthew G. Knepley E*/
2649371c9d4SSatish Balay typedef enum {
2659371c9d4SSatish Balay   PETSC_UNIT_LENGTH,
2669371c9d4SSatish Balay   PETSC_UNIT_MASS,
2679371c9d4SSatish Balay   PETSC_UNIT_TIME,
2689371c9d4SSatish Balay   PETSC_UNIT_CURRENT,
2699371c9d4SSatish Balay   PETSC_UNIT_TEMPERATURE,
2709371c9d4SSatish Balay   PETSC_UNIT_AMOUNT,
2719371c9d4SSatish Balay   PETSC_UNIT_LUMINOSITY,
2729371c9d4SSatish Balay   NUM_PETSC_UNITS
2739371c9d4SSatish Balay } PetscUnit;
2749dc85fa5SMatthew G. Knepley 
275b2b58855SToby Isaac /*S
276b2b58855SToby Isaac     DMField - PETSc object for defining a field on a mesh topology
277b2b58855SToby Isaac 
278b2b58855SToby Isaac     Level: intermediate
279b2b58855SToby Isaac S*/
280b2b58855SToby Isaac typedef struct _p_DMField *DMField;
281b2b58855SToby Isaac 
2820fdc7489SMatthew Knepley /*S
28387497f52SBarry Smith     DMUniversalLabel - A label that encodes a set of `DMLabel`s, bijectively
2840fdc7489SMatthew Knepley 
2850fdc7489SMatthew Knepley     Level: developer
2860fdc7489SMatthew Knepley S*/
2870fdc7489SMatthew Knepley typedef struct _p_UniversalLabel *DMUniversalLabel;
2880fdc7489SMatthew Knepley 
289c0517cd5SMatthew G. Knepley typedef struct _n_DMGeneratorFunctionList *DMGeneratorFunctionList;
290c0517cd5SMatthew G. Knepley 
2911e25c274SJed Brown #endif
292