Actual source code: discimpl.h

  1: /* $Id: discimpl.h,v 1.17 2000/07/16 23:20:00 knepley Exp $ */
  2: /*
  3:   This file includes the definition of structures used in PETSc for 
  4:   discretizations. This should not be included in users' code.
  5: */

  7: #ifndef __DISCIMPL_H

 10:  #include grid.h
 11:  #include discretization.h

 13: /* This structure defines the interface to registered operators. */
 14: struct _Operator {
 15:   Discretization      test;       /* This is the space of test functions */
 16:   OperatorFunction    opFunc;     /* Integrals of operators which depend on J */
 17:   ALEOperatorFunction ALEOpFunc;  /* Also incorporates an ALE velocity field */
 18:   PetscScalar        *precompInt; /* Integrals of operators on shape functions over the standard element */
 19: };

 21: /* For boundary discretizations, we need to call the assembly functions with edge information.
 22:    The current edge structure does not contain the midnode, so we pass it in a context wrapper
 23:    BoundaryContext. This will extend smoothly to 3D. */
 24: typedef struct _EdgeContext {
 25:   int   midnode; /* The midnode of an edge */
 26:   void *ctx;     /* User context */
 27: } EdgeContext;

 29: /* Discretization operations - These are the operations associated with all 
 30:     type of discretizations. Note that some of these often actually work on
 31:     the associated vectors and matrices.

 33:     The funcVal  array is just the value for every component of the function.
 34:     The fieldVal array is
 35:       f^0, df^0/dx^0, df^0/dx^1, ... , df^0/dx^d, f^1, df^1/dx^0, ... , f^c, ... , df^c/dx^d
 36: */
 37: struct _DiscretizationOps {
 38:       /* Generic Operations */
 39:   int (*setup)(Discretization),
 40:       (*setupoperators)(Discretization),
 41:       (*setfromoptions)(Discretization),
 42:       (*view)(Discretization, PetscViewer),
 43:       (*destroy)(Discretization),
 44:       /* Evaluation Operations */
 45:       (*evaluatefunctiongalerkin)(Discretization, Mesh, PointFunction, PetscScalar, int, PetscScalar *, void *),
 46:       (*evaluateoperatorgalerkin)(Discretization, Mesh, int, int, int, int, PetscScalar, int, PetscScalar *, PetscScalar *, void *),
 47:       (*evaluatealeoperatorgalerkin)(Discretization, Mesh, int, int, int, int, PetscScalar, int, PetscScalar *, PetscScalar *, PetscScalar *, void *),
 48:       (*evaluatenonlinearoperatorgalerkin)(Discretization, Mesh, NonlinearOperator, PetscScalar, int, int, PetscScalar **, PetscScalar *, void *),
 49:       (*evaluatenonlinearaleoperatorgalerkin)(Discretization, Mesh, NonlinearOperator, PetscScalar, int, int, PetscScalar **, PetscScalar *, PetscScalar *, void *),
 50:       /* Interpolation Operations */
 51:       (*interpolatefield)(Discretization, Mesh, int, double, double, double, PetscScalar *, PetscScalar *, InterpolationType),
 52:       (*interpolateelemvec)(Discretization, ElementVec, Discretization, ElementVec);
 53: };

 55: struct _Discretization {
 56:   PETSCHEADER(struct _DiscretizationOps)
 57:   void          *data;               /*          Discretization specific information and storage */
 58:   int            dim;                /* D:       The dimension */
 59:   int            funcs;              /* F:       The number of shape functions per element per component */
 60:   int            comp;               /* C:       The number of components in this field */
 61:   int            size;               /* S:       The number of shape functions per element = funcs*comp */
 62:   int            field;              /*          This is the field in the Grid using this discretization (or -1) */
 63:   /* Boundary integration */
 64:   Discretization bdDisc;             /*          The lower dimensional analog along the boundary */
 65:   /* Quadrature will be made a separate object */
 66:   int            numQuadPoints;      /* P:       Number of points used for Gaussian quadrature on the standard element */
 67:   double        *quadPoints;         /* [P*D]:   Quadrature points: x_0, y_0, x_1, y_1, etc. */
 68:   double        *quadWeights;        /* [P]:     Integration weights for each quadtrature point */
 69:   double        *quadShapeFuncs;     /* [P*F]:   Shape functions evaluated at the quadrature points */
 70:   double        *quadShapeFuncDers;  /* [P*F*D]: Shape function derivatives evaluated at the quadrature points */
 71:   /* Operators */
 72:   int            numOps;             /* O:       The number of registered operators. They are numbered sequentially */
 73:   int            maxOps;             /*          Maximum number of registered operators */
 74:   Operator      *operators;          /* [O]:     The registered operators */
 75:   /* Storage */
 76:   PetscScalar   *funcVal;            /* [C]:          Used for evaluation of function to avoid allocation during assembly */
 77:   PetscScalar  **fieldVal;           /* [2][C*(D+1)]: Used for evaluation of fields to avoid allocation during assembly */
 78: };

 80: #endif /* DISCIMPL_H */