xref: /petsc/include/petscdualspace.h (revision 28d911a8077fa87b9a4a10e8dbd0dedeb22bcffd)
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   Values:
42 +  `IDENTITY_TRANSFORM` - make no changes in the function
43 .  `COVARIANT_PIOLA_TRANSFORM` - Covariant Piola: $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$
44 -  `CONTRAVARIANT_PIOLA_TRANSFORM` - Contravariant Piola: $\sigma^*(F) = 1/|J| J F \circ \phi^{-1)$
45 
46   Level: intermediate
47 
48   Note:
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 - PETSCDUALSPACESUM       - a dual space that is a sum of other dual spaces
75 
76   Level: beginner
77 
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"
85 #define PETSCDUALSPACESUM      "sum"
86 
87 /*MC
88   PETSCDUALSPACEBDM = "bdm" - A `PetscDualSpace` object that encapsulates a dual space for Brezzi-Douglas-Marini elements
89 
90   Level: intermediate
91 
92   Note:
93   This type is a constructor alias of `PETSCDUALSPACELAGRANGE`.  During
94   `PetscDualSpaceSetUp()`, the correct value of `PetscDualSpaceSetFormDegree()` is
95   set for H-div conforming spaces. The type of the dual space is then changed to
96   to `PETSCDUALSPACELAGRANGE`.
97 
98 .seealso: `PetscDualSpaceType`, `PetscDualSpaceCreate()`, `PetscDualSpaceSetType()`, `PETSCDUALSPACELAGRANGE`, `PetscDualSpaceSetFormDegree()`
99 M*/
100 
101 PETSC_EXTERN PetscFunctionList PetscDualSpaceList;
102 PETSC_EXTERN PetscErrorCode    PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
103 PETSC_EXTERN PetscErrorCode    PetscDualSpaceDestroy(PetscDualSpace *);
104 PETSC_EXTERN PetscErrorCode    PetscDualSpaceDuplicate(PetscDualSpace, PetscDualSpace *);
105 PETSC_EXTERN PetscErrorCode    PetscDualSpaceSetType(PetscDualSpace, PetscDualSpaceType);
106 PETSC_EXTERN PetscErrorCode    PetscDualSpaceGetType(PetscDualSpace, PetscDualSpaceType *);
107 PETSC_EXTERN PetscErrorCode    PetscDualSpaceGetUniform(PetscDualSpace, PetscBool *);
108 PETSC_EXTERN PetscErrorCode    PetscDualSpaceGetNumDof(PetscDualSpace, const PetscInt **);
109 PETSC_EXTERN PetscErrorCode    PetscDualSpaceGetSection(PetscDualSpace, PetscSection *);
110 PETSC_EXTERN PetscErrorCode    PetscDualSpaceGetInteriorSection(PetscDualSpace, PetscSection *);
111 PETSC_EXTERN PetscErrorCode    PetscDualSpaceSetUp(PetscDualSpace);
112 PETSC_EXTERN PetscErrorCode    PetscDualSpaceSetFromOptions(PetscDualSpace);
113 PETSC_EXTERN PetscErrorCode    PetscDualSpaceViewFromOptions(PetscDualSpace, PetscObject, const char[]);
114 
115 PETSC_EXTERN PetscErrorCode PetscDualSpaceView(PetscDualSpace, PetscViewer);
116 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegister(const char[], PetscErrorCode (*)(PetscDualSpace));
117 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterDestroy(void);
118 
119 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace, PetscInt *);
120 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetInteriorDimension(PetscDualSpace, PetscInt *);
121 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace, PetscInt);
122 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace, PetscInt *);
123 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace, PetscInt);
124 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace, PetscInt *);
125 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace, DM);
126 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace, DM *);
127 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace, PetscInt, PetscQuadrature *);
128 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace, const PetscInt ****, const PetscScalar ****);
129 
130 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetAllData(PetscDualSpace, PetscQuadrature *, Mat *);
131 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateAllDataDefault(PetscDualSpace, PetscQuadrature *, Mat *);
132 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetInteriorData(PetscDualSpace, PetscQuadrature *, Mat *);
133 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateInteriorDataDefault(PetscDualSpace, PetscQuadrature *, Mat *);
134 PETSC_EXTERN PetscErrorCode PetscDualSpaceEqual(PetscDualSpace, PetscDualSpace, PetscBool *);
135 
136 PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace, const PetscScalar *, PetscScalar *);
137 PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace, const PetscScalar *, PetscScalar *);
138 PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyInterior(PetscDualSpace, const PetscScalar *, PetscScalar *);
139 PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyInteriorDefault(PetscDualSpace, const PetscScalar *, PetscScalar *);
140 
141 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFormDegree(PetscDualSpace, PetscInt *);
142 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFormDegree(PetscDualSpace, PetscInt);
143 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace, PetscInt *);
144 
145 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace, PetscBool *);
146 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace, PetscBool);
147 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace, PetscBool *);
148 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace, PetscBool);
149 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetTrimmed(PetscDualSpace, PetscBool *);
150 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetTrimmed(PetscDualSpace, PetscBool);
151 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetNodeType(PetscDualSpace, PetscDTNodeType *, PetscBool *, PetscReal *);
152 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetNodeType(PetscDualSpace, PetscDTNodeType, PetscBool, PetscReal);
153 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetUseMoments(PetscDualSpace, PetscBool *);
154 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetUseMoments(PetscDualSpace, PetscBool);
155 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetMomentOrder(PetscDualSpace, PetscInt *);
156 PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetMomentOrder(PetscDualSpace, PetscInt);
157 
158 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace, PetscInt, PetscDualSpace *);
159 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace, PetscInt, PetscDualSpace *);
160 PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace, PetscInt);
161 PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace, PetscInt, PetscQuadrature);
162 
163 PETSC_EXTERN PetscErrorCode PetscDualSpaceRefinedSetCellSpaces(PetscDualSpace, const PetscDualSpace[]);
164 PETSC_EXTERN PetscErrorCode PetscSpaceCreateSubspace(PetscSpace, PetscDualSpace, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscCopyMode, PetscSpace *);
165 
166 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumSetNumSubspaces(PetscDualSpace, PetscInt);
167 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumGetNumSubspaces(PetscDualSpace, PetscInt *);
168 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumSetSubspace(PetscDualSpace, PetscInt, PetscDualSpace);
169 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumGetSubspace(PetscDualSpace, PetscInt, PetscDualSpace *);
170 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumSetConcatenate(PetscDualSpace, PetscBool);
171 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumGetConcatenate(PetscDualSpace, PetscBool *);
172 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumSetInterleave(PetscDualSpace, PetscBool, PetscBool);
173 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumGetInterleave(PetscDualSpace, PetscBool *, PetscBool *);
174 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateSum(PetscInt, const PetscDualSpace[], PetscBool, PetscDualSpace *);
175