xref: /petsc/include/petscfe.h (revision a496304597bacff3545e802853d69e8765312868)
1dbe77d9eSMatthew G. Knepley /*
2dbe77d9eSMatthew G. Knepley       Objects which encapsulate finite element spaces and operations
3dbe77d9eSMatthew G. Knepley */
4*a4963045SJacob Faibussowitsch #pragma once
5dbe77d9eSMatthew G. Knepley #include <petscdm.h>
60ddb9b0bSMatthew G. Knepley #include <petscdt.h>
7dbe77d9eSMatthew G. Knepley #include <petscfetypes.h>
82764a2aaSMatthew G. Knepley #include <petscdstypes.h>
9660d4ad9SBarry Smith #include <petscspace.h>
10660d4ad9SBarry Smith #include <petscdualspace.h>
11dbe77d9eSMatthew G. Knepley 
12ac09b921SBarry Smith /* SUBMANSEC = FE */
13ac09b921SBarry Smith 
14dce8aebaSBarry Smith /*MC
15dce8aebaSBarry Smith       PetscFEGeom - Structure for geometric information for `PetscFE`
16dce8aebaSBarry Smith 
17dce8aebaSBarry Smith     Level: intermediate
18dce8aebaSBarry Smith 
1916a05f60SBarry Smith .seealso: `PetscFE`, `PetscFEGeomCreate()`, `PetscFEGeomDestroy()`, `PetscFEGeomGetChunk()`, `PetscFEGeomRestoreChunk()`, `PetscFEGeomGetPoint()`, `PetscFEGeomGetCellPoint()`,
20b24fb147SBarry Smith           `PetscFEGeomComplete()`, `PetscSpace`, `PetscDualSpace`
21dce8aebaSBarry Smith M*/
224129dba9SToby Isaac typedef struct _n_PetscFEGeom {
234129dba9SToby Isaac   const PetscReal *xi;
2441418987SMatthew G. Knepley   PetscReal       *v;     /* v[Nc*Np*dE]:           The first point in each each in real coordinates */
2541418987SMatthew G. Knepley   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) */
2641418987SMatthew G. Knepley   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) */
2741418987SMatthew G. Knepley   PetscReal       *detJ;  /* detJ[Nc*Np]:           The determinant of J, and if it is non-square its the volume change */
28f15274beSMatthew Knepley   PetscReal       *n;     /* n[Nc*Np*dE]:           For faces, the normal to the face in real coordinates, outward for the first supporting cell */
290e18dc48SMatthew G. Knepley   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 */
3041418987SMatthew G. Knepley   PetscReal *suppJ[2];    /* sJ[s][Nc*Np*dE*dE]:    For faces, the Jacobian for each supporting cell s */
3141418987SMatthew G. Knepley   PetscReal *suppInvJ[2]; /* sInvJ[s][Nc*Np*dE*dE]: For faces, the inverse Jacobian for each supporting cell s */
3241418987SMatthew G. Knepley   PetscReal *suppDetJ[2]; /* sInvJ[s][Nc*Np*dE*dE]: For faces, the Jacobian determinant for each supporting cell s */
3362bd480fSMatthew G. Knepley   PetscInt   dim;         /* dim: Topological dimension */
3462bd480fSMatthew G. Knepley   PetscInt   dimEmbed;    /* dE:  coordinate dimension */
3562bd480fSMatthew G. Knepley   PetscInt   numCells;    /* Nc:  Number of mesh points represented in the arrays */
3662bd480fSMatthew G. Knepley   PetscInt   numPoints;   /* Np:  Number of evaluation points represented in the arrays */
3741418987SMatthew G. Knepley   PetscBool  isAffine;    /* Flag for affine transforms */
385fedec97SMatthew G. Knepley   PetscBool  isCohesive;  /* Flag for a cohesive cell */
394129dba9SToby Isaac } PetscFEGeom;
404129dba9SToby Isaac 
41dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEInitializePackage(void);
42dbe77d9eSMatthew G. Knepley 
434129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomCreate(PetscQuadrature, PetscInt, PetscInt, PetscBool, PetscFEGeom **);
444129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom **);
454129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom **);
466587ee25SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetPoint(PetscFEGeom *, PetscInt, PetscInt, const PetscReal[], PetscFEGeom *);
476587ee25SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetCellPoint(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom *);
484129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomComplete(PetscFEGeom *);
494129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **);
504129dba9SToby Isaac 
51c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, PetscScalar *);
52c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, PetscScalar *);
53ebac44aeSMatthew G. Knepley 
544bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransform(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
55f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
56f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformHessian(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
572edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePullback(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
582edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforward(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
592edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
60f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardHessian(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]);
614bee2e38SMatthew G. Knepley 
62dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCFE_CLASSID;
63dbe77d9eSMatthew G. Knepley 
640483ade4SMatthew G. Knepley /*J
650483ade4SMatthew G. Knepley   PetscFEType - String with the name of a PETSc finite element space
660483ade4SMatthew G. Knepley 
670483ade4SMatthew G. Knepley   Level: beginner
680483ade4SMatthew G. Knepley 
6916a05f60SBarry Smith   Note:
7016a05f60SBarry Smith   Currently, the classes are concerned with the implementation of element integration
710483ade4SMatthew G. Knepley 
72db781477SPatrick Sanan .seealso: `PetscFESetType()`, `PetscFE`
730483ade4SMatthew G. Knepley J*/
740483ade4SMatthew G. Knepley typedef const char *PetscFEType;
750483ade4SMatthew G. Knepley #define PETSCFEBASIC     "basic"
760483ade4SMatthew G. Knepley #define PETSCFEOPENCL    "opencl"
77aaf1837cSMatthew G. Knepley #define PETSCFECOMPOSITE "composite"
780483ade4SMatthew G. Knepley 
790483ade4SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscFEList;
80568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscFECreate(MPI_Comm, PetscFE *);
81568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscFEDestroy(PetscFE *);
820483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscFESetType(PetscFE, PetscFEType);
830483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscFEGetType(PetscFE, PetscFEType *);
840483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscFESetUp(PetscFE);
850483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscFESetFromOptions(PetscFE);
86fe2efc57SMark PETSC_EXTERN PetscErrorCode    PetscFEViewFromOptions(PetscFE, PetscObject, const char[]);
873f6b16c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode    PetscFESetName(PetscFE, const char[]);
888aec7d55SBarry Smith 
89fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscFEView(PetscFE, PetscViewer);
900483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegister(const char[], PetscErrorCode (*)(PetscFE));
910483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void);
92e8d98e54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateDefault(MPI_Comm, PetscInt, PetscInt, PetscBool, const char[], PetscInt, PetscFE *);
932df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, const char[], PetscInt, PetscFE *);
94e703855dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateLagrange(MPI_Comm, PetscInt, PetscInt, PetscBool, PetscInt, PetscInt, PetscFE *);
952df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateLagrangeByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscFE *);
967c48043bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateFromSpaces(PetscSpace, PetscDualSpace, PetscQuadrature, PetscQuadrature, PetscFE *);
979a1a3eb8SMatthew G. Knepley 
989a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *);
990ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *);
1009a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt);
1019a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *);
1029a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
1039a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt);
1049a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace);
1059a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *);
1069a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace);
1079a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *);
108bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature);
109bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *);
1101bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFESetFaceQuadrature(PetscFE, PetscQuadrature);
1111bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFEGetFaceQuadrature(PetscFE, PetscQuadrature *);
1125dc5c000SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECopyQuadrature(PetscFE, PetscFE);
1130ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **);
114ef0bb6c7SMatthew G. Knepley 
115ef0bb6c7SMatthew G. Knepley /* TODO: Need a function to reuse the memory when retabulating the same FE at different points */
116f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetCellTabulation(PetscFE, PetscInt, PetscTabulation *);
117f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceTabulation(PetscFE, PetscInt, PetscTabulation *);
118ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE, PetscTabulation *);
119ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateTabulation(PetscFE, PetscInt, PetscInt, const PetscReal[], PetscInt, PetscTabulation *);
120ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEComputeTabulation(PetscFE, PetscInt, const PetscReal[], PetscInt, PetscTabulation);
121ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscTabulationDestroy(PetscTabulation *);
122ef0bb6c7SMatthew G. Knepley 
123bceba477SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERefine(PetscFE, PetscFE *);
124228cfb18SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetHeightSubspace(PetscFE, PetscInt, PetscFE *);
125a0845e3aSMatthew G. Knepley 
126c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateCellGeometry(PetscFE, PetscQuadrature, PetscFEGeom *);
127c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroyCellGeometry(PetscFE, PetscFEGeom *);
1282edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforward(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]);
1292edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforwardGradient(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]);
130f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEPushforwardHessian(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]);
1314bee2e38SMatthew G. Knepley 
1324bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrate(PetscDS, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
1339371c9d4SSatish Balay 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[]);
13406ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscDS, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
13506ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscDS, PetscWeakForm, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
13607218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridResidual(PetscDS, PetscDS, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]);
13706ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
13806ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS, PetscWeakForm, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
13907218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridJacobian(PetscDS, PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]);
140a0845e3aSMatthew G. Knepley 
14189710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECompositeGetMapping(PetscFE, PetscInt *, const PetscReal *[], const PetscReal *[], const PetscReal *[]);
14289710940SMatthew G. Knepley 
1432f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreateHeightTrace(PetscFE, PetscInt, PetscFE *);
1442f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE, PetscInt, PetscFE *);
1452f5fb066SToby Isaac 
146855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType);
147855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *);
148