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