xref: /petsc/include/petscdt.h (revision a0845e3a928e8c3de76a3c8bfba8b69f6cc922fe)
137045ce4SJed Brown /*
237045ce4SJed Brown   Common tools for constructing discretizations
337045ce4SJed Brown */
437045ce4SJed Brown #if !defined(__PETSCDT_H)
537045ce4SJed Brown #define __PETSCDT_H
637045ce4SJed Brown 
737045ce4SJed Brown #include <petscsys.h>
837045ce4SJed Brown 
9*a0845e3aSMatthew G. Knepley typedef struct {
10*a0845e3aSMatthew G. Knepley   PetscInt         numQuadPoints; /* The number of quadrature points on an element */
11*a0845e3aSMatthew G. Knepley   const PetscReal *quadPoints;    /* The quadrature point coordinates */
12*a0845e3aSMatthew G. Knepley   const PetscReal *quadWeights;   /* The quadrature weights */
13*a0845e3aSMatthew G. Knepley   PetscInt         numBasisFuncs; /* The number of finite element basis functions on an element */
14*a0845e3aSMatthew G. Knepley   PetscInt         numComponents; /* The number of components for each basis function */
15*a0845e3aSMatthew G. Knepley   const PetscReal *basis;         /* The basis functions tabulated at the quadrature points */
16*a0845e3aSMatthew G. Knepley   const PetscReal *basisDer;      /* The basis function derivatives tabulated at the quadrature points */
17*a0845e3aSMatthew G. Knepley } PetscQuadrature;
18*a0845e3aSMatthew G. Knepley 
19*a0845e3aSMatthew G. Knepley typedef struct {
20*a0845e3aSMatthew G. Knepley   PetscReal *v0, *J, *invJ, *detJ;
21*a0845e3aSMatthew G. Knepley } PetscCellGeometry;
22*a0845e3aSMatthew G. Knepley 
2337045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
2437045ce4SJed Brown PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*);
25194825f6SJed Brown PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*);
26*a0845e3aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
27*a0845e3aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature*);
2837045ce4SJed Brown 
2937045ce4SJed Brown #endif
30