1dbe77d9eSMatthew G. Knepley /* 2dbe77d9eSMatthew G. Knepley Objects which encapsulate finite element spaces and operations 3dbe77d9eSMatthew G. Knepley */ 46524c165SJacob Faibussowitsch #ifndef PETSCFE_H 526bd1501SBarry Smith #define PETSCFE_H 6dbe77d9eSMatthew G. Knepley #include <petscdm.h> 70ddb9b0bSMatthew G. Knepley #include <petscdt.h> 8dbe77d9eSMatthew G. Knepley #include <petscfetypes.h> 92764a2aaSMatthew G. Knepley #include <petscdstypes.h> 10*660d4ad9SBarry Smith #include <petscspace.h> 11*660d4ad9SBarry Smith #include <petscdualspace.h> 12dbe77d9eSMatthew G. Knepley 13ac09b921SBarry Smith /* SUBMANSEC = FE */ 14ac09b921SBarry Smith 15dce8aebaSBarry Smith /*MC 16dce8aebaSBarry Smith PetscFEGeom - Structure for geometric information for `PetscFE` 17dce8aebaSBarry Smith 18dce8aebaSBarry Smith Level: intermediate 19dce8aebaSBarry Smith 2016a05f60SBarry Smith .seealso: `PetscFE`, `PetscFEGeomCreate()`, `PetscFEGeomDestroy()`, `PetscFEGeomGetChunk()`, `PetscFEGeomRestoreChunk()`, `PetscFEGeomGetPoint()`, `PetscFEGeomGetCellPoint()`, 21dce8aebaSBarry Smith `PetscFEGeomComplete()` 22dce8aebaSBarry Smith M*/ 234129dba9SToby Isaac typedef struct _n_PetscFEGeom { 244129dba9SToby Isaac const PetscReal *xi; 2541418987SMatthew G. Knepley PetscReal *v; /* v[Nc*Np*dE]: The first point in each each in real coordinates */ 2641418987SMatthew 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) */ 2741418987SMatthew 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) */ 2841418987SMatthew G. Knepley PetscReal *detJ; /* detJ[Nc*Np]: The determinant of J, and if it is non-square its the volume change */ 29f15274beSMatthew Knepley PetscReal *n; /* n[Nc*Np*dE]: For faces, the normal to the face in real coordinates, outward for the first supporting cell */ 300e18dc48SMatthew 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 */ 3141418987SMatthew G. Knepley PetscReal *suppJ[2]; /* sJ[s][Nc*Np*dE*dE]: For faces, the Jacobian for each supporting cell s */ 3241418987SMatthew G. Knepley PetscReal *suppInvJ[2]; /* sInvJ[s][Nc*Np*dE*dE]: For faces, the inverse Jacobian for each supporting cell s */ 3341418987SMatthew G. Knepley PetscReal *suppDetJ[2]; /* sInvJ[s][Nc*Np*dE*dE]: For faces, the Jacobian determinant for each supporting cell s */ 3462bd480fSMatthew G. Knepley PetscInt dim; /* dim: Topological dimension */ 3562bd480fSMatthew G. Knepley PetscInt dimEmbed; /* dE: coordinate dimension */ 3662bd480fSMatthew G. Knepley PetscInt numCells; /* Nc: Number of mesh points represented in the arrays */ 3762bd480fSMatthew G. Knepley PetscInt numPoints; /* Np: Number of evaluation points represented in the arrays */ 3841418987SMatthew G. Knepley PetscBool isAffine; /* Flag for affine transforms */ 395fedec97SMatthew G. Knepley PetscBool isCohesive; /* Flag for a cohesive cell */ 404129dba9SToby Isaac } PetscFEGeom; 414129dba9SToby Isaac 42dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEInitializePackage(void); 43dbe77d9eSMatthew G. Knepley 44dbe77d9eSMatthew G. Knepley 454129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomCreate(PetscQuadrature, PetscInt, PetscInt, PetscBool, PetscFEGeom **); 4662bd480fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetQuadrature(PetscFEGeom *, PetscQuadrature *); 4762bd480fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomSetQuadrature(PetscFEGeom *, PetscQuadrature); 484129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomGetChunk(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom **); 494129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomRestoreChunk(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom **); 506587ee25SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetPoint(PetscFEGeom *, PetscInt, PetscInt, const PetscReal[], PetscFEGeom *); 516587ee25SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGeomGetCellPoint(PetscFEGeom *, PetscInt, PetscInt, PetscFEGeom *); 524129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomComplete(PetscFEGeom *); 534129dba9SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEGeomDestroy(PetscFEGeom **); 544129dba9SToby Isaac 55c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, PetscScalar *); 56c330f8ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace, PetscInt, PetscReal, PetscFEGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar *, void *), void *, PetscScalar *); 57ebac44aeSMatthew G. Knepley 584bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransform(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 59f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 60f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceTransformHessian(PetscDualSpace, PetscDualSpaceTransformType, PetscBool, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 612edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePullback(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 622edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforward(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 632edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 64f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpacePushforwardHessian(PetscDualSpace, PetscFEGeom *, PetscInt, PetscInt, PetscScalar[]); 654bee2e38SMatthew G. Knepley 66dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCFE_CLASSID; 67dbe77d9eSMatthew G. Knepley 680483ade4SMatthew G. Knepley /*J 690483ade4SMatthew G. Knepley PetscFEType - String with the name of a PETSc finite element space 700483ade4SMatthew G. Knepley 710483ade4SMatthew G. Knepley Level: beginner 720483ade4SMatthew G. Knepley 7316a05f60SBarry Smith Note: 7416a05f60SBarry Smith Currently, the classes are concerned with the implementation of element integration 750483ade4SMatthew G. Knepley 76db781477SPatrick Sanan .seealso: `PetscFESetType()`, `PetscFE` 770483ade4SMatthew G. Knepley J*/ 780483ade4SMatthew G. Knepley typedef const char *PetscFEType; 790483ade4SMatthew G. Knepley #define PETSCFEBASIC "basic" 800483ade4SMatthew G. Knepley #define PETSCFEOPENCL "opencl" 81aaf1837cSMatthew G. Knepley #define PETSCFECOMPOSITE "composite" 820483ade4SMatthew G. Knepley 830483ade4SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscFEList; 84568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreate(MPI_Comm, PetscFE *); 85568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroy(PetscFE *); 860483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetType(PetscFE, PetscFEType); 870483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetType(PetscFE, PetscFEType *); 880483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetUp(PetscFE); 890483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetFromOptions(PetscFE); 90fe2efc57SMark PETSC_EXTERN PetscErrorCode PetscFEViewFromOptions(PetscFE, PetscObject, const char[]); 913f6b16c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetName(PetscFE, const char[]); 928aec7d55SBarry Smith 93fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscFEView(PetscFE, PetscViewer); 940483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegister(const char[], PetscErrorCode (*)(PetscFE)); 950483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void); 96e8d98e54SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateDefault(MPI_Comm, PetscInt, PetscInt, PetscBool, const char[], PetscInt, PetscFE *); 972df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, const char[], PetscInt, PetscFE *); 98e703855dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateLagrange(MPI_Comm, PetscInt, PetscInt, PetscBool, PetscInt, PetscInt, PetscFE *); 992df84da0SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateLagrangeByCell(MPI_Comm, PetscInt, PetscInt, DMPolytopeType, PetscInt, PetscInt, PetscFE *); 1007c48043bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateFromSpaces(PetscSpace, PetscDualSpace, PetscQuadrature, PetscQuadrature, PetscFE *); 1019a1a3eb8SMatthew G. Knepley 1029a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *); 1030ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *); 1049a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt); 1059a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *); 1069a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 1079a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt); 1089a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace); 1099a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *); 1109a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace); 1119a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *); 112bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature); 113bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *); 1141bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFESetFaceQuadrature(PetscFE, PetscQuadrature); 1151bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFEGetFaceQuadrature(PetscFE, PetscQuadrature *); 1165dc5c000SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECopyQuadrature(PetscFE, PetscFE); 1170ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **); 118ef0bb6c7SMatthew G. Knepley 119ef0bb6c7SMatthew G. Knepley /* TODO: Need a function to reuse the memory when retabulating the same FE at different points */ 120f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetCellTabulation(PetscFE, PetscInt, PetscTabulation *); 121f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceTabulation(PetscFE, PetscInt, PetscTabulation *); 122ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE, PetscTabulation *); 123ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateTabulation(PetscFE, PetscInt, PetscInt, const PetscReal[], PetscInt, PetscTabulation *); 124ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEComputeTabulation(PetscFE, PetscInt, const PetscReal[], PetscInt, PetscTabulation); 125ef0bb6c7SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscTabulationDestroy(PetscTabulation *); 126ef0bb6c7SMatthew G. Knepley 127bceba477SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERefine(PetscFE, PetscFE *); 128228cfb18SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetHeightSubspace(PetscFE, PetscInt, PetscFE *); 129a0845e3aSMatthew G. Knepley 130c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateCellGeometry(PetscFE, PetscQuadrature, PetscFEGeom *); 131c9ba7969SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroyCellGeometry(PetscFE, PetscFEGeom *); 1322edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforward(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]); 1332edcad52SToby Isaac PETSC_EXTERN PetscErrorCode PetscFEPushforwardGradient(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]); 134f9244615SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEPushforwardHessian(PetscFE, PetscFEGeom *, PetscInt, PetscScalar[]); 1354bee2e38SMatthew G. Knepley 1364bee2e38SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrate(PetscDS, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]); 1379371c9d4SSatish 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[]); 13806ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscDS, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]); 13906ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscDS, PetscWeakForm, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]); 14007218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridResidual(PetscDS, PetscDS, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]); 14106ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]); 14206ad1575SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdJacobian(PetscDS, PetscWeakForm, PetscFormKey, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]); 14307218a29SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateHybridJacobian(PetscDS, PetscDS, PetscFEJacobianType, PetscFormKey, PetscInt, PetscInt, PetscFEGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]); 144a0845e3aSMatthew G. Knepley 14589710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECompositeGetMapping(PetscFE, PetscInt *, const PetscReal *[], const PetscReal *[], const PetscReal *[]); 14689710940SMatthew G. Knepley 1472f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreateHeightTrace(PetscFE, PetscInt, PetscFE *); 1482f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE, PetscInt, PetscFE *); 1492f5fb066SToby Isaac 150855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType); 151855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *); 152855cd083SMatthew G. Knepley 153dbe77d9eSMatthew G. Knepley #endif 154