Actual source code: petscfvimpl.h

petsc-3.9.4 2018-09-11
Report Typos and Errors
  1: #if !defined(_PETSCFVIMPL_H)
  2: #define _PETSCFVIMPL_H

  4:  #include <petscfv.h>
  5:  #include <petsc/private/petscimpl.h>

  7: PETSC_EXTERN PetscBool PetscLimiterRegisterAllCalled;
  8: PETSC_EXTERN PetscBool PetscFVRegisterAllCalled;
  9: PETSC_EXTERN PetscErrorCode PetscLimiterRegisterAll(void);
 10: PETSC_EXTERN PetscErrorCode PetscFVRegisterAll(void);

 12: typedef struct _PetscLimiterOps *PetscLimiterOps;
 13: struct _PetscLimiterOps {
 14:   PetscErrorCode (*setfromoptions)(PetscLimiter);
 15:   PetscErrorCode (*setup)(PetscLimiter);
 16:   PetscErrorCode (*view)(PetscLimiter,PetscViewer);
 17:   PetscErrorCode (*destroy)(PetscLimiter);
 18:   PetscErrorCode (*limit)(PetscLimiter, PetscReal, PetscReal *);
 19: };

 21: struct _p_PetscLimiter {
 22:   PETSCHEADER(struct _PetscLimiterOps);
 23:   void           *data;             /* Implementation object */
 24: };

 26: typedef struct {
 27:   PetscInt dummy;
 28: } PetscLimiter_Sin;

 30: typedef struct {
 31:   PetscInt dummy;
 32: } PetscLimiter_Zero;

 34: typedef struct {
 35:   PetscInt dummy;
 36: } PetscLimiter_None;

 38: typedef struct {
 39:   PetscInt dummy;
 40: } PetscLimiter_Minmod;

 42: typedef struct {
 43:   PetscInt dummy;
 44: } PetscLimiter_VanLeer;

 46: typedef struct {
 47:   PetscInt dummy;
 48: } PetscLimiter_VanAlbada;

 50: typedef struct {
 51:   PetscInt dummy;
 52: } PetscLimiter_Superbee;

 54: typedef struct {
 55:   PetscInt dummy;
 56: } PetscLimiter_MC;

 58: typedef struct _PetscFVOps *PetscFVOps;
 59: struct _PetscFVOps {
 60:   PetscErrorCode (*setfromoptions)(PetscFV);
 61:   PetscErrorCode (*setup)(PetscFV);
 62:   PetscErrorCode (*view)(PetscFV,PetscViewer);
 63:   PetscErrorCode (*destroy)(PetscFV);
 64:   PetscErrorCode (*computegradient)(PetscFV, PetscInt, const PetscScalar[], PetscScalar []);
 65:   PetscErrorCode (*integraterhsfunction)(PetscFV, PetscDS, PetscInt, PetscInt, PetscFVFaceGeom *, PetscReal *, PetscScalar[], PetscScalar[], PetscScalar[], PetscScalar[]);
 66: };

 68: struct _p_PetscFV {
 69:   PETSCHEADER(struct _PetscFVOps);
 70:   void           *data;             /* Implementation object */
 71:   PetscLimiter    limiter;          /* The slope limiter */
 72:   PetscDualSpace  dualSpace;        /* The dual space P', usually simple */
 73:   PetscInt        numComponents;    /* The number of field components */
 74:   PetscInt        dim;              /* The spatial dimension */
 75:   PetscBool       computeGradients; /* Flag for gradient computation */
 76:   PetscScalar    *fluxWork;         /* The work array for flux calculation */
 77:   PetscQuadrature quadrature;       /* Suitable quadrature on the volume */
 78:   PetscReal      *B, *D, *H;        /* Tabulation of pseudo-basis and derivatives at quadrature points */
 79:   char          **componentNames;   /* Names of the component fields */
 80: };

 82: typedef struct {
 83:   PetscInt cellType;
 84: } PetscFV_Upwind;

 86: typedef struct {
 87:   PetscInt     maxFaces, workSize;
 88:   PetscScalar *B, *Binv, *tau, *work;
 89: } PetscFV_LeastSquares;

 91: PETSC_STATIC_INLINE PetscErrorCode PetscFVInterpolate_Static(PetscFV fv, const PetscScalar x[], PetscInt q, PetscScalar interpolant[])
 92: {
 93:   PetscInt       Nc, fc;

 97:   PetscFVGetNumComponents(fv, &Nc);
 98:   for (fc = 0; fc < Nc; ++fc) {interpolant[fc] = x[fc];}
 99:   return(0);
100: }

102: #endif