Actual source code: petscfetypes.h

  1: #pragma once

  3: /* MANSEC = DM */
  4: /* SUBMANSEC = FE */

  6: /*S
  7:   PetscFE - PETSc object that manages a finite element space, e.g. the P_1 Lagrange element

  9:   Level: beginner

 11: .seealso: `PetscFECreate()`, `PetscSpace`, `PetscDualSpace`, `PetscSpaceCreate()`, `PetscDualSpaceCreate()`, `PetscFESetType()`, `PetscFEType`
 12: S*/
 13: typedef struct _p_PetscFE *PetscFE;

 15: /*MC
 16:   PetscFEJacobianType - indicates which pointwise functions should be used to fill the Jacobian matrix

 18:   Level: beginner

 20: .seealso: `PetscFEIntegrateJacobian()`
 21: M*/
 22: typedef enum {
 23:   PETSCFE_JACOBIAN,
 24:   PETSCFE_JACOBIAN_PRE,
 25:   PETSCFE_JACOBIAN_DYN
 26: } PetscFEJacobianType;

 28: /*E
 29:   PetscFEGeomMode - Describes the type of geometry being encoded.

 31:   Values:
 32: + `PETSC_FEGEOM_BASIC`    - These are normal dim-cells, with dim == dE, and only bulk data is stored.
 33: . `PETSC_FEGEOM_EMBEDDED` - These are dim-cells embedded in a higher dimension, as an embedded manifold, where dim < dE and only bulk data is stored.
 34: . `PETSC_FEGEOM_BOUNDARY` - These are dim-cells on the boundary of a dE-mesh, so that dim < dE, and both bulk and s = 1 face data are stored.
 35: - `PETSC_FEGEOM_COHESIVE` - These are dim-cells in the interior of a dE-mesh, so that dim < dE, and both bulk and s = 2 face data are stored.

 37:   Level: beginner

 39:   Note:
 40:   .vb
 41:   dim - The topological dimension and reference coordinate dimension
 42:   dE  - The real coordinate dimension
 43:   s   - The number of supporting cells for a face
 44:   .ve

 46: .seealso: [](ch_dmbase), `PetscFEGeom`, `DM`, `DMPLEX`, `PetscFEGeomCreate()`
 47: E*/
 48: typedef enum {
 49:   PETSC_FEGEOM_BASIC,
 50:   PETSC_FEGEOM_EMBEDDED,
 51:   PETSC_FEGEOM_BOUNDARY,
 52:   PETSC_FEGEOM_COHESIVE
 53: } PetscFEGeomMode;