Actual source code: petscdualspace.h

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

 11: /* MANSEC = DM */
 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:   Values:
 43: +  `IDENTITY_TRANSFORM`            - make no changes in the function
 44: .  `COVARIANT_PIOLA_TRANSFORM`     - Covariant Piola: $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$
 45: -  `CONTRAVARIANT_PIOLA_TRANSFORM` - Contravariant Piola: $\sigma^*(F) = 1/|J| J F \circ \phi^{-1)$

 47:   Level: intermediate

 49:   Note:
 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.
 53:   {cite}`rogneskirbylogg2009`

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

 63: PETSC_EXTERN PetscClassId PETSCDUALSPACE_CLASSID;

 65: /*J
 66:   PetscDualSpaceType - String with the name of a PETSc dual space

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

 75:   Level: beginner

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

 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    PetscDualSpaceGetInteriorSection(PetscDualSpace, PetscSection *);
110: PETSC_EXTERN PetscErrorCode    PetscDualSpaceSetUp(PetscDualSpace);
111: PETSC_EXTERN PetscErrorCode    PetscDualSpaceSetFromOptions(PetscDualSpace);
112: PETSC_EXTERN PetscErrorCode    PetscDualSpaceViewFromOptions(PetscDualSpace, PetscObject, const char[]);

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

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

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

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

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

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

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

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

165: PETSC_EXTERN PetscErrorCode PetscDualSpaceSumSetNumSubspaces(PetscDualSpace, PetscInt);
166: PETSC_EXTERN PetscErrorCode PetscDualSpaceSumGetNumSubspaces(PetscDualSpace, PetscInt *);
167: PETSC_EXTERN PetscErrorCode PetscDualSpaceSumSetSubspace(PetscDualSpace, PetscInt, PetscDualSpace);
168: PETSC_EXTERN PetscErrorCode PetscDualSpaceSumGetSubspace(PetscDualSpace, PetscInt, PetscDualSpace *);
169: PETSC_EXTERN PetscErrorCode PetscDualSpaceSumSetConcatenate(PetscDualSpace, PetscBool);
170: PETSC_EXTERN PetscErrorCode PetscDualSpaceSumGetConcatenate(PetscDualSpace, PetscBool *);
171: PETSC_EXTERN PetscErrorCode PetscDualSpaceSumSetInterleave(PetscDualSpace, PetscBool, PetscBool);
172: PETSC_EXTERN PetscErrorCode PetscDualSpaceSumGetInterleave(PetscDualSpace, PetscBool *, PetscBool *);
173: PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateSum(PetscInt, const PetscDualSpace[], PetscBool, PetscDualSpace *);