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 #define PETSCSPACEDG "dg" 24 25 PETSC_EXTERN PetscFunctionList PetscSpaceList; 26 PETSC_EXTERN PetscBool PetscSpaceRegisterAllCalled; 27 PETSC_EXTERN PetscErrorCode PetscSpaceCreate(MPI_Comm, PetscSpace *); 28 PETSC_EXTERN PetscErrorCode PetscSpaceDestroy(PetscSpace *); 29 PETSC_EXTERN PetscErrorCode PetscSpaceSetType(PetscSpace, PetscSpaceType); 30 PETSC_EXTERN PetscErrorCode PetscSpaceGetType(PetscSpace, PetscSpaceType *); 31 PETSC_EXTERN PetscErrorCode PetscSpaceSetUp(PetscSpace); 32 PETSC_EXTERN PetscErrorCode PetscSpaceSetFromOptions(PetscSpace); 33 PETSC_EXTERN PetscErrorCode PetscSpaceRegister(const char [], PetscErrorCode (*)(PetscSpace)); 34 PETSC_EXTERN PetscErrorCode PetscSpaceRegisterAll(void); 35 PETSC_EXTERN PetscErrorCode PetscSpaceRegisterDestroy(void); 36 37 PETSC_EXTERN PetscErrorCode PetscSpaceGetDimension(PetscSpace, PetscInt *); 38 PETSC_EXTERN PetscErrorCode PetscSpaceSetOrder(PetscSpace, PetscInt); 39 PETSC_EXTERN PetscErrorCode PetscSpaceGetOrder(PetscSpace, PetscInt *); 40 PETSC_EXTERN PetscErrorCode PetscSpaceEvaluate(PetscSpace, PetscInt, const PetscReal[], PetscReal[], PetscReal[], PetscReal[]); 41 42 PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetNumVariables(PetscSpace, PetscInt); 43 PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetNumVariables(PetscSpace, PetscInt *); 44 PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace, PetscBool); 45 PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace, PetscBool *); 46 47 PETSC_EXTERN PetscErrorCode PetscSpaceDGSetQuadrature(PetscSpace, PetscQuadrature); 48 PETSC_EXTERN PetscErrorCode PetscSpaceDGGetQuadrature(PetscSpace, PetscQuadrature *); 49 50 PETSC_EXTERN PetscClassId PETSCDUALSPACE_CLASSID; 51 52 /*J 53 PetscDualSpaceType - String with the name of a PETSc dual space 54 55 Level: beginner 56 57 .seealso: PetscDualSpaceSetType(), PetscDualSpace 58 J*/ 59 typedef const char *PetscDualSpaceType; 60 #define PETSCDUALSPACELAGRANGE "lagrange" 61 62 PETSC_EXTERN PetscFunctionList PetscDualSpaceList; 63 PETSC_EXTERN PetscBool PetscDualSpaceRegisterAllCalled; 64 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *); 65 PETSC_EXTERN PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *); 66 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetType(PetscDualSpace, PetscDualSpaceType); 67 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetType(PetscDualSpace, PetscDualSpaceType *); 68 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace); 69 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace); 70 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegister(const char [], PetscErrorCode (*)(PetscDualSpace)); 71 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterAll(void); 72 PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterDestroy(void); 73 74 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace, PetscInt *); 75 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace, PetscInt); 76 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace, PetscInt *); 77 PETSC_EXTERN PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace, DM); 78 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace, DM *); 79 PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace, PetscInt, PetscQuadrature *); 80 PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace, PetscInt, PetscBool, DM *); 81 82 PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscCellGeometry, PetscInt, void (*)(const PetscReal [], PetscScalar *, void *), void *, PetscScalar *); 83 84 PETSC_EXTERN PetscClassId PETSCFE_CLASSID; 85 86 /*J 87 PetscFEType - String with the name of a PETSc finite element space 88 89 Level: beginner 90 91 Note: Currently, the classes are concerned with the implementation of element integration 92 93 .seealso: PetscFESetType(), PetscFE 94 J*/ 95 typedef const char *PetscFEType; 96 #define PETSCFEBASIC "basic" 97 #define PETSCFENONAFFINE "nonaffine" 98 #define PETSCFEOPENCL "opencl" 99 #define PETSCFECOMPOSITE "composite" 100 101 PETSC_EXTERN PetscFunctionList PetscFEList; 102 PETSC_EXTERN PetscBool PetscFERegisterAllCalled; 103 PETSC_EXTERN PetscErrorCode PetscFECreate(MPI_Comm, PetscFE *); 104 PETSC_EXTERN PetscErrorCode PetscFEDestroy(PetscFE *); 105 PETSC_EXTERN PetscErrorCode PetscFESetType(PetscFE, PetscFEType); 106 PETSC_EXTERN PetscErrorCode PetscFEGetType(PetscFE, PetscFEType *); 107 PETSC_EXTERN PetscErrorCode PetscFESetUp(PetscFE); 108 PETSC_EXTERN PetscErrorCode PetscFESetFromOptions(PetscFE); 109 PETSC_EXTERN PetscErrorCode PetscFERegister(const char [], PetscErrorCode (*)(PetscFE)); 110 PETSC_EXTERN PetscErrorCode PetscFERegisterAll(void); 111 PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void); 112 PETSC_EXTERN PetscErrorCode PetscFECreateDefault(DM, PetscInt, PetscInt, PetscBool, const char [], PetscInt, PetscFE *); 113 114 PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *); 115 PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *); 116 PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt); 117 PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *); 118 PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *); 119 PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt); 120 PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace); 121 PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *); 122 PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace); 123 PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *); 124 PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature); 125 PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *); 126 PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **); 127 PETSC_EXTERN PetscErrorCode PetscFEGetDefaultTabulation(PetscFE, PetscReal **, PetscReal **, PetscReal **); 128 PETSC_EXTERN PetscErrorCode PetscFEGetTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **); 129 PETSC_EXTERN PetscErrorCode PetscFERestoreTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **); 130 PETSC_EXTERN PetscErrorCode PetscFERefine(PetscFE, PetscFE *); 131 132 PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[], 133 PetscInt, PetscFE[], const PetscScalar[], 134 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 135 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 136 PetscScalar[]); 137 PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[], 138 PetscInt, PetscFE[], const PetscScalar[], 139 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]), 140 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]), 141 PetscScalar[]); 142 PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobianAction(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[], const PetscScalar[], 143 void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 144 void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 145 void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 146 void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 147 PetscScalar[]); 148 PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscInt, PetscCellGeometry, const PetscScalar[], 149 PetscInt, PetscFE[], const PetscScalar[], 150 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 151 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 152 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 153 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 154 PetscScalar[]); 155 PETSC_EXTERN PetscErrorCode PetscFEIntegrateIFunction(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[], const PetscScalar[], 156 PetscInt, PetscFE[], const PetscScalar[], 157 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 158 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]), 159 PetscScalar[]); 160 PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdIFunction(PetscFE, PetscInt, PetscInt, PetscFE[], PetscInt, PetscCellGeometry, const PetscScalar[], const PetscScalar[], 161 PetscInt, PetscFE[], const PetscScalar[], 162 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]), 163 void (*)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]), 164 PetscScalar[]); 165 166 PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType); 167 PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *); 168 169 /* TODO: Should be moved inside DM */ 170 typedef struct { 171 PetscFE *fe; 172 PetscFE *feAux; 173 void (**f0Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */ 174 void (**f1Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */ 175 void (**f0IFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */ 176 void (**f1IFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */ 177 void (**g0Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_0 functions for each field pair */ 178 void (**g1Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_1 functions for each field pair */ 179 void (**g2Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_2 functions for each field pair */ 180 void (**g3Funcs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], PetscScalar[]); /* The g_3 functions for each field pair */ 181 void (**bcFuncs)(const PetscReal[], PetscScalar *, void *); /* The boundary condition function for each field component */ 182 void **bcCtxs; /* Contexts for each boundary condition function, or null if all contexts are null */ 183 PetscFE *feBd; 184 void (**f0BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */ 185 void (**f1BdFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */ 186 void (**f0BdIFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_0 functions for each field */ 187 void (**f1BdIFuncs)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal[], const PetscReal[], PetscScalar[]); /* The f_1 functions for each field */ 188 } PetscFEM; 189 190 #endif 191