1dbe77d9eSMatthew G. Knepley /* 2dbe77d9eSMatthew G. Knepley Objects which encapsulate finite element spaces and operations 3dbe77d9eSMatthew G. Knepley */ 4dbe77d9eSMatthew G. Knepley #if !defined(__PETSCFE_H) 5dbe77d9eSMatthew G. Knepley #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> 10dbe77d9eSMatthew G. Knepley 11c0d900a5SMatthew G. Knepley /* Assuming dim <= 3 */ 12c0d900a5SMatthew G. Knepley typedef struct { 13c0d900a5SMatthew G. Knepley PetscReal v0[3]; 14c0d900a5SMatthew G. Knepley PetscReal J[9]; 15c0d900a5SMatthew G. Knepley PetscReal invJ[9]; 16c0d900a5SMatthew G. Knepley PetscReal detJ; 177a1a1bd4SToby Isaac PetscInt dim; 187a1a1bd4SToby Isaac PetscInt dimEmbed; 19c0d900a5SMatthew G. Knepley } PetscFECellGeom; 20c0d900a5SMatthew G. Knepley 214d0b9603SSander Arens typedef struct { 224d0b9603SSander Arens PetscReal v0[3]; 234d0b9603SSander Arens PetscReal J[9]; /* Size could be reduced to 6 */ 244d0b9603SSander Arens PetscReal invJ[2][9]; /* invJ in the adjacent cells */ 254d0b9603SSander Arens PetscReal detJ; /* Not really the determinant of J, but \sqrt{\det{J^T J}} */ 264d0b9603SSander Arens PetscReal n[3]; 274d0b9603SSander Arens PetscInt dim; 284d0b9603SSander Arens PetscInt dimEmbed; 294d0b9603SSander Arens PetscInt face[2]; /* The local face numbers in the adjacent cells */ 304d0b9603SSander Arens } PetscFEFaceGeom; 314d0b9603SSander Arens 32dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEInitializePackage(void); 33dbe77d9eSMatthew G. Knepley 34dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCSPACE_CLASSID; 35dbe77d9eSMatthew G. Knepley 36dbe77d9eSMatthew G. Knepley /*J 37dbe77d9eSMatthew G. Knepley PetscSpaceType - String with the name of a PETSc linear space 38dbe77d9eSMatthew G. Knepley 39dbe77d9eSMatthew G. Knepley Level: beginner 40dbe77d9eSMatthew G. Knepley 41dbe77d9eSMatthew G. Knepley .seealso: PetscSpaceSetType(), PetscSpace 42dbe77d9eSMatthew G. Knepley J*/ 43dbe77d9eSMatthew G. Knepley typedef const char *PetscSpaceType; 44dbe77d9eSMatthew G. Knepley #define PETSCSPACEPOLYNOMIAL "poly" 459c3cf19fSMatthew G. Knepley #define PETSCSPACEPOINT "point" 46*2f5fb066SToby Isaac #define PETSCSPACESUBSPACE "subspace" 47dbe77d9eSMatthew G. Knepley 48dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscSpaceList; 49dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceCreate(MPI_Comm, PetscSpace *); 50568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceDestroy(PetscSpace *); 51dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetType(PetscSpace, PetscSpaceType); 52dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetType(PetscSpace, PetscSpaceType *); 539a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetUp(PetscSpace); 54dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetFromOptions(PetscSpace); 558aec7d55SBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscSpaceViewFromOptions(PetscSpace A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);} 568aec7d55SBarry Smith 57fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscSpaceView(PetscSpace,PetscViewer); 58dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceRegister(const char [], PetscErrorCode (*)(PetscSpace)); 59dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceRegisterDestroy(void); 60dbe77d9eSMatthew G. Knepley 619a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetDimension(PetscSpace, PetscInt *); 629c3cf19fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetNumComponents(PetscSpace, PetscInt); 639c3cf19fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetNumComponents(PetscSpace, PetscInt *); 64157782e2SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceSetNumVariables(PetscSpace, PetscInt); 65157782e2SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceGetNumVariables(PetscSpace, PetscInt *); 669a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceSetOrder(PetscSpace, PetscInt); 679a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetOrder(PetscSpace, PetscInt *); 682bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceEvaluate(PetscSpace, PetscInt, const PetscReal[], PetscReal[], PetscReal[], PetscReal[]); 690c16f5adSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpaceGetHeightSubspace(PetscSpace, PetscInt, PetscSpace *); 709a1a3eb8SMatthew G. Knepley 712bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace, PetscBool); 722bdb15eaSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace, PetscBool *); 73c43d9e4aSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace, PetscBool); 74c43d9e4aSToby Isaac PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace, PetscBool *); 752bdb15eaSMatthew G. Knepley 768049c7f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePointGetPoints(PetscSpace, PetscQuadrature *); 778049c7f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscSpacePointSetPoints(PetscSpace, PetscQuadrature); 789a1a3eb8SMatthew G. Knepley 79*2f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscSpaceCreateSubspace(PetscSpace, PetscDualSpace, PetscReal *, PetscReal *, PetscReal *, PetscReal *, PetscCopyMode, PetscSpace *); 80*2f5fb066SToby Isaac 81dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCDUALSPACE_CLASSID; 82dbe77d9eSMatthew G. Knepley 83dbe77d9eSMatthew G. Knepley /*J 84dbe77d9eSMatthew G. Knepley PetscDualSpaceType - String with the name of a PETSc dual space 85dbe77d9eSMatthew G. Knepley 86dbe77d9eSMatthew G. Knepley Level: beginner 87dbe77d9eSMatthew G. Knepley 88dbe77d9eSMatthew G. Knepley .seealso: PetscDualSpaceSetType(), PetscDualSpace 89dbe77d9eSMatthew G. Knepley J*/ 90dbe77d9eSMatthew G. Knepley typedef const char *PetscDualSpaceType; 91dbe77d9eSMatthew G. Knepley #define PETSCDUALSPACELAGRANGE "lagrange" 92c2765ee2SMatthew G. Knepley #define PETSCDUALSPACESIMPLE "simple" 93dbe77d9eSMatthew G. Knepley 94dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscDualSpaceList; 95dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *); 96568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *); 97907df5ccSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace, PetscDualSpace *); 98dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetType(PetscDualSpace, PetscDualSpaceType); 99dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetType(PetscDualSpace, PetscDualSpaceType *); 1008d2f55e7SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace, const PetscInt **); 101*2f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateSection(PetscDualSpace, PetscSection *); 1029a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace); 103dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace); 1048aec7d55SBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);} 1058aec7d55SBarry Smith 106fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscDualSpaceView(PetscDualSpace,PetscViewer); 107dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceRegister(const char [], PetscErrorCode (*)(PetscDualSpace)); 108dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterDestroy(void); 109dbe77d9eSMatthew G. Knepley 110dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace, PetscInt *); 1119c3cf19fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace, PetscInt); 1129c3cf19fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace, PetscInt *); 1139a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace, PetscInt); 1149a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace, PetscInt *); 1159a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace, DM); 1169a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace, DM *); 117ebac44aeSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace, PetscInt, PetscQuadrature *); 1180ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace, PetscInt, PetscBool, DM *); 119fffc2ff8SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace, const PetscInt ****, const PetscScalar ****); 120dbe77d9eSMatthew G. Knepley 121476ae45aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscReal, PetscFECellGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *, PetscScalar *); 12255b905e8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace, PetscInt, PetscReal, PetscFECellGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *, PetscScalar *); 123ebac44aeSMatthew G. Knepley 12419ded50bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace, PetscBool *); 12519ded50bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace, PetscBool); 12689cad2ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace, PetscBool *); 12789cad2ffSToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace, PetscBool); 12819ded50bSMatthew G. Knepley 1296b594bd0SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace,PetscInt,PetscDualSpace *); 130*2f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace,PetscInt,PetscDualSpace *); 131c2765ee2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace, PetscInt); 132c2765ee2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace, PetscInt, PetscQuadrature); 1336b594bd0SToby Isaac 134dbe77d9eSMatthew G. Knepley PETSC_EXTERN PetscClassId PETSCFE_CLASSID; 135dbe77d9eSMatthew G. Knepley 1360483ade4SMatthew G. Knepley /*J 1370483ade4SMatthew G. Knepley PetscFEType - String with the name of a PETSc finite element space 1380483ade4SMatthew G. Knepley 1390483ade4SMatthew G. Knepley Level: beginner 1400483ade4SMatthew G. Knepley 1410483ade4SMatthew G. Knepley Note: Currently, the classes are concerned with the implementation of element integration 1420483ade4SMatthew G. Knepley 1430483ade4SMatthew G. Knepley .seealso: PetscFESetType(), PetscFE 1440483ade4SMatthew G. Knepley J*/ 1450483ade4SMatthew G. Knepley typedef const char *PetscFEType; 1460483ade4SMatthew G. Knepley #define PETSCFEBASIC "basic" 1473dff1f76SMatthew G. Knepley #define PETSCFENONAFFINE "nonaffine" 1480483ade4SMatthew G. Knepley #define PETSCFEOPENCL "opencl" 149aaf1837cSMatthew G. Knepley #define PETSCFECOMPOSITE "composite" 1500483ade4SMatthew G. Knepley 1510483ade4SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PetscFEList; 152568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreate(MPI_Comm, PetscFE *); 153568c31deSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEDestroy(PetscFE *); 1540483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetType(PetscFE, PetscFEType); 1550483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetType(PetscFE, PetscFEType *); 1560483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetUp(PetscFE); 1570483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetFromOptions(PetscFE); 1588aec7d55SBarry Smith PETSC_STATIC_INLINE PetscErrorCode PetscFEViewFromOptions(PetscFE A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);} 1598aec7d55SBarry Smith 160fd41f350SJed Brown PETSC_EXTERN PetscErrorCode PetscFEView(PetscFE,PetscViewer); 1610483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegister(const char [], PetscErrorCode (*)(PetscFE)); 1620483ade4SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void); 163a8e86d1aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECreateDefault(DM, PetscInt, PetscInt, PetscBool, const char [], PetscInt, PetscFE *); 1649a1a3eb8SMatthew G. Knepley 1659a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *); 1660ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *); 1679a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt); 1689a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *); 1699a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 1709a769b6bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt); 1719a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace); 1729a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *); 1739a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace); 1749a1a3eb8SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *); 175bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature); 176bfa639d9SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *); 1771bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFESetFaceQuadrature(PetscFE, PetscQuadrature); 1781bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFEGetFaceQuadrature(PetscFE, PetscQuadrature *); 1790ddb9b0bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **); 180a319912fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetDefaultTabulation(PetscFE, PetscReal **, PetscReal **, PetscReal **); 1811bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFEGetFaceTabulation(PetscFE, PetscReal **, PetscReal **, PetscReal **); 1821bafc85dSSander Arens PETSC_EXTERN PetscErrorCode PetscFEGetFaceCentroidTabulation(PetscFE, PetscReal **); 183a319912fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **); 184a319912fSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERestoreTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **); 185bceba477SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFERefine(PetscFE, PetscFE *); 186228cfb18SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEGetHeightSubspace(PetscFE, PetscInt, PetscFE *); 187a0845e3aSMatthew G. Knepley 188822ee3b2SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrate(PetscFE, PetscDS, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]); 18911dd639bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscFE, PetscDS, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]); 1904d0b9603SSander Arens PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscFE, PetscDS, PetscInt, PetscInt, PetscFEFaceGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscScalar[]); 19111dd639bSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscFE, PetscDS, PetscFEJacobianType, PetscInt, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]); 1924d0b9603SSander Arens PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdJacobian(PetscFE, PetscDS, PetscInt, PetscInt, PetscInt, PetscFEFaceGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscReal, PetscReal, PetscScalar[]); 193a0845e3aSMatthew G. Knepley 19489710940SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFECompositeGetMapping(PetscFE, PetscInt *, const PetscReal *[], const PetscReal *[], const PetscReal *[]); 19589710940SMatthew G. Knepley 196*2f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreateHeightTrace(PetscFE, PetscInt, PetscFE *); 197*2f5fb066SToby Isaac PETSC_EXTERN PetscErrorCode PetscFECreatePointTrace(PetscFE, PetscInt, PetscFE *); 198*2f5fb066SToby Isaac 199855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType); 200855cd083SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *); 201855cd083SMatthew G. Knepley 202dbe77d9eSMatthew G. Knepley #endif 203