xref: /petsc/include/petscdt.h (revision 6a6dd048690cb4dc28d91d1ce2d9f1de8cbcadd6)
1 /*
2   Common tools for constructing discretizations
3 */
4 #if !defined(__PETSCDT_H)
5 #define __PETSCDT_H
6 
7 #include <petscsys.h>
8 
9 /*S
10   PetscQuadrature - Quadrature rule for integration.
11 
12   Level: developer
13 
14 .seealso:  PetscQuadratureCreate(), PetscQuadratureDestroy()
15 S*/
16 typedef struct _p_PetscQuadrature *PetscQuadrature;
17 
18 /*E
19   PetscGLLCreateType - algorithm used to compute the Gauss-Lobatto-Legendre nodes and weights
20 
21   Level: intermediate
22 
23 $  PETSCGLL_VIA_LINEARALGEBRA - compute the nodes via linear algebra
24 $  PETSCGLL_VIA_NEWTON - compute the nodes by solving a nonlinear equation with Newton's method
25 
26 E*/
27 typedef enum {PETSCGLL_VIA_LINEARALGEBRA,PETSCGLL_VIA_NEWTON} PetscGLLCreateType;
28 
29 PETSC_EXTERN PetscErrorCode PetscQuadratureCreate(MPI_Comm, PetscQuadrature *);
30 PETSC_EXTERN PetscErrorCode PetscQuadratureDuplicate(PetscQuadrature, PetscQuadrature *);
31 PETSC_EXTERN PetscErrorCode PetscQuadratureGetOrder(PetscQuadrature, PetscInt*);
32 PETSC_EXTERN PetscErrorCode PetscQuadratureSetOrder(PetscQuadrature, PetscInt);
33 PETSC_EXTERN PetscErrorCode PetscQuadratureGetNumComponents(PetscQuadrature, PetscInt*);
34 PETSC_EXTERN PetscErrorCode PetscQuadratureSetNumComponents(PetscQuadrature, PetscInt);
35 PETSC_EXTERN PetscErrorCode PetscQuadratureGetData(PetscQuadrature, PetscInt*, PetscInt*, PetscInt*, const PetscReal *[], const PetscReal *[]);
36 PETSC_EXTERN PetscErrorCode PetscQuadratureSetData(PetscQuadrature, PetscInt, PetscInt, PetscInt, const PetscReal [], const PetscReal []);
37 PETSC_EXTERN PetscErrorCode PetscQuadratureView(PetscQuadrature, PetscViewer);
38 PETSC_EXTERN PetscErrorCode PetscQuadratureDestroy(PetscQuadrature *);
39 
40 PETSC_EXTERN PetscErrorCode PetscQuadratureExpandComposite(PetscQuadrature, PetscInt, const PetscReal[], const PetscReal[], PetscQuadrature *);
41 
42 PETSC_EXTERN PetscErrorCode PetscDTLegendreEval(PetscInt,const PetscReal*,PetscInt,const PetscInt*,PetscReal*,PetscReal*,PetscReal*);
43 PETSC_EXTERN PetscErrorCode PetscDTGaussQuadrature(PetscInt,PetscReal,PetscReal,PetscReal*,PetscReal*);
44 PETSC_EXTERN PetscErrorCode PetscDTGaussLobattoLegendreQuadrature(PetscInt,PetscGLLCreateType,PetscReal*,PetscReal*);
45 PETSC_EXTERN PetscErrorCode PetscDTReconstructPoly(PetscInt,PetscInt,const PetscReal*,PetscInt,const PetscReal*,PetscReal*);
46 PETSC_EXTERN PetscErrorCode PetscDTGaussTensorQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
47 PETSC_EXTERN PetscErrorCode PetscDTGaussJacobiQuadrature(PetscInt,PetscInt,PetscInt,PetscReal,PetscReal,PetscQuadrature*);
48 
49 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhTensorQuadrature(PetscInt, PetscInt, PetscReal, PetscReal, PetscQuadrature *);
50 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrate(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *);
51 PETSC_EXTERN PetscErrorCode PetscDTTanhSinhIntegrateMPFR(void (*)(PetscReal, PetscReal *), PetscReal, PetscReal, PetscInt, PetscReal *);
52 
53 #endif
54