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