Actual source code: petscdmtypes.h

petsc-3.11.4 2019-09-28
Report Typos and Errors
  1: #if !defined(PETSCDMTYPES_H)
  2: #define PETSCDMTYPES_H

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

  7:    Level: intermediate

  9:   Concepts: grids, grid refinement

 11:    Notes:
 12:     The DMDACreate() based object and the DMCompositeCreate() based object are examples of DMs

 14: .seealso:  DMCompositeCreate(), DMDACreate(), DMSetType(), DMType
 15: S*/
 16: typedef struct _p_DM* DM;

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

 21:   Level: beginner

 23:   A boundary may be of type DM_BOUNDARY_NONE (no ghost nodes), DM_BOUNDARY_GHOSTED (ghost vertices/cells
 24:   exist but aren't filled; you can put values into them and then apply a stencil that uses those ghost locations),
 25:   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; 
 26:   not yet implemented for 3d), DM_BOUNDARY_PERIODIC (ghost vertices/cells filled by the opposite
 27:   edge of the domain), or DM_BOUNDARY_TWIST (like periodic, only glued backwards like a Mobius strip).

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

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

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

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

 42: .seealso: DMDASetBoundaryType(), DMDACreate1d(), DMDACreate2d(), DMDACreate3d(), DMDACreate()
 43: E*/
 44: typedef enum {DM_BOUNDARY_NONE, DM_BOUNDARY_GHOSTED, DM_BOUNDARY_MIRROR, DM_BOUNDARY_PERIODIC, DM_BOUNDARY_TWIST} DMBoundaryType;
 45: /*E
 46:   DMBoundaryConditionType - indicates what type of boundary condition is to be imposed

 48:   Note: This flag indicates the type of function which will define the condition:
 49: $ DM_BC_ESSENTIAL       - A Dirichlet condition using a function of the coordinates
 50: $ DM_BC_ESSENTIAL_FIELD - A Dirichlet condition using a function of the coordinates and auxiliary field data
 51: $ DM_BC_NATURAL         - A Neumann condition using a function of the coordinates
 52: $ DM_BC_NATURAL_FIELD   - A Dirichlet condition using a function of the coordinates and auxiliary field data
 53: $ DM_BC_NATURAL_RIEMANN - A flux condition which determines the state in ghost cells
 54: The user can check whether a boundary condition is essential using (type & DM_BC_ESSENTIAL), and similarly for
 55: natural conditions (type & DM_BC_NATURAL)

 57:   Level: beginner

 59: .seealso: DMAddBoundary(), DMGetBoundary()
 60: E*/
 61: typedef enum {DM_BC_ESSENTIAL = 1, DM_BC_ESSENTIAL_FIELD = 5, DM_BC_NATURAL = 2, DM_BC_NATURAL_FIELD = 6, DM_BC_NATURAL_RIEMANN = 10} DMBoundaryConditionType;

 63: /*E
 64:   DMPointLocationType - Describes the method to handle point location failure

 66:   Level: beginner

 68:   If a search using DM_POINTLOCATION_NONE fails, the failure is signaled with a negative cell number. On the
 69:   other hand, if DM_POINTLOCATION_NEAREST is used, on failure, the (approximate) nearest point in the mesh is
 70:   used, replacing the given point in the input vector. DM_POINTLOCATION_REMOVE returns values only for points
 71:   which were located.

 73: .seealso: DMLocatePoints()
 74: E*/
 75: typedef enum {DM_POINTLOCATION_NONE, DM_POINTLOCATION_NEAREST, DM_POINTLOCATION_REMOVE} DMPointLocationType;

 77: /*E
 78:   DMAdaptationStrategy - Describes the strategy used for adaptive solves

 80:   Level: beginner

 82:   DM_ADAPTATION_INITIAL will refine a mesh based on an initial guess. DM_ADAPTATION_SEQUENTIAL will refine the
 83:   mesh based on a sequence of solves, much like grid sequencing. DM_ADAPTATION_MULTILEVEL will use the sequence
 84:   of constructed meshes in a multilevel solve, much like the Systematic Upscaling of Brandt.

 86: .seealso: DMAdaptorSolve()
 87: E*/
 88: typedef enum {DM_ADAPTATION_INITIAL, DM_ADAPTATION_SEQUENTIAL, DM_ADAPTATION_MULTILEVEL} DMAdaptationStrategy;

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

 93:   Level: beginner

 95:   DM_ADAPTATION_REFINE will uniformly refine a mesh, much like grid sequencing. DM_ADAPTATION_LABEL will adapt
 96:   the mesh based upon a label of the cells filled with DMAdaptFlag markers. DM_ADAPTATION_METRIC will try to
 97:   mesh the manifold described by the input metric tensor uniformly. PETSc can also construct such a metric based
 98:   upon an input primal or a gradient field.

100: .seealso: DMAdaptorSolve()
101: E*/
102: typedef enum {DM_ADAPTATION_NONE, DM_ADAPTATION_REFINE, DM_ADAPTATION_LABEL, DM_ADAPTATION_METRIC} DMAdaptationCriterion;

104: /*E
105:   DMAdaptFlag - Marker in the label prescribing adaptation

107:   Level: beginner

109: .seealso: DMAdaptLabel()
110: E*/
111: 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;

113: /*S
114:   PetscPartitioner - PETSc object that manages a graph partitioner

116:   Level: intermediate

118:   Concepts: partition, mesh

120: .seealso: PetscPartitionerCreate(), PetscPartitionerSetType(), PetscPartitionerType
121: S*/
122: typedef struct _p_PetscPartitioner *PetscPartitioner;

124: /*E
125:   PetscUnit - The seven fundamental SI units

127:   Level: beginner

129: .seealso: DMPlexGetScale(), DMPlexSetScale()
130: E*/
131: 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;

133: /*S
134:     DMField - PETSc object for defining a field on a mesh topology

136:     Level: intermediate
137: S*/
138: typedef struct _p_DMField* DMField;

140: #endif