16524c165SJacob Faibussowitsch #ifndef PETSCDMTYPES_H 226bd1501SBarry Smith #define PETSCDMTYPES_H 31e25c274SJed Brown 4ac09b921SBarry Smith /* SUBMANSEC = DM */ 5ac09b921SBarry Smith 61e25c274SJed Brown /*S 71e25c274SJed Brown DM - Abstract PETSc object that manages an abstract grid object and its interactions with the algebraic solvers 81e25c274SJed Brown 91e25c274SJed Brown Level: intermediate 101e25c274SJed Brown 11dce8aebaSBarry Smith .seealso: `DMType`, `DMGetType()`, `DMCompositeCreate()`, `DMDACreate()`, `DMSetType()`, `DMType`, `DMDA`, `DMPLEX` 121e25c274SJed Brown S*/ 131e25c274SJed Brown typedef struct _p_DM *DM; 141e25c274SJed Brown 15bff4a2f0SMatthew G. Knepley /*E 16bff4a2f0SMatthew G. Knepley DMBoundaryType - Describes the choice for fill of ghost cells on physical domain boundaries. 17bff4a2f0SMatthew G. Knepley 18bff4a2f0SMatthew G. Knepley Level: beginner 19bff4a2f0SMatthew G. Knepley 2087497f52SBarry Smith A boundary may be of type `DM_BOUNDARY_NONE` (no ghost nodes), `DM_BOUNDARY_GHOSTED` (ghost vertices/cells 21dbb368e6SPatrick Sanan exist but aren't filled; you can put values into them and then apply a stencil that uses those ghost locations), 2287497f52SBarry 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 the ghost point value; 2387497f52SBarry Smith not yet implemented for 3d), `DM_BOUNDARY_PERIODIC` (ghost vertices/cells filled by the opposite 2487497f52SBarry Smith edge of the domain), or `DM_BOUNDARY_TWIST` (like periodic, only glued backwards like a Mobius strip). 25bff4a2f0SMatthew G. Knepley 26dbb368e6SPatrick Sanan Notes: 27dbb368e6SPatrick Sanan This is information for the boundary of the __PHYSICAL__ domain. It has nothing to do with boundaries between 2887497f52SBarry Smith processes. That width is always determined by the stencil width; see `DMDASetStencilWidth()`. 29bff4a2f0SMatthew G. Knepley 3087497f52SBarry 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 . 31288e7d53SBarry Smith 3295452b02SPatrick Sanan Developer Notes: 3387497f52SBarry 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 34288e7d53SBarry Smith as the 0th grid point where the physical boundary serves as the mirror? 35288e7d53SBarry Smith 36dbb368e6SPatrick Sanan References: 37606c0280SSatish Balay . * - https://scicomp.stackexchange.com/questions/5355/writing-the-poisson-equation-finite-difference-matrix-with-neumann-boundary-cond 38288e7d53SBarry Smith 39db781477SPatrick Sanan .seealso: `DMDASetBoundaryType()`, `DMDACreate1d()`, `DMDACreate2d()`, `DMDACreate3d()`, `DMDACreate()` 40bff4a2f0SMatthew G. Knepley E*/ 419371c9d4SSatish Balay typedef enum { 429371c9d4SSatish Balay DM_BOUNDARY_NONE, 439371c9d4SSatish Balay DM_BOUNDARY_GHOSTED, 449371c9d4SSatish Balay DM_BOUNDARY_MIRROR, 459371c9d4SSatish Balay DM_BOUNDARY_PERIODIC, 469371c9d4SSatish Balay DM_BOUNDARY_TWIST 479371c9d4SSatish Balay } DMBoundaryType; 4862a38674SMatthew G. Knepley /*E 499dc85fa5SMatthew G. Knepley DMBoundaryConditionType - indicates what type of boundary condition is to be imposed 509dc85fa5SMatthew G. Knepley 519dc85fa5SMatthew G. Knepley Note: This flag indicates the type of function which will define the condition: 529dc85fa5SMatthew G. Knepley $ DM_BC_ESSENTIAL - A Dirichlet condition using a function of the coordinates 539dc85fa5SMatthew G. Knepley $ DM_BC_ESSENTIAL_FIELD - A Dirichlet condition using a function of the coordinates and auxiliary field data 54b18799e0SMatthew G. Knepley $ DM_BC_ESSENTIAL_BD_FIELD - A Dirichlet condition using a function of the coordinates, facet normal, and auxiliary field data 559dc85fa5SMatthew G. Knepley $ DM_BC_NATURAL - A Neumann condition using a function of the coordinates 56b18799e0SMatthew G. Knepley $ DM_BC_NATURAL_FIELD - A Neumann condition using a function of the coordinates and auxiliary field data 579dc85fa5SMatthew G. Knepley $ DM_BC_NATURAL_RIEMANN - A flux condition which determines the state in ghost cells 589dc85fa5SMatthew G. Knepley The user can check whether a boundary condition is essential using (type & DM_BC_ESSENTIAL), and similarly for 599dc85fa5SMatthew G. Knepley natural conditions (type & DM_BC_NATURAL) 609dc85fa5SMatthew G. Knepley 619dc85fa5SMatthew G. Knepley Level: beginner 629dc85fa5SMatthew G. Knepley 63db781477SPatrick Sanan .seealso: `DMAddBoundary()`, `DSAddBoundary()`, `DSGetBoundary()` 649dc85fa5SMatthew G. Knepley E*/ 659371c9d4SSatish Balay typedef enum { 669371c9d4SSatish Balay DM_BC_ESSENTIAL = 1, 679371c9d4SSatish Balay DM_BC_ESSENTIAL_FIELD = 5, 689371c9d4SSatish Balay DM_BC_NATURAL = 2, 699371c9d4SSatish Balay DM_BC_NATURAL_FIELD = 6, 709371c9d4SSatish Balay DM_BC_ESSENTIAL_BD_FIELD = 9, 719371c9d4SSatish Balay DM_BC_NATURAL_RIEMANN = 10 729371c9d4SSatish Balay } DMBoundaryConditionType; 739dc85fa5SMatthew G. Knepley 749dc85fa5SMatthew G. Knepley /*E 7562a38674SMatthew G. Knepley DMPointLocationType - Describes the method to handle point location failure 7662a38674SMatthew G. Knepley 7762a38674SMatthew G. Knepley Level: beginner 7862a38674SMatthew G. Knepley 7987497f52SBarry Smith If a search using `DM_POINTLOCATION_NONE` fails, the failure is signaled with a negative cell number. On the 8087497f52SBarry Smith other hand, if `DM_POINTLOCATION_NEAREST` is used, on failure, the (approximate) nearest point in the mesh is 8187497f52SBarry Smith used, replacing the given point in the input vector. `DM_POINTLOCATION_REMOVE` returns values only for points 822d1fa6caSMatthew G. Knepley which were located. 8362a38674SMatthew G. Knepley 84db781477SPatrick Sanan .seealso: `DMLocatePoints()` 8562a38674SMatthew G. Knepley E*/ 869371c9d4SSatish Balay typedef enum { 879371c9d4SSatish Balay DM_POINTLOCATION_NONE, 889371c9d4SSatish Balay DM_POINTLOCATION_NEAREST, 899371c9d4SSatish Balay DM_POINTLOCATION_REMOVE 909371c9d4SSatish Balay } DMPointLocationType; 9162a38674SMatthew G. Knepley 925675c177SMatthew G. Knepley /*E 93*863027abSJed Brown DMBlockingType - Describes how to choose variable block sizes 94*863027abSJed Brown 95*863027abSJed Brown Level: intermediate 96*863027abSJed Brown 97*863027abSJed Brown When using `PCVPBJACOBI`, one can choose to block by topological point (all fields at a cell center, at a face, etc.) 98*863027abSJed Brown or by field nodes (using number of components per field to identify "nodes"). Field nodes lead to smaller blocks, but 99*863027abSJed Brown may converge more slowly. For example, a cubic Lagrange hexahedron will have one node at vertices, two at edges, four 100*863027abSJed Brown at faces, and eight at cell centers. If using point blocking, the `PCVPBJACOBI` preconditioner will work with block 101*863027abSJed Brown sizes up to 8 Lagrange nodes. For 5-component CFD, this produces matrices up to 40x40, which increases memory 102*863027abSJed Brown footprint and may harm performance. With field node blocking, the max block size will correspond to one Lagrange node, 103*863027abSJed Brown or 5x5 blocks for the CFD example. 104*863027abSJed Brown 105*863027abSJed Brown .seealso: `PCVPBJACOBI`, `MatSetVariableBlockSizes()`, `DMSetBlockingType()` 106*863027abSJed Brown E*/ 107*863027abSJed Brown typedef enum { 108*863027abSJed Brown DM_BLOCKING_POINT, 109*863027abSJed Brown DM_BLOCKING_FIELD_NODE 110*863027abSJed Brown } DMBlockingType; 111*863027abSJed Brown 112*863027abSJed Brown /*E 113174e7490SMatthew G. Knepley DMAdaptationStrategy - Describes the strategy used for adaptive solves 1145675c177SMatthew G. Knepley 1155675c177SMatthew G. Knepley Level: beginner 1165675c177SMatthew G. Knepley 11759b28e79SMatthew G. Knepley DM_ADAPTATION_INITIAL will refine a mesh based on an initial guess. DM_ADAPTATION_SEQUENTIAL will refine the 11859b28e79SMatthew G. Knepley mesh based on a sequence of solves, much like grid sequencing. DM_ADAPTATION_MULTILEVEL will use the sequence 11959b28e79SMatthew G. Knepley of constructed meshes in a multilevel solve, much like the Systematic Upscaling of Brandt. 1205675c177SMatthew G. Knepley 121db781477SPatrick Sanan .seealso: `DMAdaptorSolve()` 1225675c177SMatthew G. Knepley E*/ 1239371c9d4SSatish Balay typedef enum { 1249371c9d4SSatish Balay DM_ADAPTATION_INITIAL, 1259371c9d4SSatish Balay DM_ADAPTATION_SEQUENTIAL, 1269371c9d4SSatish Balay DM_ADAPTATION_MULTILEVEL 1279371c9d4SSatish Balay } DMAdaptationStrategy; 128174e7490SMatthew G. Knepley 129174e7490SMatthew G. Knepley /*E 130174e7490SMatthew G. Knepley DMAdaptationCriterion - Describes the test used to decide whether to coarsen or refine parts of the mesh 131174e7490SMatthew G. Knepley 132174e7490SMatthew G. Knepley Level: beginner 133174e7490SMatthew G. Knepley 13487497f52SBarry Smith `DM_ADAPTATION_REFINE` will uniformly refine a mesh, much like grid sequencing. `DM_ADAPTATION_LABEL` will adapt 13587497f52SBarry Smith the mesh based upon a label of the cells filled with `DMAdaptFlag` markers. `DM_ADAPTATION_METRIC` will try to 136174e7490SMatthew G. Knepley mesh the manifold described by the input metric tensor uniformly. PETSc can also construct such a metric based 137174e7490SMatthew G. Knepley upon an input primal or a gradient field. 138174e7490SMatthew G. Knepley 139db781477SPatrick Sanan .seealso: `DMAdaptorSolve()` 140174e7490SMatthew G. Knepley E*/ 1419371c9d4SSatish Balay typedef enum { 1429371c9d4SSatish Balay DM_ADAPTATION_NONE, 1439371c9d4SSatish Balay DM_ADAPTATION_REFINE, 1449371c9d4SSatish Balay DM_ADAPTATION_LABEL, 1459371c9d4SSatish Balay DM_ADAPTATION_METRIC 1469371c9d4SSatish Balay } DMAdaptationCriterion; 1475675c177SMatthew G. Knepley 1489dc85fa5SMatthew G. Knepley /*E 1499dc85fa5SMatthew G. Knepley DMAdaptFlag - Marker in the label prescribing adaptation 1509dc85fa5SMatthew G. Knepley 1519dc85fa5SMatthew G. Knepley Level: beginner 1529dc85fa5SMatthew G. Knepley 153db781477SPatrick Sanan .seealso: `DMAdaptLabel()` 1549dc85fa5SMatthew G. Knepley E*/ 1559371c9d4SSatish Balay typedef enum { 1569371c9d4SSatish Balay DM_ADAPT_DETERMINE = PETSC_DETERMINE, 1579371c9d4SSatish Balay DM_ADAPT_KEEP = 0, 1589371c9d4SSatish Balay DM_ADAPT_REFINE, 1599371c9d4SSatish Balay DM_ADAPT_COARSEN, 1609371c9d4SSatish Balay DM_ADAPT_COARSEN_LAST, 1619371c9d4SSatish Balay DM_ADAPT_RESERVED_COUNT 1629371c9d4SSatish Balay } DMAdaptFlag; 1633ee9839eSMatthew G. Knepley 1643ee9839eSMatthew G. Knepley /*E 1653ee9839eSMatthew G. Knepley DMDirection - Indicates a coordinate direction 1663ee9839eSMatthew G. Knepley 1673ee9839eSMatthew G. Knepley Level: beginner 1683ee9839eSMatthew G. Knepley 169db781477SPatrick Sanan .seealso: `DMDAGetRay()`, `DMDAGetProcessorSubset()`, `DMPlexShearGeometry()` 1703ee9839eSMatthew G. Knepley E*/ 1719371c9d4SSatish Balay typedef enum { 1729371c9d4SSatish Balay DM_X, 1739371c9d4SSatish Balay DM_Y, 1749371c9d4SSatish Balay DM_Z 1759371c9d4SSatish Balay } DMDirection; 1769dc85fa5SMatthew G. Knepley 177a6e0b375SMatthew G. Knepley /*E 17887497f52SBarry Smith DMEnclosureType - The type of enclosure relation between one `DM` and another 179a6e0b375SMatthew G. Knepley 180a6e0b375SMatthew G. Knepley Level: beginner 181a6e0b375SMatthew G. Knepley 18287497f52SBarry Smith Notes: 18387497f52SBarry Smith For example, one `DM` dmA may be the boundary of another dmB, in which case it would be labeled `DM_ENC_SUBMESH`. 18487497f52SBarry Smith 18587497f52SBarry Smith If the situation is reversed, and dmA has boundary dmB, it would be labeled `DM_ENC_SUPERMESH`. 18687497f52SBarry Smith 18787497f52SBarry Smith Likewise, if dmA was a subregion of dmB, it would be labeled `DM_ENC_SUBMESH`. 18887497f52SBarry Smith 18987497f52SBarry Smith If no relation can be determined, `DM_ENC_NONE` is used. 19087497f52SBarry Smith 19187497f52SBarry Smith If a relation is not yet known, then `DM_ENC_UNKNOWN` is used. 192a6e0b375SMatthew G. Knepley 193db781477SPatrick Sanan .seealso: `DMGetEnclosureRelation()` 194a6e0b375SMatthew G. Knepley E*/ 1959371c9d4SSatish Balay typedef enum { 1969371c9d4SSatish Balay DM_ENC_EQUALITY, 1979371c9d4SSatish Balay DM_ENC_SUPERMESH, 1989371c9d4SSatish Balay DM_ENC_SUBMESH, 1999371c9d4SSatish Balay DM_ENC_NONE, 2009371c9d4SSatish Balay DM_ENC_UNKNOWN 2019371c9d4SSatish Balay } DMEnclosureType; 202ba2698f1SMatthew G. Knepley 203ba2698f1SMatthew G. Knepley /*E 204ba2698f1SMatthew G. Knepley DMPolytopeType - This describes the polytope represented by each cell. 205ba2698f1SMatthew G. Knepley 206ba2698f1SMatthew G. Knepley Level: beginner 207ba2698f1SMatthew G. Knepley 20887497f52SBarry Smith While most operations only need the topology information in the `DMPLEX`, we must sometimes have the 209ba2698f1SMatthew G. Knepley user specify a polytope. For instance, when interpolating from a cell-vertex mesh, the type of 21087497f52SBarry Smith polytope can be ambiguous. Also, `DMPLEX` allows different symmetries of a prism cell with the same 211ba2698f1SMatthew G. Knepley constituent points. Normally these types are autoamtically inferred and the user does not specify 212ba2698f1SMatthew G. Knepley them. 213ba2698f1SMatthew G. Knepley 214db781477SPatrick Sanan .seealso: `DMPlexComputeCellTypes()` 215ba2698f1SMatthew G. Knepley E*/ 2169371c9d4SSatish Balay typedef enum { 2179371c9d4SSatish Balay DM_POLYTOPE_POINT, 2189371c9d4SSatish Balay DM_POLYTOPE_SEGMENT, 2199371c9d4SSatish Balay DM_POLYTOPE_POINT_PRISM_TENSOR, 2209371c9d4SSatish Balay DM_POLYTOPE_TRIANGLE, 2219371c9d4SSatish Balay DM_POLYTOPE_QUADRILATERAL, 2229371c9d4SSatish Balay DM_POLYTOPE_SEG_PRISM_TENSOR, 2239371c9d4SSatish Balay DM_POLYTOPE_TETRAHEDRON, 2249371c9d4SSatish Balay DM_POLYTOPE_HEXAHEDRON, 2259371c9d4SSatish Balay DM_POLYTOPE_TRI_PRISM, 2269371c9d4SSatish Balay DM_POLYTOPE_TRI_PRISM_TENSOR, 2279371c9d4SSatish Balay DM_POLYTOPE_QUAD_PRISM_TENSOR, 2289371c9d4SSatish Balay DM_POLYTOPE_PYRAMID, 2299371c9d4SSatish Balay DM_POLYTOPE_FV_GHOST, 2309371c9d4SSatish Balay DM_POLYTOPE_INTERIOR_GHOST, 2319371c9d4SSatish Balay DM_POLYTOPE_UNKNOWN, 2329371c9d4SSatish Balay DM_NUM_POLYTOPES 2339371c9d4SSatish Balay } DMPolytopeType; 234ba2698f1SMatthew G. Knepley PETSC_EXTERN const char *const DMPolytopeTypes[]; 235a6e0b375SMatthew G. Knepley 2369dc85fa5SMatthew G. Knepley /*E 2379dc85fa5SMatthew G. Knepley PetscUnit - The seven fundamental SI units 2389dc85fa5SMatthew G. Knepley 2399dc85fa5SMatthew G. Knepley Level: beginner 2409dc85fa5SMatthew G. Knepley 241db781477SPatrick Sanan .seealso: `DMPlexGetScale()`, `DMPlexSetScale()` 2429dc85fa5SMatthew G. Knepley E*/ 2439371c9d4SSatish Balay typedef enum { 2449371c9d4SSatish Balay PETSC_UNIT_LENGTH, 2459371c9d4SSatish Balay PETSC_UNIT_MASS, 2469371c9d4SSatish Balay PETSC_UNIT_TIME, 2479371c9d4SSatish Balay PETSC_UNIT_CURRENT, 2489371c9d4SSatish Balay PETSC_UNIT_TEMPERATURE, 2499371c9d4SSatish Balay PETSC_UNIT_AMOUNT, 2509371c9d4SSatish Balay PETSC_UNIT_LUMINOSITY, 2519371c9d4SSatish Balay NUM_PETSC_UNITS 2529371c9d4SSatish Balay } PetscUnit; 2539dc85fa5SMatthew G. Knepley 254b2b58855SToby Isaac /*S 255b2b58855SToby Isaac DMField - PETSc object for defining a field on a mesh topology 256b2b58855SToby Isaac 257b2b58855SToby Isaac Level: intermediate 258b2b58855SToby Isaac S*/ 259b2b58855SToby Isaac typedef struct _p_DMField *DMField; 260b2b58855SToby Isaac 2610fdc7489SMatthew Knepley /*S 26287497f52SBarry Smith DMUniversalLabel - A label that encodes a set of `DMLabel`s, bijectively 2630fdc7489SMatthew Knepley 2640fdc7489SMatthew Knepley Level: developer 2650fdc7489SMatthew Knepley S*/ 2660fdc7489SMatthew Knepley typedef struct _p_UniversalLabel *DMUniversalLabel; 2670fdc7489SMatthew Knepley 268c0517cd5SMatthew G. Knepley typedef struct _n_DMGeneratorFunctionList *DMGeneratorFunctionList; 269c0517cd5SMatthew G. Knepley 2701e25c274SJed Brown #endif 271