Actual source code: petscfe.h

petsc-3.7.3 2016-08-01
Report Typos and Errors
  1: /*
  2:       Objects which encapsulate finite element spaces and operations
  3: */
  6: #include <petscdm.h>
  7: #include <petscdt.h>
  8: #include <petscfetypes.h>
  9: #include <petscdstypes.h>

 11: /* Assuming dim <= 3 */
 12: typedef struct {
 13:   PetscReal v0[3];
 14:   PetscReal J[9];
 15:   PetscReal invJ[9];
 16:   PetscReal detJ;
 17:   PetscReal n[3];
 18:   PetscInt  dim;
 19:   PetscInt  dimEmbed;
 20: } PetscFECellGeom;

 22: PETSC_EXTERN PetscErrorCode PetscFEInitializePackage(void);

 24: PETSC_EXTERN PetscClassId PETSCSPACE_CLASSID;

 26: /*J
 27:   PetscSpaceType - String with the name of a PETSc linear space

 29:   Level: beginner

 31: .seealso: PetscSpaceSetType(), PetscSpace
 32: J*/
 33: typedef const char *PetscSpaceType;
 34: #define PETSCSPACEPOLYNOMIAL "poly"
 35: #define PETSCSPACEDG         "dg"

 37: PETSC_EXTERN PetscFunctionList PetscSpaceList;
 38: PETSC_EXTERN PetscErrorCode PetscSpaceCreate(MPI_Comm, PetscSpace *);
 39: PETSC_EXTERN PetscErrorCode PetscSpaceDestroy(PetscSpace *);
 40: PETSC_EXTERN PetscErrorCode PetscSpaceSetType(PetscSpace, PetscSpaceType);
 41: PETSC_EXTERN PetscErrorCode PetscSpaceGetType(PetscSpace, PetscSpaceType *);
 42: PETSC_EXTERN PetscErrorCode PetscSpaceSetUp(PetscSpace);
 43: PETSC_EXTERN PetscErrorCode PetscSpaceSetFromOptions(PetscSpace);
 44: PETSC_STATIC_INLINE PetscErrorCode PetscSpaceViewFromOptions(PetscSpace A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);}

 46: PETSC_EXTERN PetscErrorCode PetscSpaceView(PetscSpace,PetscViewer);
 47: PETSC_EXTERN PetscErrorCode PetscSpaceRegister(const char [], PetscErrorCode (*)(PetscSpace));
 48: PETSC_EXTERN PetscErrorCode PetscSpaceRegisterDestroy(void);

 50: PETSC_EXTERN PetscErrorCode PetscSpaceGetDimension(PetscSpace, PetscInt *);
 51: PETSC_EXTERN PetscErrorCode PetscSpaceSetOrder(PetscSpace, PetscInt);
 52: PETSC_EXTERN PetscErrorCode PetscSpaceGetOrder(PetscSpace, PetscInt *);
 53: PETSC_EXTERN PetscErrorCode PetscSpaceEvaluate(PetscSpace, PetscInt, const PetscReal[], PetscReal[], PetscReal[], PetscReal[]);

 55: PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetNumVariables(PetscSpace, PetscInt);
 56: PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetNumVariables(PetscSpace, PetscInt *);
 57: PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetSymmetric(PetscSpace, PetscBool);
 58: PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetSymmetric(PetscSpace, PetscBool *);
 59: PETSC_EXTERN PetscErrorCode PetscSpacePolynomialSetTensor(PetscSpace, PetscBool);
 60: PETSC_EXTERN PetscErrorCode PetscSpacePolynomialGetTensor(PetscSpace, PetscBool *);

 62: PETSC_EXTERN PetscErrorCode PetscSpaceDGSetQuadrature(PetscSpace, PetscQuadrature);
 63: PETSC_EXTERN PetscErrorCode PetscSpaceDGGetQuadrature(PetscSpace, PetscQuadrature *);

 65: PETSC_EXTERN PetscClassId PETSCDUALSPACE_CLASSID;

 67: /*J
 68:   PetscDualSpaceType - String with the name of a PETSc dual space

 70:   Level: beginner

 72: .seealso: PetscDualSpaceSetType(), PetscDualSpace
 73: J*/
 74: typedef const char *PetscDualSpaceType;
 75: #define PETSCDUALSPACELAGRANGE "lagrange"
 76: #define PETSCDUALSPACESIMPLE   "simple"

 78: PETSC_EXTERN PetscFunctionList PetscDualSpaceList;
 79: PETSC_EXTERN PetscErrorCode PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *);
 80: PETSC_EXTERN PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *);
 81: PETSC_EXTERN PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace, PetscDualSpace *);
 82: PETSC_EXTERN PetscErrorCode PetscDualSpaceSetType(PetscDualSpace, PetscDualSpaceType);
 83: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetType(PetscDualSpace, PetscDualSpaceType *);
 84: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace, const PetscInt **);
 85: PETSC_EXTERN PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace);
 86: PETSC_EXTERN PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace);
 87: PETSC_STATIC_INLINE PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);}

 89: PETSC_EXTERN PetscErrorCode PetscDualSpaceView(PetscDualSpace,PetscViewer);
 90: PETSC_EXTERN PetscErrorCode PetscDualSpaceRegister(const char [], PetscErrorCode (*)(PetscDualSpace));
 91: PETSC_EXTERN PetscErrorCode PetscDualSpaceRegisterDestroy(void);

 93: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace, PetscInt *);
 94: PETSC_EXTERN PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace, PetscInt);
 95: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace, PetscInt *);
 96: PETSC_EXTERN PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace, DM);
 97: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace, DM *);
 98: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace, PetscInt, PetscQuadrature *);
 99: PETSC_EXTERN PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace, PetscInt, PetscBool, DM *);

101: PETSC_EXTERN PetscErrorCode PetscDualSpaceApply(PetscDualSpace, PetscInt, PetscReal, PetscFECellGeom *, PetscInt, PetscErrorCode (*)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *, PetscScalar *);

103: PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeGetContinuity(PetscDualSpace, PetscBool *);
104: PETSC_EXTERN PetscErrorCode PetscDualSpaceLagrangeSetContinuity(PetscDualSpace, PetscBool);

106: PETSC_EXTERN PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace,PetscInt,PetscDualSpace *);
107: PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetDimension(PetscDualSpace, PetscInt);
108: PETSC_EXTERN PetscErrorCode PetscDualSpaceSimpleSetFunctional(PetscDualSpace, PetscInt, PetscQuadrature);

110: PETSC_EXTERN PetscClassId PETSCFE_CLASSID;

112: /*J
113:   PetscFEType - String with the name of a PETSc finite element space

115:   Level: beginner

117:   Note: Currently, the classes are concerned with the implementation of element integration

119: .seealso: PetscFESetType(), PetscFE
120: J*/
121: typedef const char *PetscFEType;
122: #define PETSCFEBASIC     "basic"
123: #define PETSCFENONAFFINE "nonaffine"
124: #define PETSCFEOPENCL    "opencl"
125: #define PETSCFECOMPOSITE "composite"

127: PETSC_EXTERN PetscFunctionList PetscFEList;
128: PETSC_EXTERN PetscErrorCode PetscFECreate(MPI_Comm, PetscFE *);
129: PETSC_EXTERN PetscErrorCode PetscFEDestroy(PetscFE *);
130: PETSC_EXTERN PetscErrorCode PetscFESetType(PetscFE, PetscFEType);
131: PETSC_EXTERN PetscErrorCode PetscFEGetType(PetscFE, PetscFEType *);
132: PETSC_EXTERN PetscErrorCode PetscFESetUp(PetscFE);
133: PETSC_EXTERN PetscErrorCode PetscFESetFromOptions(PetscFE);
134: PETSC_STATIC_INLINE PetscErrorCode PetscFEViewFromOptions(PetscFE A,PetscObject B,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,B,name);}

136: PETSC_EXTERN PetscErrorCode PetscFEView(PetscFE,PetscViewer);
137: PETSC_EXTERN PetscErrorCode PetscFERegister(const char [], PetscErrorCode (*)(PetscFE));
138: PETSC_EXTERN PetscErrorCode PetscFERegisterDestroy(void);
139: PETSC_EXTERN PetscErrorCode PetscFECreateDefault(DM, PetscInt, PetscInt, PetscBool, const char [], PetscInt, PetscFE *);

141: PETSC_EXTERN PetscErrorCode PetscFEGetDimension(PetscFE, PetscInt *);
142: PETSC_EXTERN PetscErrorCode PetscFEGetSpatialDimension(PetscFE, PetscInt *);
143: PETSC_EXTERN PetscErrorCode PetscFESetNumComponents(PetscFE, PetscInt);
144: PETSC_EXTERN PetscErrorCode PetscFEGetNumComponents(PetscFE, PetscInt *);
145: PETSC_EXTERN PetscErrorCode PetscFEGetTileSizes(PetscFE, PetscInt *, PetscInt *, PetscInt *, PetscInt *);
146: PETSC_EXTERN PetscErrorCode PetscFESetTileSizes(PetscFE, PetscInt, PetscInt, PetscInt, PetscInt);
147: PETSC_EXTERN PetscErrorCode PetscFESetBasisSpace(PetscFE, PetscSpace);
148: PETSC_EXTERN PetscErrorCode PetscFEGetBasisSpace(PetscFE, PetscSpace *);
149: PETSC_EXTERN PetscErrorCode PetscFESetDualSpace(PetscFE, PetscDualSpace);
150: PETSC_EXTERN PetscErrorCode PetscFEGetDualSpace(PetscFE, PetscDualSpace *);
151: PETSC_EXTERN PetscErrorCode PetscFESetQuadrature(PetscFE, PetscQuadrature);
152: PETSC_EXTERN PetscErrorCode PetscFEGetQuadrature(PetscFE, PetscQuadrature *);
153: PETSC_EXTERN PetscErrorCode PetscFEGetNumDof(PetscFE, const PetscInt **);
154: PETSC_EXTERN PetscErrorCode PetscFEGetDefaultTabulation(PetscFE, PetscReal **, PetscReal **, PetscReal **);
155: PETSC_EXTERN PetscErrorCode PetscFEGetFaceTabulation(PetscFE, PetscReal **);
156: PETSC_EXTERN PetscErrorCode PetscFEGetTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **);
157: PETSC_EXTERN PetscErrorCode PetscFERestoreTabulation(PetscFE, PetscInt, const PetscReal[], PetscReal **, PetscReal **, PetscReal **);
158: PETSC_EXTERN PetscErrorCode PetscFERefine(PetscFE, PetscFE *);

160: PETSC_EXTERN PetscErrorCode PetscFEIntegrate(PetscFE, PetscDS, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], PetscDS, const PetscScalar[], PetscReal[]);
161: PETSC_EXTERN PetscErrorCode PetscFEIntegrateResidual(PetscFE, PetscDS, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
162: PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdResidual(PetscFE, PetscDS, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
163: PETSC_EXTERN PetscErrorCode PetscFEIntegrateJacobian(PetscFE, PetscDS, PetscFEJacobianType, PetscInt, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);
164: PETSC_EXTERN PetscErrorCode PetscFEIntegrateBdJacobian(PetscFE, PetscDS, PetscInt, PetscInt, PetscInt, PetscFECellGeom *, const PetscScalar[], const PetscScalar[], PetscDS, const PetscScalar[], PetscScalar[]);

166: PETSC_EXTERN PetscErrorCode PetscFECompositeGetMapping(PetscFE, PetscInt *, const PetscReal *[], const PetscReal *[], const PetscReal *[]);

168: PETSC_EXTERN PetscErrorCode PetscFEOpenCLSetRealType(PetscFE, PetscDataType);
169: PETSC_EXTERN PetscErrorCode PetscFEOpenCLGetRealType(PetscFE, PetscDataType *);

171: #endif