xref: /petsc/include/petsc/private/petscfvimpl.h (revision 4796d426df8bce7e71eb92de021c7fe5f7327dc3)
1 #ifndef PETSCFVIMPL_H
2 #define PETSCFVIMPL_H
3 
4 #include <petscfv.h>
5 #ifdef PETSC_HAVE_LIBCEED
6   #include <petscfvceed.h>
7 #endif
8 #include <petsc/private/petscimpl.h>
9 #include <petsc/private/dmimpl.h>
10 
11 PETSC_EXTERN PetscBool      PetscLimiterRegisterAllCalled;
12 PETSC_EXTERN PetscBool      PetscFVRegisterAllCalled;
13 PETSC_EXTERN PetscErrorCode PetscLimiterRegisterAll(void);
14 PETSC_EXTERN PetscErrorCode PetscFVRegisterAll(void);
15 
16 typedef struct _PetscLimiterOps *PetscLimiterOps;
17 struct _PetscLimiterOps {
18   PetscErrorCode (*setfromoptions)(PetscLimiter);
19   PetscErrorCode (*setup)(PetscLimiter);
20   PetscErrorCode (*view)(PetscLimiter, PetscViewer);
21   PetscErrorCode (*destroy)(PetscLimiter);
22   PetscErrorCode (*limit)(PetscLimiter, PetscReal, PetscReal *);
23 };
24 
25 struct _p_PetscLimiter {
26   PETSCHEADER(struct _PetscLimiterOps);
27   void *data; /* Implementation object */
28 };
29 
30 typedef struct {
31   PetscInt dummy;
32 } PetscLimiter_Sin;
33 
34 typedef struct {
35   PetscInt dummy;
36 } PetscLimiter_Zero;
37 
38 typedef struct {
39   PetscInt dummy;
40 } PetscLimiter_None;
41 
42 typedef struct {
43   PetscInt dummy;
44 } PetscLimiter_Minmod;
45 
46 typedef struct {
47   PetscInt dummy;
48 } PetscLimiter_VanLeer;
49 
50 typedef struct {
51   PetscInt dummy;
52 } PetscLimiter_VanAlbada;
53 
54 typedef struct {
55   PetscInt dummy;
56 } PetscLimiter_Superbee;
57 
58 typedef struct {
59   PetscInt dummy;
60 } PetscLimiter_MC;
61 
62 typedef struct _PetscFVOps *PetscFVOps;
63 struct _PetscFVOps {
64   PetscErrorCode (*setfromoptions)(PetscFV);
65   PetscErrorCode (*setup)(PetscFV);
66   PetscErrorCode (*view)(PetscFV, PetscViewer);
67   PetscErrorCode (*destroy)(PetscFV);
68   PetscErrorCode (*computegradient)(PetscFV, PetscInt, const PetscScalar[], PetscScalar[]);
69   PetscErrorCode (*integraterhsfunction)(PetscFV, PetscDS, PetscInt, PetscInt, PetscFVFaceGeom *, PetscReal *, PetscScalar[], PetscScalar[], PetscScalar[], PetscScalar[]);
70 };
71 
72 struct _p_PetscFV {
73   PETSCHEADER(struct _PetscFVOps);
74   void           *data;             /* Implementation object */
75   PetscLimiter    limiter;          /* The slope limiter */
76   PetscDualSpace  dualSpace;        /* The dual space P', usually simple */
77   PetscInt        numComponents;    /* The number of field components */
78   PetscInt        dim;              /* The spatial dimension */
79   PetscBool       computeGradients; /* Flag for gradient computation */
80   PetscScalar    *fluxWork;         /* The work array for flux calculation */
81   PetscQuadrature quadrature;       /* Suitable quadrature on the volume */
82   PetscTabulation T;                /* Tabulation of pseudo-basis and derivatives at quadrature points */
83   char          **componentNames;   /* Names of the component fields */
84 #ifdef PETSC_HAVE_LIBCEED
85   Ceed      ceed;      /* The LibCEED context, usually set by the DM */
86   CeedBasis ceedBasis; /* Basis for libCEED matching this element */
87 #endif
88 };
89 
90 typedef struct {
91   PetscInt cellType;
92 } PetscFV_Upwind;
93 
94 typedef struct {
95   PetscInt     maxFaces, workSize;
96   PetscScalar *B, *Binv, *tau, *work;
97 } PetscFV_LeastSquares;
98 
99 static inline PetscErrorCode PetscFVInterpolate_Static(PetscFV fv, const PetscScalar x[], PetscInt q, PetscScalar interpolant[])
100 {
101   PetscInt Nc;
102 
103   PetscFunctionBeginHot;
104   PetscCall(PetscFVGetNumComponents(fv, &Nc));
105   PetscCall(PetscArraycpy(interpolant, x, Nc));
106   PetscFunctionReturn(PETSC_SUCCESS);
107 }
108 
109 #endif
110