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> 10 11 /* SUBMANSEC = DUALSPACE */ 12 13 /*S 14 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 15 16 Level: beginner 17 18 .seealso: `PetscDualSpaceCreate()`, `PetscSpace`, `PetscSpaceCreate()`, `PetscDualSpaceSetType()`, `PetscDualSpaceType` 19 S*/ 20 typedef struct _p_PetscDualSpace *PetscDualSpace; 21 22 /*MC 23 PetscDualSpaceReferenceCell - The type of reference cell 24 25 Level: beginner 26 27 Note: 28 This is used only for automatic creation of reference cells. A `PetscDualSpace` can accept an arbitrary `DM` for a reference cell. 29 30 .seealso: `PetscSpace`, `PetscDualSpaceCreate()`, `PetscDualSpaceType` 31 M*/ 32 typedef enum { 33 PETSCDUALSPACE_REFCELL_SIMPLEX, 34 PETSCDUALSPACE_REFCELL_TENSOR 35 } PetscDualSpaceReferenceCell; 36 PETSC_EXTERN const char *const PetscDualSpaceReferenceCells[]; 37 38 /*MC 39 PetscDualSpaceTransformType - The type of function transform 40 41 Level: intermediate 42 43 Values: 44 + `IDENTITY_TRANSFORM` - make no changes in the function 45 . `COVARIANT_PIOLA_TRANSFORM` - Covariant Piola: $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ 46 - `CONTRAVARIANT_PIOLA_TRANSFORM` - Contravariant Piola: $\sigma^*(F) = 1/|J| J F \circ \phi^{-1)$ 47 48 Notes: 49 These transforms, and their inverses, are used to move functions and functionals between the reference element and real space. 50 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, 51 so that it acts on real space, we use the pushforward transform $\sigma^*$. The pullback $\sigma_*$ is the inverse transform. 52 53 References: 54 . Rognes, Kirby, and Logg, Efficient Assembly of Hdiv and Hrot Conforming Finite Elements, SISC, 31(6), 4130-4151, arXiv 1205.3085, 2010 55 56 .seealso: `PetscDualSpaceGetDeRahm()` 57 M*/ 58 typedef enum { 59 IDENTITY_TRANSFORM, 60 COVARIANT_PIOLA_TRANSFORM, 61 CONTRAVARIANT_PIOLA_TRANSFORM 62 } PetscDualSpaceTransformType; 63 64 PETSC_EXTERN PetscClassId PETSCDUALSPACE_CLASSID; 65 66 /*J 67 PetscDualSpaceType - String with the name of a PETSc dual space 68 69 Values: 70 + PETSCDUALSPACELAGRANGE - a dual space of pointwise evaluation functionals 71 . PETSCDUALSPACESIMPLE - a dual space defined by functionals provided with `PetscDualSpaceSimpleSetFunctional()` 72 . PETSCDUALSPACEREFINED - the joint dual space defined by a group of cells, usually refined from one larger cell 73 - PETSCDUALSPACEBDM - a dual space for Brezzi-Douglas-Marini elements 74 75 Level: beginner 76 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 85 /*MC 86 PETSCDUALSPACEBDM = "bdm" - A `PetscDualSpace` object that encapsulates a dual space for Brezzi-Douglas-Marini elements 87 88 Level: intermediate 89 90 Note: 91 This type is a constructor alias of `PETSCDUALSPACELAGRANGE`. During 92 `PetscDualSpaceSetUp()`, the correct value of `PetscDualSpaceSetFormDegree()` is 93 set for H-div conforming spaces. The type of the dual space is then changed to 94 to `PETSCDUALSPACELAGRANGE`. 95 96 .seealso: `PetscDualSpaceType`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`, `PETSCDUALSPACELAGRANGE`, `PetscDualSpaceSetFormDegree()` 97 M*/ 98 99 PETSC_EXTERN PetscFunctionList PetscDualSpaceList; 100 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *); 101 PETSC_EXTERN PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *); 102 PETSC_EXTERN PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace, PetscDualSpace *); 103 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetType(PetscDualSpace, PetscDualSpaceType); 104 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetType(PetscDualSpace, PetscDualSpaceType *); 105 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetUniform(PetscDualSpace, PetscBool *); 106 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace, const PetscInt **); 107 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetSection(PetscDualSpace, PetscSection *); 108 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace); 109 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace); 110 PETSC_EXTERN PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace, PetscObject, const char[]); 111 112 PETSC_EXTERN PetscErrorCode PetscDualSpaceView(PetscDualSpace, PetscViewer); 113 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegister(const char[], PetscErrorCode (*)(PetscDualSpace)); 114 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterDestroy(void); 115 116 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace, PetscInt *); 117 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetInteriorDimension(PetscDualSpace, PetscInt *); 118 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace, PetscInt); 119 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace, PetscInt *); 120 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace, PetscInt); 121 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace, PetscInt *); 122 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace, DM); 123 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace, DM *); 124 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace, PetscInt, PetscQuadrature *); 125 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace, const PetscInt ****, const PetscScalar ****); 126 127 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetAllData(PetscDualSpace, PetscQuadrature *, Mat *); 128 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateAllDataDefault(PetscDualSpace, PetscQuadrature *, Mat *); 129 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetInteriorData(PetscDualSpace, PetscQuadrature *, Mat *); 130 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateInteriorDataDefault(PetscDualSpace, PetscQuadrature *, Mat *); 131 PETSC_EXTERN PetscErrorCode PetscDualSpaceEqual(PetscDualSpace, PetscDualSpace, PetscBool *); 132 133 PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace, const PetscScalar *, PetscScalar *); 134 PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace, const PetscScalar *, PetscScalar *); 135 PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyInterior(PetscDualSpace, const PetscScalar *, PetscScalar *); 136 PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyInteriorDefault(PetscDualSpace, const PetscScalar *, PetscScalar *); 137 138 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFormDegree(PetscDualSpace, PetscInt *); 139 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFormDegree(PetscDualSpace, PetscInt); 140 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace, PetscInt *); 141 142 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace, PetscBool *); 143 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace, PetscBool); 144 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace, PetscBool *); 145 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace, PetscBool); 146 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetTrimmed(PetscDualSpace, PetscBool *); 147 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetTrimmed(PetscDualSpace, PetscBool); 148 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetNodeType(PetscDualSpace, PetscDTNodeType *, PetscBool *, PetscReal *); 149 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetNodeType(PetscDualSpace, PetscDTNodeType, PetscBool, PetscReal); 150 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetUseMoments(PetscDualSpace, PetscBool *); 151 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetUseMoments(PetscDualSpace, PetscBool); 152 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetMomentOrder(PetscDualSpace, PetscInt *); 153 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetMomentOrder(PetscDualSpace, PetscInt); 154 155 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace, PetscInt, PetscDualSpace *); 156 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace, PetscInt, PetscDualSpace *); 157 PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace, PetscInt); 158 PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace, PetscInt, PetscQuadrature); 159 160 PETSC_EXTERN PetscErrorCode PetscDualSpaceRefinedSetCellSpaces(PetscDualSpace, const PetscDualSpace[]); 161 PETSC_EXTERN PetscErrorCode PetscSpaceCreateSubspace(PetscSpace, PetscDualSpace, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscCopyMode, PetscSpace *); 162