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