1 /* 2 Objects which encapsulate finite element spaces and operations 3 */ 4 #if !defined(__PETSCFE_H) 5 #define __PETSCFE_H 6 #include <petscdm.h> 7 #include <petscdt.h> 8 #include <petscfetypes.h> 9 10 PETSC_EXTERN PetscErrorCode PetscFEInitializePackage(void); 11 12 PETSC_EXTERN PetscClassId PETSCSPACE_CLASSID; 13 14 /*J 15 PetscSpaceType - String with the name of a PETSc linear space 16 17 Level: beginner 18 19 .seealso: PetscSpaceSetType(), PetscSpace 20 J*/ 21 typedef const char *PetscSpaceType; 22 #define PETSCSPACEPOLYNOMIAL "poly" 23 24 PETSC_EXTERN PetscFunctionList PetscSpaceList; 25 PETSC_EXTERN PetscBool PetscSpaceRegisterAllCalled; 26 PETSC_EXTERN PetscErrorCode PetscSpaceCreate(MPI_Comm, PetscSpace *); 27 PETSC_EXTERN PetscErrorCode PetscSpaceDestroy(PetscSpace *); 28 PETSC_EXTERN PetscErrorCode PetscSpaceSetType(PetscSpace, PetscSpaceType); 29 PETSC_EXTERN PetscErrorCode PetscSpaceGetType(PetscSpace, PetscSpaceType *); 30 PETSC_EXTERN PetscErrorCode PetscSpaceSetUp(PetscSpace); 31 PETSC_EXTERN PetscErrorCode PetscSpaceSetFromOptions(PetscSpace); 32 PETSC_EXTERN PetscErrorCode PetscSpaceRegister(const char [], PetscErrorCode (*)(PetscSpace)); 33 PETSC_EXTERN PetscErrorCode PetscSpaceRegisterAll(void); 34 PETSC_EXTERN PetscErrorCode PetscSpaceRegisterDestroy(void); 35 36 PETSC_EXTERN PetscErrorCode PetscSpaceGetDimension(PetscSpace, PetscInt *); 37 PETSC_EXTERN PetscErrorCode PetscSpaceSetOrder(PetscSpace, PetscInt); 38 PETSC_EXTERN PetscErrorCode PetscSpaceGetOrder(PetscSpace, PetscInt *); 39 40 PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetNumVariables(PetscSpace, PetscInt); 41 PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetNumVariables(PetscSpace, PetscInt *); 42 43 PETSC_EXTERN PetscClassId PETSCDUALSPACE_CLASSID; 44 45 /*J 46 PetscDualSpaceType - String with the name of a PETSc dual space 47 48 Level: beginner 49 50 .seealso: PetscDualSpaceSetType(), PetscDualSpace 51 J*/ 52 typedef const char *PetscDualSpaceType; 53 #define PETSCDUALSPACELAGRANGE "lagrange" 54 55 PETSC_EXTERN PetscFunctionList PetscDualSpaceList; 56 PETSC_EXTERN PetscBool PetscDualSpaceRegisterAllCalled; 57 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *); 58 PETSC_EXTERN PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *); 59 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetType(PetscDualSpace, PetscDualSpaceType); 60 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetType(PetscDualSpace, PetscDualSpaceType *); 61 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace); 62 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace); 63 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegister(const char [], PetscErrorCode (*)(PetscDualSpace)); 64 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterAll(void); 65 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterDestroy(void); 66 67 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace, PetscInt *); 68 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace, PetscInt); 69 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace, PetscInt *); 70 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace, DM); 71 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace, DM *); 72 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace, PetscInt, PetscBool, DM *); 73 74 PETSC_EXTERN PetscClassId PETSCFE_CLASSID; 75 76 /*J 77 PetscFEType - String with the name of a PETSc finite element space 78 79 Level: beginner 80 81 Note: Currently, the classes are concerned with the implementation of element integration 82 83 .seealso: PetscFESetType(), PetscFE 84 J*/ 85 typedef const char *PetscFEType; 86 #define PETSCFEBASIC "basic" 87 #define PETSCFEOPENCL "opencl" 88 89 PETSC_EXTERN PetscFunctionList PetscFEList; 90 PETSC_EXTERN PetscBool PetscFERegisterAllCalled; 91 PETSC_EXTERN PetscErrorCode PetscFECreate(MPI_Comm, PetscFE *); 92 PETSC_EXTERN PetscErrorCode PetscFEDestroy(PetscFE *); 93 PETSC_EXTERN PetscErrorCode PetscFESetType(PetscFE, PetscFEType); 94 PETSC_EXTERN PetscErrorCode PetscFEGetType(PetscFE, PetscFEType *); 95 PETSC_EXTERN PetscErrorCode PetscFESetUp(PetscFE); 96 PETSC_EXTERN PetscErrorCode PetscFESetFromOptions(PetscFE); 97 PETSC_EXTERN PetscErrorCode PetscFERegister(const char [], PetscErrorCode (*)(PetscFE)); 98 PETSC_EXTERN PetscErrorCode PetscFERegisterAll(void); 99 PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void); 100 101 PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *); 102 PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *); 103 PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt); 104 PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *); 105 PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 106 PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt); 107 PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace); 108 PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *); 109 PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace); 110 PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *); 111 PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature); 112 PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *); 113 PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **); 114 PETSC_EXTERN PetscErrorCode PetscFEGetDefaultTabulation(PetscFE, PetscReal **, PetscReal **, PetscReal **); 115 PETSC_EXTERN PetscErrorCode PetscFEGetTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **); 116 PETSC_EXTERN PetscErrorCode PetscFERestoreTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **); 117 118 PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[], 119 PetscInt, PetscFE[], const PetscScalar[], 120 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 121 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 122 PetscScalar[]); 123 PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[], 124 PetscInt, PetscFE[], const PetscScalar[], 125 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]), 126 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]), 127 PetscScalar[]); 128 PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobianAction(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[], const PetscScalar[], 129 void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 130 void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 131 void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 132 void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 133 PetscScalar[]); 134 PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscInt, PetscCellGeometry, const PetscScalar[], 135 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 136 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 137 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 138 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 139 PetscScalar[]); 140 141 PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType); 142 PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *); 143 144 /* TODO: Should be moved inside DM */ 145 typedef struct { 146 PetscFE *fe; 147 PetscFE *feAux; 148 void (**f0Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */ 149 void (**f1Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */ 150 void (**g0Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_0 functions for each field pair */ 151 void (**g1Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_1 functions for each field pair */ 152 void (**g2Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_2 functions for each field pair */ 153 void (**g3Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_3 functions for each field pair */ 154 void (**bcFuncs)(const PetscReal[], PetscScalar *); /* The boundary condition function for each field component */ 155 PetscFE *feBd; 156 void (**f0BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */ 157 void (**f1BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */ 158 } PetscFEM; 159 160 #endif 161