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