Actual source code: petscdsimpl.h

petsc-3.13.2 2020-06-02
Report Typos and Errors
  1: #if !defined(PETSCDSIMPL_H)
  2: #define PETSCDSIMPL_H

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

  7: PETSC_EXTERN PetscBool      PetscDSRegisterAllCalled;
  8: PETSC_EXTERN PetscErrorCode PetscDSRegisterAll(void);

 10: typedef struct _n_DSBoundary *DSBoundary;

 12: struct _n_DSBoundary {
 13:   const char *name;
 14:   const char *labelname;
 15:   DMBoundaryConditionType type;
 16:   PetscInt    field;
 17:   PetscInt    numcomps;
 18:   PetscInt   *comps;
 19:   void      (*func)(void);
 20:   PetscInt    numids;
 21:   PetscInt   *ids;
 22:   void       *ctx;
 23:   DSBoundary  next;
 24: };

 26: typedef struct _PetscDSOps *PetscDSOps;
 27: struct _PetscDSOps {
 28:   PetscErrorCode (*setfromoptions)(PetscDS);
 29:   PetscErrorCode (*setup)(PetscDS);
 30:   PetscErrorCode (*view)(PetscDS,PetscViewer);
 31:   PetscErrorCode (*destroy)(PetscDS);
 32: };

 34: struct _p_PetscDS {
 35:   PETSCHEADER(struct _PetscDSOps);
 36:   void        *data;              /* Implementation object */
 37:   PetscDS     *subprobs;          /* The subspaces for each dimension */
 38:   PetscBool    setup;             /* Flag for setup */
 39:   PetscBool    isHybrid;          /* Flag for hybrid cell (this is crappy, but the only thing I can see to do now) */
 40:   PetscInt     dimEmbed;          /* The real space coordinate dimension */
 41:   PetscInt     Nf;                /* The number of solution fields */
 42:   PetscObject *disc;              /* The discretization for each solution field (PetscFE, PetscFV, etc.) */
 43:   /* Equations */
 44:   DSBoundary            boundary;      /* Linked list of boundary conditions */
 45:   PetscBool             useJacPre;     /* Flag for using the Jacobian preconditioner */
 46:   PetscBool            *implicit;      /* Flag for implicit or explicit solve for each field */
 47:   PetscPointFunc       *obj;           /* Scalar integral (like an objective function) */
 48:   PetscPointFunc       *f;             /* Weak form integrands for F, f_0, f_1 */
 49:   PetscPointJac        *g;             /* Weak form integrands for J = dF/du, g_0, g_1, g_2, g_3 */
 50:   PetscPointJac        *gp;            /* Weak form integrands for preconditioner for J, g_0, g_1, g_2, g_3 */
 51:   PetscPointJac        *gt;            /* Weak form integrands for dF/du_t, g_0, g_1, g_2, g_3 */
 52:   PetscBdPointFunc     *fBd;           /* Weak form boundary integrands F_bd, f_0, f_1 */
 53:   PetscBdPointJac      *gBd;           /* Weak form boundary integrands J_bd = dF_bd/du, g_0, g_1, g_2, g_3 */
 54:   PetscRiemannFunc     *r;             /* Riemann solvers */
 55:   PetscPointFunc       *update;        /* Direct update of field coefficients */
 56:   PetscSimplePointFunc *exactSol;      /* Exact solutions for each field */
 57:   void                **exactCtx;      /* Contexts for the exact solution functions */
 58:   PetscInt              numConstants;  /* Number of constants passed to point functions */
 59:   PetscScalar          *constants;     /* Array of constants passed to point functions */
 60:   void                 **ctx;          /* User contexts for each field */
 61:   /* Computed sizes */
 62:   PetscInt         totDim;             /* Total system dimension */
 63:   PetscInt         totComp;            /* Total field components */
 64:   PetscInt        *Nc;                 /* Number of components for each field */
 65:   PetscInt        *Nb;                 /* Number of basis functions for each field */
 66:   PetscInt        *off;                /* Offsets for each field */
 67:   PetscInt        *offDer;             /* Derivative offsets for each field */
 68:   PetscTabulation *T;                  /* Basis function and derivative tabulation for each field */
 69:   PetscTabulation *Tf;                 /* Basis function and derivative tabulation for each local face and field */
 70:   /* Work space */
 71:   PetscScalar *u;                      /* Field evaluation */
 72:   PetscScalar *u_t;                    /* Field time derivative evaluation */
 73:   PetscScalar *u_x;                    /* Field gradient evaluation */
 74:   PetscScalar *basisReal;              /* Workspace for pushforward */
 75:   PetscScalar *basisDerReal;           /* Workspace for derivative pushforward */
 76:   PetscScalar *testReal;               /* Workspace for pushforward */
 77:   PetscScalar *testDerReal;            /* Workspace for derivative pushforward */
 78:   PetscReal   *x;                      /* Workspace for computing real coordinates */
 79:   PetscScalar *f0, *f1;                /* Point evaluations of weak form residual integrands */
 80:   PetscScalar *g0, *g1, *g2, *g3;      /* Point evaluations of weak form Jacobian integrands */
 81: };

 83: typedef struct {
 84:   PetscInt dummy; /* */
 85: } PetscDS_Basic;

 87: PETSC_INTERN PetscErrorCode PetscDSIsFE_Internal(PetscDS, PetscInt, PetscBool *);

 89: #endif