Actual source code: petscdualspace.h

  1: /*
  2:       Objects which encapsulate finite element spaces
  3: */
  4: #ifndef PETSCDUALSPACE_H
  5: #define PETSCDUALSPACE_H
  6: #include <petscdm.h>
  7: #include <petscdt.h>
  8: #include <petscfetypes.h>
  9: #include <petscdstypes.h>
 10: #include <petscspace.h>

 12: /* SUBMANSEC = DUALSPACE */

 14: /*S
 15:   PetscDualSpace - PETSc object that manages the dual space to a linear space, e.g. the space of evaluation functionals at the vertices of a triangle

 17:   Level: beginner

 19: .seealso: `PetscDualSpaceCreate()`, `PetscSpace`, `PetscSpaceCreate()`, `PetscDualSpaceSetType()`, `PetscDualSpaceType`
 20: S*/
 21: typedef struct _p_PetscDualSpace *PetscDualSpace;

 23: /*MC
 24:   PetscDualSpaceReferenceCell - The type of reference cell

 26:   Level: beginner

 28:   Note:
 29:   This is used only for automatic creation of reference cells. A `PetscDualSpace` can accept an arbitrary `DM` for a reference cell.

 31: .seealso: `PetscSpace`, `PetscDualSpaceCreate()`, `PetscDualSpaceType`
 32: M*/
 33: typedef enum {
 34:   PETSCDUALSPACE_REFCELL_SIMPLEX,
 35:   PETSCDUALSPACE_REFCELL_TENSOR
 36: } PetscDualSpaceReferenceCell;
 37: PETSC_EXTERN const char *const PetscDualSpaceReferenceCells[];

 39: /*MC
 40:   PetscDualSpaceTransformType - The type of function transform

 42:   Level: intermediate

 44:   Values:
 45: +  `IDENTITY_TRANSFORM` - make no changes in the function
 46: .  `COVARIANT_PIOLA_TRANSFORM` - Covariant Piola: $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$
 47: -  `CONTRAVARIANT_PIOLA_TRANSFORM` - Contravariant Piola: $\sigma^*(F) = 1/|J| J F \circ \phi^{-1)$

 49:   Notes:
 50:   These transforms, and their inverses, are used to move functions and functionals between the reference element and real space.
 51:   Suppose that we have a mapping $\phi$ which maps the reference cell to real space, and its Jacobian $J$. If we want to transform function $F$ on the reference element,
 52:   so that it acts on real space, we use the pushforward transform $\sigma^*$. The pullback $\sigma_*$ is the inverse transform.

 54:   References:
 55: .    Rognes, Kirby, and Logg, Efficient Assembly of Hdiv and Hrot Conforming Finite Elements, SISC, 31(6), 4130-4151, arXiv 1205.3085, 2010

 57: .seealso: `PetscDualSpaceGetDeRahm()`
 58: M*/
 59: typedef enum {
 60:   IDENTITY_TRANSFORM,
 61:   COVARIANT_PIOLA_TRANSFORM,
 62:   CONTRAVARIANT_PIOLA_TRANSFORM
 63: } PetscDualSpaceTransformType;

 65: PETSC_EXTERN PetscClassId PETSCDUALSPACE_CLASSID;

 67: /*J
 68:   PetscDualSpaceType - String with the name of a PETSc dual space

 70:   Values:
 71: + PETSCDUALSPACELAGRANGE  - a dual space of pointwise evaluation functionals
 72: . PETSCDUALSPACESIMPLE -  a dual space defined by functionals provided with `PetscDualSpaceSimpleSetFunctional()`
 73: . PETSCDUALSPACEREFINED - the joint dual space defined by a group of cells, usually refined from one larger cell
 74: - PETSCDUALSPACEBDM - a dual space for Brezzi-Douglas-Marini elements

 76:   Level: beginner

 78: .seealso: `PetscDualSpaceSetType()`, `PetscDualSpace`, `PetscSpace`
 79: J*/
 80: typedef const char *PetscDualSpaceType;
 81: #define PETSCDUALSPACELAGRANGE "lagrange"
 82: #define PETSCDUALSPACESIMPLE   "simple"
 83: #define PETSCDUALSPACEREFINED  "refined"
 84: #define PETSCDUALSPACEBDM      "bdm"

 86: /*MC
 87:   PETSCDUALSPACEBDM = "bdm" - A `PetscDualSpace` object that encapsulates a dual space for Brezzi-Douglas-Marini elements

 89:   Level: intermediate

 91:   Note:
 92:   This type is a constructor alias of `PETSCDUALSPACELAGRANGE`.  During
 93:   `PetscDualSpaceSetUp()`, the correct value of `PetscDualSpaceSetFormDegree()` is
 94:   set for H-div conforming spaces. The type of the dual space is then changed to
 95:   to `PETSCDUALSPACELAGRANGE`.

 97: .seealso: `PetscDualSpaceType`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`, `PETSCDUALSPACELAGRANGE`, `PetscDualSpaceSetFormDegree()`
 98: M*/

100: PETSC_EXTERN PetscFunctionList PetscDualSpaceList;
101: PETSC_EXTERN PetscErrorCode    PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
102: PETSC_EXTERN PetscErrorCode    PetscDualSpaceDestroy(PetscDualSpace *);
103: PETSC_EXTERN PetscErrorCode    PetscDualSpaceDuplicate(PetscDualSpace, PetscDualSpace *);
104: PETSC_EXTERN PetscErrorCode    PetscDualSpaceSetType(PetscDualSpace, PetscDualSpaceType);
105: PETSC_EXTERN PetscErrorCode    PetscDualSpaceGetType(PetscDualSpace, PetscDualSpaceType *);
106: PETSC_EXTERN PetscErrorCode    PetscDualSpaceGetUniform(PetscDualSpace, PetscBool *);
107: PETSC_EXTERN PetscErrorCode    PetscDualSpaceGetNumDof(PetscDualSpace, const PetscInt **);
108: PETSC_EXTERN PetscErrorCode    PetscDualSpaceGetSection(PetscDualSpace, PetscSection *);
109: PETSC_EXTERN PetscErrorCode    PetscDualSpaceSetUp(PetscDualSpace);
110: PETSC_EXTERN PetscErrorCode    PetscDualSpaceSetFromOptions(PetscDualSpace);
111: PETSC_EXTERN PetscErrorCode    PetscDualSpaceViewFromOptions(PetscDualSpace, PetscObject, const char[]);

113: PETSC_EXTERN PetscErrorCode PetscDualSpaceView(PetscDualSpace, PetscViewer);
114: PETSC_EXTERN PetscErrorCode PetscDualSpaceRegister(const char[], PetscErrorCode (*)(PetscDualSpace));
115: PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterDestroy(void);

117: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace, PetscInt *);
118: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetInteriorDimension(PetscDualSpace, PetscInt *);
119: PETSC_EXTERN PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace, PetscInt);
120: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace, PetscInt *);
121: PETSC_EXTERN PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace, PetscInt);
122: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace, PetscInt *);
123: PETSC_EXTERN PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace, DM);
124: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace, DM *);
125: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace, PetscInt, PetscQuadrature *);
126: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace, const PetscInt ****, const PetscScalar ****);

128: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetAllData(PetscDualSpace, PetscQuadrature *, Mat *);
129: PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateAllDataDefault(PetscDualSpace, PetscQuadrature *, Mat *);
130: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetInteriorData(PetscDualSpace, PetscQuadrature *, Mat *);
131: PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateInteriorDataDefault(PetscDualSpace, PetscQuadrature *, Mat *);
132: PETSC_EXTERN PetscErrorCode PetscDualSpaceEqual(PetscDualSpace, PetscDualSpace, PetscBool *);

134: PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace, const PetscScalar *, PetscScalar *);
135: PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace, const PetscScalar *, PetscScalar *);
136: PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyInterior(PetscDualSpace, const PetscScalar *, PetscScalar *);
137: PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyInteriorDefault(PetscDualSpace, const PetscScalar *, PetscScalar *);

139: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFormDegree(PetscDualSpace, PetscInt *);
140: PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFormDegree(PetscDualSpace, PetscInt);
141: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace, PetscInt *);

143: PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace, PetscBool *);
144: PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace, PetscBool);
145: PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace, PetscBool *);
146: PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace, PetscBool);
147: PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetTrimmed(PetscDualSpace, PetscBool *);
148: PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetTrimmed(PetscDualSpace, PetscBool);
149: PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetNodeType(PetscDualSpace, PetscDTNodeType *, PetscBool *, PetscReal *);
150: PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetNodeType(PetscDualSpace, PetscDTNodeType, PetscBool, PetscReal);
151: PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetUseMoments(PetscDualSpace, PetscBool *);
152: PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetUseMoments(PetscDualSpace, PetscBool);
153: PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetMomentOrder(PetscDualSpace, PetscInt *);
154: PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetMomentOrder(PetscDualSpace, PetscInt);

156: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace, PetscInt, PetscDualSpace *);
157: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace, PetscInt, PetscDualSpace *);
158: PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace, PetscInt);
159: PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace, PetscInt, PetscQuadrature);

161: PETSC_EXTERN PetscErrorCode PetscDualSpaceRefinedSetCellSpaces(PetscDualSpace, const PetscDualSpace[]);
162: PETSC_EXTERN PetscErrorCode PetscSpaceCreateSubspace(PetscSpace, PetscDualSpace, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscCopyMode, PetscSpace *);

164: #endif