Actual source code: petscdsimpl.h

petsc-3.7.6 2017-04-24
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 _PetscDSOps *PetscDSOps;
 11: struct _PetscDSOps {
 12:   PetscErrorCode (*setfromoptions)(PetscDS);
 13:   PetscErrorCode (*setup)(PetscDS);
 14:   PetscErrorCode (*view)(PetscDS,PetscViewer);
 15:   PetscErrorCode (*destroy)(PetscDS);
 16: };

 18: struct _p_PetscDS {
 19:   PETSCHEADER(struct _PetscDSOps);
 20:   void        *data;      /* Implementation object */
 21:   PetscBool    setup;     /* Flag for setup */
 22:   PetscInt     Nf;        /* The number of solution fields */
 23:   PetscBool   *implicit;  /* Flag for implicit or explicit solve */
 24:   PetscBool   *adjacency; /* Flag for variable influence */
 25:   PetscObject *disc;      /* The discretization for each solution field (PetscFE, PetscFV, etc.) */
 26:   PetscObject *discBd;    /* The boundary discretization for each solution field (PetscFE, PetscFV, etc.) */
 27:   PetscPointFunc   *obj;  /* Scalar integral (like an objective function) */
 28:   PetscPointFunc   *f;    /* Weak form integrands for F, f_0, f_1 */
 29:   PetscPointJac    *g;    /* Weak form integrands for J = dF/du, g_0, g_1, g_2, g_3 */
 30:   PetscPointJac    *gp;   /* Weak form integrands for preconditioner for J, g_0, g_1, g_2, g_3 */
 31:   PetscPointJac    *gt;   /* Weak form integrands for dF/du_t, g_0, g_1, g_2, g_3 */
 32:   PetscBdPointFunc *fBd;  /* Weak form boundary integrands F_bd, f_0, f_1 */
 33:   PetscBdPointJac  *gBd;  /* Weak form boundary integrands J_bd = dF_bd/du, g_0, g_1, g_2, g_3 */
 34:   PetscRiemannFunc *r;    /* Riemann solvers */
 35:   void       **ctx;       /* User contexts for each field */
 36:   PetscInt     dim;       /* The spatial dimension */
 37:   /* Computed sizes */
 38:   PetscInt     totDim, totDimBd;       /* Total system dimension */
 39:   PetscInt     totComp;                /* Total field components */
 40:   /* Work space */
 41:   PetscInt    *off,       *offBd;      /* Offsets for each field */
 42:   PetscInt    *offDer,    *offDerBd;   /* Derivative offsets for each field */
 43:   PetscReal  **basis,    **basisBd;    /* Default basis tabulation for each field */
 44:   PetscReal  **basisDer, **basisDerBd; /* Default basis derivative tabulation for each field */
 45:   PetscScalar *u;                      /* Field evaluation */
 46:   PetscScalar *u_t;                    /* Field time derivative evaluation */
 47:   PetscScalar *u_x;                    /* Field gradient evaluation */
 48:   PetscScalar *refSpaceDer;            /* Workspace for computing derivative in the reference coordinates */
 49:   PetscReal   *x;                      /* Workspace for computing real coordinates */
 50:   PetscScalar *f0, *f1;                /* Point evaluations of weak form residual integrands */
 51:   PetscScalar *g0, *g1, *g2, *g3;      /* Point evaluations of weak form Jacobian integrands */
 52: };

 54: typedef struct {
 55:   PetscInt dummy; /* */
 56: } PetscDS_Basic;

 58: #endif