xref: /petsc/include/petscfe.h (revision 3a2a065f407e7da3bb6fc7eaf619dbab56034fde)
1 /*
2       Objects which encapsulate finite element spaces and operations
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 #include <petscdualspace.h>
11 
12 /* SUBMANSEC = FE */
13 
14 /*MC
15     PetscFEGeom - Structure for geometric information for `PetscFE`
16 
17     Level: intermediate
18 
19     Note:
20     This is a struct, not a `PetscObject`
21 
22 .seealso: `PetscFE`, `PetscFEGeomCreate()`, `PetscFEGeomDestroy()`, `PetscFEGeomGetChunk()`, `PetscFEGeomRestoreChunk()`, `PetscFEGeomGetPoint()`, `PetscFEGeomGetCellPoint()`,
23           `PetscFEGeomComplete()`, `PetscSpace`, `PetscDualSpace`
24 M*/
25 typedef struct _n_PetscFEGeom {
26   const PetscReal *xi;
27   PetscReal       *v;     /* v[Nc*Np*dE]:           The first point in each each in real coordinates */
28   PetscReal       *J;     /* J[Nc*Np*dE*dE]:        The Jacobian of the map from reference to real coordinates (if nonsquare it is completed with orthogonal columns) */
29   PetscReal       *invJ;  /* invJ[Nc*Np*dE*dE]:     The inverse of the Jacobian of the map from reference to real coordinates (if nonsquare it is completed with orthogonal columns) */
30   PetscReal       *detJ;  /* detJ[Nc*Np]:           The determinant of J, and if it is non-square its the volume change */
31   PetscReal       *n;     /* n[Nc*Np*dE]:           For faces, the normal to the face in real coordinates, outward for the first supporting cell */
32   PetscInt (*face)[4];    /* face[Nc][s*2]:         For faces, the local face number (cone index) and orientation for this face in each supporting cell s */
33   PetscReal *suppJ[2];    /* sJ[s][Nc*Np*dE*dE]:    For faces, the Jacobian for each supporting cell s */
34   PetscReal *suppInvJ[2]; /* sInvJ[s][Nc*Np*dE*dE]: For faces, the inverse Jacobian for each supporting cell s */
35   PetscReal *suppDetJ[2]; /* sInvJ[s][Nc*Np*dE*dE]: For faces, the Jacobian determinant for each supporting cell s */
36   PetscInt   dim;         /* dim: Topological dimension */
37   PetscInt   dimEmbed;    /* dE:  coordinate dimension */
38   PetscInt   numCells;    /* Nc:  Number of mesh points represented in the arrays */
39   PetscInt   numPoints;   /* Np:  Number of evaluation points represented in the arrays */
40   PetscBool  isAffine;    /* Flag for affine transforms */
41   PetscBool  isCohesive;  /* Flag for a cohesive cell */
42 } PetscFEGeom;
43 
44 PETSC_EXTERN PetscErrorCode PetscFEInitializePackage(void);
45 
46 PETSC_EXTERN PetscErrorCode PetscFEGeomCreate(PetscQuadrature, PetscInt, PetscInt, PetscBool, PetscFEGeom **);
47 PETSC_EXTERN PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom **);
48 PETSC_EXTERN PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom **);
49 PETSC_EXTERN PetscErrorCode PetscFEGeomGetPoint(PetscFEGeom *, PetscInt, PetscInt, const PetscReal[], PetscFEGeom *);
50 PETSC_EXTERN PetscErrorCode PetscFEGeomGetCellPoint(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom *);
51 PETSC_EXTERN PetscErrorCode PetscFEGeomComplete(PetscFEGeom *);
52 PETSC_EXTERN PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **);
53 
54 PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, PetscScalar *);
55 PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, PetscScalar *);
56 
57 PETSC_EXTERN PetscErrorCode PetscDualSpaceTransform(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
58 PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
59 PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformHessian(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
60 PETSC_EXTERN PetscErrorCode PetscDualSpacePullback(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
61 PETSC_EXTERN PetscErrorCode PetscDualSpacePushforward(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
62 PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
63 PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardHessian(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
64 
65 PETSC_EXTERN PetscClassId PETSCFE_CLASSID;
66 
67 /*J
68   PetscFEType - String with the name of a PETSc finite element space
69 
70   Level: beginner
71 
72   Note:
73   Currently, the classes are concerned with the implementation of element integration
74 
75 .seealso: `PetscFESetType()`, `PetscFE`
76 J*/
77 typedef const char *PetscFEType;
78 #define PETSCFEBASIC     "basic"
79 #define PETSCFEOPENCL    "opencl"
80 #define PETSCFECOMPOSITE "composite"
81 #define PETSCFEVECTOR    "vector"
82 
83 PETSC_EXTERN PetscFunctionList PetscFEList;
84 PETSC_EXTERN PetscErrorCode    PetscFECreate(MPI_Comm, PetscFE *);
85 PETSC_EXTERN PetscErrorCode    PetscFEDestroy(PetscFE *);
86 PETSC_EXTERN PetscErrorCode    PetscFESetType(PetscFE, PetscFEType);
87 PETSC_EXTERN PetscErrorCode    PetscFEGetType(PetscFE, PetscFEType *);
88 PETSC_EXTERN PetscErrorCode    PetscFESetUp(PetscFE);
89 PETSC_EXTERN PetscErrorCode    PetscFESetFromOptions(PetscFE);
90 PETSC_EXTERN PetscErrorCode    PetscFEViewFromOptions(PetscFE, PetscObject, const char[]);
91 PETSC_EXTERN PetscErrorCode    PetscFESetName(PetscFE, const char[]);
92 PETSC_EXTERN PetscErrorCode    PetscFECreateVector(PetscFE, PetscInt, PetscBool, PetscBool, PetscFE *);
93 
94 PETSC_EXTERN PetscErrorCode PetscFEView(PetscFE, PetscViewer);
95 PETSC_EXTERN PetscErrorCode PetscFERegister(const char[], PetscErrorCode (*)(PetscFE));
96 PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void);
97 PETSC_EXTERN PetscErrorCode PetscFECreateDefault(MPI_Comm, PetscInt, PetscInt, PetscBool, const char[], PetscInt, PetscFE *);
98 PETSC_EXTERN PetscErrorCode PetscFECreateByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, const char[], PetscInt, PetscFE *);
99 PETSC_EXTERN PetscErrorCode PetscFECreateLagrange(MPI_Comm, PetscInt, PetscInt, PetscBool, PetscInt, PetscInt, PetscFE *);
100 PETSC_EXTERN PetscErrorCode PetscFECreateLagrangeByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscFE *);
101 PETSC_EXTERN PetscErrorCode PetscFECreateFromSpaces(PetscSpace, PetscDualSpace, PetscQuadrature, PetscQuadrature, PetscFE *);
102 
103 PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *);
104 PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *);
105 PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt);
106 PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *);
107 PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
108 PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt);
109 PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace);
110 PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *);
111 PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace);
112 PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *);
113 PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature);
114 PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *);
115 PETSC_EXTERN PetscErrorCode PetscFESetFaceQuadrature(PetscFE, PetscQuadrature);
116 PETSC_EXTERN PetscErrorCode PetscFEGetFaceQuadrature(PetscFE, PetscQuadrature *);
117 PETSC_EXTERN PetscErrorCode PetscFECopyQuadrature(PetscFE, PetscFE);
118 PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **);
119 
120 /* TODO: Need a function to reuse the memory when retabulating the same FE at different points */
121 PETSC_EXTERN PetscErrorCode PetscFEGetCellTabulation(PetscFE, PetscInt, PetscTabulation *);
122 PETSC_EXTERN PetscErrorCode PetscFEGetFaceTabulation(PetscFE, PetscInt, PetscTabulation *);
123 PETSC_EXTERN PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE, PetscTabulation *);
124 PETSC_EXTERN PetscErrorCode PetscFECreateTabulation(PetscFE, PetscInt, PetscInt, const PetscReal[], PetscInt, PetscTabulation *);
125 PETSC_EXTERN PetscErrorCode PetscFEComputeTabulation(PetscFE, PetscInt, const PetscReal[], PetscInt, PetscTabulation);
126 PETSC_EXTERN PetscErrorCode PetscTabulationDestroy(PetscTabulation *);
127 
128 PETSC_EXTERN PetscErrorCode PetscFERefine(PetscFE, PetscFE *);
129 PETSC_EXTERN PetscErrorCode PetscFEGetHeightSubspace(PetscFE, PetscInt, PetscFE *);
130 
131 PETSC_EXTERN PetscErrorCode PetscFECreateCellGeometry(PetscFE, PetscQuadrature, PetscFEGeom *);
132 PETSC_EXTERN PetscErrorCode PetscFEDestroyCellGeometry(PetscFE, PetscFEGeom *);
133 PETSC_EXTERN PetscErrorCode PetscFEPushforward(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]);
134 PETSC_EXTERN PetscErrorCode PetscFEPushforwardGradient(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]);
135 PETSC_EXTERN PetscErrorCode PetscFEPushforwardHessian(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]);
136 
137 PETSC_EXTERN PetscErrorCode PetscFEIntegrate(PetscDS, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
138 PETSC_EXTERN PetscErrorCode PetscFEIntegrateBd(PetscDS, PetscInt, void (*)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
139 PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscDS, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
140 PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscDS, PetscWeakForm, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
141 PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridResidual(PetscDS, PetscDS, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
142 PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
143 PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS, PetscWeakForm, PetscFEJacobianType, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
144 PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridJacobian(PetscDS, PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
145 
146 PETSC_EXTERN PetscErrorCode PetscFECompositeGetMapping(PetscFE, PetscInt *, const PetscReal *[], const PetscReal *[], const PetscReal *[]);
147 
148 PETSC_EXTERN PetscErrorCode PetscFECreateHeightTrace(PetscFE, PetscInt, PetscFE *);
149 PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE, PetscInt, PetscFE *);
150 
151 PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType);
152 PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *);
153 
154 #ifdef PETSC_HAVE_LIBCEED
155 
156   #ifndef PLEXFE_QFUNCTION
157     #define PLEXFE_QFUNCTION(fname, f0_name, f1_name) \
158       CEED_QFUNCTION(PlexQFunction##fname)(void *ctx, const CeedInt Q, const CeedScalar *const *in, CeedScalar *const *out) \
159       { \
160         const CeedScalar *u = in[0], *du = in[1], *qdata = in[2]; \
161         CeedScalar       *v = out[0], *dv = out[1]; \
162         const PetscInt    Nc   = 1; \
163         const PetscInt    cdim = 2; \
164 \
165         CeedPragmaSIMD for (CeedInt i = 0; i < Q; ++i) \
166         { \
167           const PetscInt   uOff[2]    = {0, Nc}; \
168           const PetscInt   uOff_x[2]  = {0, Nc * cdim}; \
169           const CeedScalar x[2]       = {qdata[i + Q * 1], qdata[i + Q * 2]}; \
170           const CeedScalar invJ[2][2] = { \
171             {qdata[i + Q * 3], qdata[i + Q * 5]}, \
172             {qdata[i + Q * 4], qdata[i + Q * 6]} \
173           }; \
174           const CeedScalar u_x[2] = {invJ[0][0] * du[i + Q * 0] + invJ[1][0] * du[i + Q * 1], invJ[0][1] * du[i + Q * 0] + invJ[1][1] * du[i + Q * 1]}; \
175           PetscScalar      f0[Nc]; \
176           PetscScalar      f1[Nc * cdim]; \
177 \
178           for (PetscInt k = 0; k < Nc; ++k) f0[k] = 0; \
179           for (PetscInt k = 0; k < Nc * cdim; ++k) f1[k] = 0; \
180           f0_name(2, 1, 0, uOff, uOff_x, u, NULL, u_x, NULL, NULL, NULL, NULL, NULL, 0.0, x, 0, NULL, f0); \
181           f1_name(2, 1, 0, uOff, uOff_x, u, NULL, u_x, NULL, NULL, NULL, NULL, NULL, 0.0, x, 0, NULL, f1); \
182 \
183           dv[i + Q * 0] = qdata[i + Q * 0] * (invJ[0][0] * f1[0] + invJ[0][1] * f1[1]); \
184           dv[i + Q * 1] = qdata[i + Q * 0] * (invJ[1][0] * f1[0] + invJ[1][1] * f1[1]); \
185           v[i]          = qdata[i + Q * 0] * f0[0]; \
186         } \
187         return CEED_ERROR_SUCCESS; \
188       }
189   #endif
190 
191 #else
192 
193   #ifndef PLEXFE_QFUNCTION
194     #define PLEXFE_QFUNCTION(fname, f0_name, f1_name)
195   #endif
196 
197 #endif
198