xref: /petsc/include/petscdualspace.h (revision 4ad8454beace47809662cdae21ee081016eaa39a)
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   {cite}`rogneskirbylogg2009`
53 
54 .seealso: `PetscDualSpaceGetDeRahm()`
55 M*/
56 typedef enum {
57   IDENTITY_TRANSFORM,
58   COVARIANT_PIOLA_TRANSFORM,
59   CONTRAVARIANT_PIOLA_TRANSFORM
60 } PetscDualSpaceTransformType;
61 
62 PETSC_EXTERN PetscClassId PETSCDUALSPACE_CLASSID;
63 
64 /*J
65   PetscDualSpaceType - String with the name of a PETSc dual space
66 
67   Values:
68 + PETSCDUALSPACELAGRANGE  - a dual space of pointwise evaluation functionals
69 . PETSCDUALSPACESIMPLE    - a dual space defined by functionals provided with `PetscDualSpaceSimpleSetFunctional()`
70 . PETSCDUALSPACEREFINED   - the joint dual space defined by a group of cells, usually refined from one larger cell
71 . PETSCDUALSPACEBDM       - a dual space for Brezzi-Douglas-Marini elements
72 - PETSCDUALSPACESUM       - a dual space that is a sum of other dual spaces
73 
74   Level: beginner
75 
76 .seealso: `PetscDualSpaceSetType()`, `PetscDualSpace`, `PetscSpace`
77 J*/
78 typedef const char *PetscDualSpaceType;
79 #define PETSCDUALSPACELAGRANGE "lagrange"
80 #define PETSCDUALSPACESIMPLE   "simple"
81 #define PETSCDUALSPACEREFINED  "refined"
82 #define PETSCDUALSPACEBDM      "bdm"
83 #define PETSCDUALSPACESUM      "sum"
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    PetscDualSpaceGetInteriorSection(PetscDualSpace, PetscSection *);
109 PETSC_EXTERN PetscErrorCode    PetscDualSpaceSetUp(PetscDualSpace);
110 PETSC_EXTERN PetscErrorCode    PetscDualSpaceSetFromOptions(PetscDualSpace);
111 PETSC_EXTERN PetscErrorCode    PetscDualSpaceViewFromOptions(PetscDualSpace, PetscObject, const char[]);
112 
113 PETSC_EXTERN PetscErrorCode PetscDualSpaceView(PetscDualSpace, PetscViewer);
114 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegister(const char[], PetscErrorCode (*)(PetscDualSpace));
115 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterDestroy(void);
116 
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 ****);
127 
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 *);
133 
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 *);
138 
139 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFormDegree(PetscDualSpace, PetscInt *);
140 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFormDegree(PetscDualSpace, PetscInt);
141 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace, PetscInt *);
142 
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);
155 
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);
160 
161 PETSC_EXTERN PetscErrorCode PetscDualSpaceRefinedSetCellSpaces(PetscDualSpace, const PetscDualSpace[]);
162 PETSC_EXTERN PetscErrorCode PetscSpaceCreateSubspace(PetscSpace, PetscDualSpace, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscCopyMode, PetscSpace *);
163 
164 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumSetNumSubspaces(PetscDualSpace, PetscInt);
165 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumGetNumSubspaces(PetscDualSpace, PetscInt *);
166 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumSetSubspace(PetscDualSpace, PetscInt, PetscDualSpace);
167 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumGetSubspace(PetscDualSpace, PetscInt, PetscDualSpace *);
168 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumSetConcatenate(PetscDualSpace, PetscBool);
169 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumGetConcatenate(PetscDualSpace, PetscBool *);
170 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumSetInterleave(PetscDualSpace, PetscBool, PetscBool);
171 PETSC_EXTERN PetscErrorCode PetscDualSpaceSumGetInterleave(PetscDualSpace, PetscBool *, PetscBool *);
172 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateSum(PetscInt, const PetscDualSpace[], PetscBool, PetscDualSpace *);
173