126bd1501SBarry Smith #if !defined(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 11*87497f52SBarry Smith .seealso: `DMType`, `DMDGetType()`, `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 20*87497f52SBarry 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), 22*87497f52SBarry 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; 23*87497f52SBarry Smith not yet implemented for 3d), `DM_BOUNDARY_PERIODIC` (ghost vertices/cells filled by the opposite 24*87497f52SBarry 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 28*87497f52SBarry Smith processes. That width is always determined by the stencil width; see `DMDASetStencilWidth()`. 29bff4a2f0SMatthew G. Knepley 30*87497f52SBarry 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: 33*87497f52SBarry 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*/ 41bff4a2f0SMatthew G. Knepley typedef enum {DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_MIRROR, DM_BOUNDARY_PERIODIC, DM_BOUNDARY_TWIST} DMBoundaryType; 4262a38674SMatthew G. Knepley /*E 439dc85fa5SMatthew G. Knepley DMBoundaryConditionType - indicates what type of boundary condition is to be imposed 449dc85fa5SMatthew G. Knepley 459dc85fa5SMatthew G. Knepley Note: This flag indicates the type of function which will define the condition: 469dc85fa5SMatthew G. Knepley $ DM_BC_ESSENTIAL - A Dirichlet condition using a function of the coordinates 479dc85fa5SMatthew G. Knepley $ DM_BC_ESSENTIAL_FIELD - A Dirichlet condition using a function of the coordinates and auxiliary field data 48b18799e0SMatthew G. Knepley $ DM_BC_ESSENTIAL_BD_FIELD - A Dirichlet condition using a function of the coordinates, facet normal, and auxiliary field data 499dc85fa5SMatthew G. Knepley $ DM_BC_NATURAL - A Neumann condition using a function of the coordinates 50b18799e0SMatthew G. Knepley $ DM_BC_NATURAL_FIELD - A Neumann condition using a function of the coordinates and auxiliary field data 519dc85fa5SMatthew G. Knepley $ DM_BC_NATURAL_RIEMANN - A flux condition which determines the state in ghost cells 529dc85fa5SMatthew G. Knepley The user can check whether a boundary condition is essential using (type & DM_BC_ESSENTIAL), and similarly for 539dc85fa5SMatthew G. Knepley natural conditions (type & DM_BC_NATURAL) 549dc85fa5SMatthew G. Knepley 559dc85fa5SMatthew G. Knepley Level: beginner 569dc85fa5SMatthew G. Knepley 57db781477SPatrick Sanan .seealso: `DMAddBoundary()`, `DSAddBoundary()`, `DSGetBoundary()` 589dc85fa5SMatthew G. Knepley E*/ 59b18799e0SMatthew G. Knepley typedef enum {DM_BC_ESSENTIAL = 1, DM_BC_ESSENTIAL_FIELD = 5, DM_BC_NATURAL = 2, DM_BC_NATURAL_FIELD = 6, DM_BC_ESSENTIAL_BD_FIELD = 9, DM_BC_NATURAL_RIEMANN = 10} DMBoundaryConditionType; 609dc85fa5SMatthew G. Knepley 619dc85fa5SMatthew G. Knepley /*E 6262a38674SMatthew G. Knepley DMPointLocationType - Describes the method to handle point location failure 6362a38674SMatthew G. Knepley 6462a38674SMatthew G. Knepley Level: beginner 6562a38674SMatthew G. Knepley 66*87497f52SBarry Smith If a search using `DM_POINTLOCATION_NONE` fails, the failure is signaled with a negative cell number. On the 67*87497f52SBarry Smith other hand, if `DM_POINTLOCATION_NEAREST` is used, on failure, the (approximate) nearest point in the mesh is 68*87497f52SBarry Smith used, replacing the given point in the input vector. `DM_POINTLOCATION_REMOVE` returns values only for points 692d1fa6caSMatthew G. Knepley which were located. 7062a38674SMatthew G. Knepley 71db781477SPatrick Sanan .seealso: `DMLocatePoints()` 7262a38674SMatthew G. Knepley E*/ 732d1fa6caSMatthew G. Knepley typedef enum {DM_POINTLOCATION_NONE, DM_POINTLOCATION_NEAREST, DM_POINTLOCATION_REMOVE} DMPointLocationType; 7462a38674SMatthew G. Knepley 755675c177SMatthew G. Knepley /*E 76174e7490SMatthew G. Knepley DMAdaptationStrategy - Describes the strategy used for adaptive solves 775675c177SMatthew G. Knepley 785675c177SMatthew G. Knepley Level: beginner 795675c177SMatthew G. Knepley 8059b28e79SMatthew G. Knepley DM_ADAPTATION_INITIAL will refine a mesh based on an initial guess. DM_ADAPTATION_SEQUENTIAL will refine the 8159b28e79SMatthew G. Knepley mesh based on a sequence of solves, much like grid sequencing. DM_ADAPTATION_MULTILEVEL will use the sequence 8259b28e79SMatthew G. Knepley of constructed meshes in a multilevel solve, much like the Systematic Upscaling of Brandt. 835675c177SMatthew G. Knepley 84db781477SPatrick Sanan .seealso: `DMAdaptorSolve()` 855675c177SMatthew G. Knepley E*/ 86174e7490SMatthew G. Knepley typedef enum {DM_ADAPTATION_INITIAL, DM_ADAPTATION_SEQUENTIAL, DM_ADAPTATION_MULTILEVEL} DMAdaptationStrategy; 87174e7490SMatthew G. Knepley 88174e7490SMatthew G. Knepley /*E 89174e7490SMatthew G. Knepley DMAdaptationCriterion - Describes the test used to decide whether to coarsen or refine parts of the mesh 90174e7490SMatthew G. Knepley 91174e7490SMatthew G. Knepley Level: beginner 92174e7490SMatthew G. Knepley 93*87497f52SBarry Smith `DM_ADAPTATION_REFINE` will uniformly refine a mesh, much like grid sequencing. `DM_ADAPTATION_LABEL` will adapt 94*87497f52SBarry Smith the mesh based upon a label of the cells filled with `DMAdaptFlag` markers. `DM_ADAPTATION_METRIC` will try to 95174e7490SMatthew G. Knepley mesh the manifold described by the input metric tensor uniformly. PETSc can also construct such a metric based 96174e7490SMatthew G. Knepley upon an input primal or a gradient field. 97174e7490SMatthew G. Knepley 98db781477SPatrick Sanan .seealso: `DMAdaptorSolve()` 99174e7490SMatthew G. Knepley E*/ 100174e7490SMatthew G. Knepley typedef enum {DM_ADAPTATION_NONE, DM_ADAPTATION_REFINE, DM_ADAPTATION_LABEL, DM_ADAPTATION_METRIC} DMAdaptationCriterion; 1015675c177SMatthew G. Knepley 1029dc85fa5SMatthew G. Knepley /*E 1039dc85fa5SMatthew G. Knepley DMAdaptFlag - Marker in the label prescribing adaptation 1049dc85fa5SMatthew G. Knepley 1059dc85fa5SMatthew G. Knepley Level: beginner 1069dc85fa5SMatthew G. Knepley 107db781477SPatrick Sanan .seealso: `DMAdaptLabel()` 1089dc85fa5SMatthew G. Knepley E*/ 109bf2d5fbbSStefano Zampini typedef enum {DM_ADAPT_DETERMINE = PETSC_DETERMINE, DM_ADAPT_KEEP = 0, DM_ADAPT_REFINE, DM_ADAPT_COARSEN, DM_ADAPT_COARSEN_LAST, DM_ADAPT_RESERVED_COUNT} DMAdaptFlag; 1103ee9839eSMatthew G. Knepley 1113ee9839eSMatthew G. Knepley /*E 1123ee9839eSMatthew G. Knepley DMDirection - Indicates a coordinate direction 1133ee9839eSMatthew G. Knepley 1143ee9839eSMatthew G. Knepley Level: beginner 1153ee9839eSMatthew G. Knepley 116db781477SPatrick Sanan .seealso: `DMDAGetRay()`, `DMDAGetProcessorSubset()`, `DMPlexShearGeometry()` 1173ee9839eSMatthew G. Knepley E*/ 1183ee9839eSMatthew G. Knepley typedef enum {DM_X, DM_Y, DM_Z} DMDirection; 1199dc85fa5SMatthew G. Knepley 120a6e0b375SMatthew G. Knepley /*E 121*87497f52SBarry Smith DMEnclosureType - The type of enclosure relation between one `DM` and another 122a6e0b375SMatthew G. Knepley 123a6e0b375SMatthew G. Knepley Level: beginner 124a6e0b375SMatthew G. Knepley 125*87497f52SBarry Smith Notes: 126*87497f52SBarry Smith For example, one `DM` dmA may be the boundary of another dmB, in which case it would be labeled `DM_ENC_SUBMESH`. 127*87497f52SBarry Smith 128*87497f52SBarry Smith If the situation is reversed, and dmA has boundary dmB, it would be labeled `DM_ENC_SUPERMESH`. 129*87497f52SBarry Smith 130*87497f52SBarry Smith Likewise, if dmA was a subregion of dmB, it would be labeled `DM_ENC_SUBMESH`. 131*87497f52SBarry Smith 132*87497f52SBarry Smith If no relation can be determined, `DM_ENC_NONE` is used. 133*87497f52SBarry Smith 134*87497f52SBarry Smith If a relation is not yet known, then `DM_ENC_UNKNOWN` is used. 135a6e0b375SMatthew G. Knepley 136db781477SPatrick Sanan .seealso: `DMGetEnclosureRelation()` 137a6e0b375SMatthew G. Knepley E*/ 138a6e0b375SMatthew G. Knepley typedef enum {DM_ENC_EQUALITY, DM_ENC_SUPERMESH, DM_ENC_SUBMESH, DM_ENC_NONE, DM_ENC_UNKNOWN} DMEnclosureType; 139ba2698f1SMatthew G. Knepley 140ba2698f1SMatthew G. Knepley /*E 141ba2698f1SMatthew G. Knepley DMPolytopeType - This describes the polytope represented by each cell. 142ba2698f1SMatthew G. Knepley 143ba2698f1SMatthew G. Knepley Level: beginner 144ba2698f1SMatthew G. Knepley 145*87497f52SBarry Smith While most operations only need the topology information in the `DMPLEX`, we must sometimes have the 146ba2698f1SMatthew G. Knepley user specify a polytope. For instance, when interpolating from a cell-vertex mesh, the type of 147*87497f52SBarry Smith polytope can be ambiguous. Also, `DMPLEX` allows different symmetries of a prism cell with the same 148ba2698f1SMatthew G. Knepley constituent points. Normally these types are autoamtically inferred and the user does not specify 149ba2698f1SMatthew G. Knepley them. 150ba2698f1SMatthew G. Knepley 151db781477SPatrick Sanan .seealso: `DMPlexComputeCellTypes()` 152ba2698f1SMatthew G. Knepley E*/ 153da9060c4SMatthew G. Knepley typedef enum {DM_POLYTOPE_POINT, DM_POLYTOPE_SEGMENT, DM_POLYTOPE_POINT_PRISM_TENSOR, DM_POLYTOPE_TRIANGLE, DM_POLYTOPE_QUADRILATERAL, DM_POLYTOPE_SEG_PRISM_TENSOR, DM_POLYTOPE_TETRAHEDRON, DM_POLYTOPE_HEXAHEDRON, DM_POLYTOPE_TRI_PRISM, DM_POLYTOPE_TRI_PRISM_TENSOR, DM_POLYTOPE_QUAD_PRISM_TENSOR, DM_POLYTOPE_PYRAMID, DM_POLYTOPE_FV_GHOST, DM_POLYTOPE_INTERIOR_GHOST, DM_POLYTOPE_UNKNOWN, DM_NUM_POLYTOPES} DMPolytopeType; 154ba2698f1SMatthew G. Knepley PETSC_EXTERN const char *const DMPolytopeTypes[]; 155a6e0b375SMatthew G. Knepley 1569dc85fa5SMatthew G. Knepley /*E 1579dc85fa5SMatthew G. Knepley PetscUnit - The seven fundamental SI units 1589dc85fa5SMatthew G. Knepley 1599dc85fa5SMatthew G. Knepley Level: beginner 1609dc85fa5SMatthew G. Knepley 161db781477SPatrick Sanan .seealso: `DMPlexGetScale()`, `DMPlexSetScale()` 1629dc85fa5SMatthew G. Knepley E*/ 1639dc85fa5SMatthew G. Knepley typedef enum {PETSC_UNIT_LENGTH, PETSC_UNIT_MASS, PETSC_UNIT_TIME, PETSC_UNIT_CURRENT, PETSC_UNIT_TEMPERATURE, PETSC_UNIT_AMOUNT, PETSC_UNIT_LUMINOSITY, NUM_PETSC_UNITS} PetscUnit; 1649dc85fa5SMatthew G. Knepley 165b2b58855SToby Isaac /*S 166b2b58855SToby Isaac DMField - PETSc object for defining a field on a mesh topology 167b2b58855SToby Isaac 168b2b58855SToby Isaac Level: intermediate 169b2b58855SToby Isaac S*/ 170b2b58855SToby Isaac typedef struct _p_DMField* DMField; 171b2b58855SToby Isaac 1720fdc7489SMatthew Knepley /*S 173*87497f52SBarry Smith DMUniversalLabel - A label that encodes a set of `DMLabel`s, bijectively 1740fdc7489SMatthew Knepley 1750fdc7489SMatthew Knepley Level: developer 1760fdc7489SMatthew Knepley S*/ 1770fdc7489SMatthew Knepley typedef struct _p_UniversalLabel* DMUniversalLabel; 1780fdc7489SMatthew Knepley 179c0517cd5SMatthew G. Knepley typedef struct _n_DMGeneratorFunctionList *DMGeneratorFunctionList; 180c0517cd5SMatthew G. Knepley 1811e25c274SJed Brown #endif 182