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