Actual source code: daimpl.h

  1: /*
  2:    Distributed arrays - communication tools for parallel, rectangular grids.
  3: */

  5: #if !defined(_DAIMPL_H)
  6: #define _DAIMPL_H

 8:  #include petscda.h

 10: /*
 11:    The DM interface is shared by DA, VecPack, and any other object that may 
 12:   be used with the DMMG class. If you change this MAKE SURE you change
 13:   struct _DAOps and struct _VecPackOps!
 14: */
 15: typedef struct _DMOps *DMOps;
 16: struct _DMOps {
 17:   PetscErrorCode (*view)(DM,PetscViewer);
 18:   PetscErrorCode (*createglobalvector)(DM,Vec*);
 19:   PetscErrorCode (*getcoloring)(DM,ISColoringType,ISColoring*);
 20:   PetscErrorCode (*getmatrix)(DM, MatType,Mat*);
 21:   PetscErrorCode (*getinterpolation)(DM,DM,Mat*,Vec*);
 22:   PetscErrorCode (*refine)(DM,MPI_Comm,DM*);
 23:   PetscErrorCode (*getinjection)(DM,DM,VecScatter*);
 24: };

 26: struct _p_DM {
 27:   PETSCHEADER(struct _DMOps);
 28: };

 30: typedef struct _DAOps *DAOps;
 31: struct _DAOps {
 32:   PetscErrorCode (*view)(DA,PetscViewer);
 33:   PetscErrorCode (*createglobalvector)(DA,Vec*);
 34:   PetscErrorCode (*getcoloring)(DA,ISColoringType,ISColoring*);
 35:   PetscErrorCode (*getmatrix)(DA, MatType,Mat*);
 36:   PetscErrorCode (*getinterpolation)(DA,DA,Mat*,Vec*);
 37:   PetscErrorCode (*refine)(DA,MPI_Comm,DA*);
 38:   PetscErrorCode (*getinjection)(DA,DA,VecScatter*);
 39:   PetscErrorCode (*getelements)(DA,PetscInt*,const PetscInt*[]);
 40:   PetscErrorCode (*restoreelements)(DA,PetscInt*,const PetscInt*[]);
 41: };

 43: struct _p_DA {
 44:   PETSCHEADER(struct _DAOps);
 45:   PetscInt            M,N,P;                 /* array dimensions */
 46:   PetscInt            m,n,p;                 /* processor layout */
 47:   PetscInt            w;                     /* degrees of freedom per node */
 48:   PetscInt            s;                     /* stencil width */
 49:   PetscInt            xs,xe,ys,ye,zs,ze;     /* range of local values */
 50:   PetscInt            Xs,Xe,Ys,Ye,Zs,Ze;     /* range including ghost values
 51:                                                    values above already scaled by w */
 52:   PetscInt            *idx,Nl;               /* local to global map */
 53:   PetscInt            base;                  /* global number of 1st local node */
 54:   DAPeriodicType      wrap;                  /* indicates type of periodic boundaries */
 55:   VecScatter          gtol,ltog,ltol;        /* scatters, see below for details */
 56:   DAStencilType       stencil_type;          /* stencil, either box or star */
 57:   PetscInt            dim;                   /* DA dimension (1,2, or 3) */
 58:   DAInterpolationType interptype;

 60:   PetscInt            nlocal,Nlocal;         /* local size of local vector and global vector */

 62:   AO                  ao;                    /* application ordering context */

 64:   ISLocalToGlobalMapping ltogmap,ltogmapb;   /* local to global mapping for associated vectors */
 65:   Vec                    coordinates;        /* coordinates (x,y,z) of local nodes, not including ghosts*/
 66:   DA                     da_coordinates;     /* da for getting ghost values of coordinates */
 67:   Vec                    ghosted_coordinates;/* coordinates with ghost nodes */
 68:   char                   **fieldname;        /* names of individual components in vectors */

 70:   PetscInt               *lx,*ly,*lz;        /* number of nodes in each partition block along 3 axis */
 71:   Vec                    natural;            /* global vector for storing items in natural order */
 72:   VecScatter             gton;               /* vector scatter from global to natural */

 74:   ISColoring             localcoloring;       /* set by DAGetColoring() */
 75:   ISColoring             ghostedcoloring;

 77:   DAElementType          elementtype;
 78:   PetscInt               ne;                  /* number of elements */
 79:   PetscInt               *e;                  /* the elements */

 81: #define DA_MAX_WORK_VECTORS 10 /* work vectors available to users  via DAVecGetArray() */
 82:   Vec                    localin[DA_MAX_WORK_VECTORS],localout[DA_MAX_WORK_VECTORS];
 83:   Vec                    globalin[DA_MAX_WORK_VECTORS],globalout[DA_MAX_WORK_VECTORS];

 85:   PetscInt               refine_x,refine_y,refine_z; /* ratio used in refining */

 87: #define DA_MAX_AD_ARRAYS 2 /* work arrays for holding derivative type data, via DAGetAdicArray() */
 88:   void                   *adarrayin[DA_MAX_AD_ARRAYS],*adarrayout[DA_MAX_AD_ARRAYS];
 89:   void                   *adarrayghostedin[DA_MAX_AD_ARRAYS],*adarrayghostedout[DA_MAX_AD_ARRAYS];
 90:   void                   *adstartin[DA_MAX_AD_ARRAYS],*adstartout[DA_MAX_AD_ARRAYS];
 91:   void                   *adstartghostedin[DA_MAX_AD_ARRAYS],*adstartghostedout[DA_MAX_AD_ARRAYS];
 92:   PetscInt                    tdof,ghostedtdof;

 94:                             /* work arrays for holding derivative type data, via DAGetAdicMFArray() */
 95:   void                   *admfarrayin[DA_MAX_AD_ARRAYS],*admfarrayout[DA_MAX_AD_ARRAYS];
 96:   void                   *admfarrayghostedin[DA_MAX_AD_ARRAYS],*admfarrayghostedout[DA_MAX_AD_ARRAYS];
 97:   void                   *admfstartin[DA_MAX_AD_ARRAYS],*admfstartout[DA_MAX_AD_ARRAYS];
 98:   void                   *admfstartghostedin[DA_MAX_AD_ARRAYS],*admfstartghostedout[DA_MAX_AD_ARRAYS];

100: #define DA_MAX_WORK_ARRAYS 2 /* work arrays for holding work via DAGetArray() */
101:   void                   *arrayin[DA_MAX_WORK_ARRAYS],*arrayout[DA_MAX_WORK_ARRAYS];
102:   void                   *arrayghostedin[DA_MAX_WORK_ARRAYS],*arrayghostedout[DA_MAX_WORK_ARRAYS];
103:   void                   *startin[DA_MAX_WORK_ARRAYS],*startout[DA_MAX_WORK_ARRAYS];
104:   void                   *startghostedin[DA_MAX_WORK_ARRAYS],*startghostedout[DA_MAX_WORK_ARRAYS];

106:   DALocalFunction1       lf;
107:   DALocalFunction1       lj;
108:   DALocalFunction1       adic_lf;
109:   DALocalFunction1       adicmf_lf;
110:   DALocalFunction1       adifor_lf;
111:   DALocalFunction1       adiformf_lf;

113:   PetscErrorCode (*lfi)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*);
114:   PetscErrorCode (*adic_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*);
115:   PetscErrorCode (*adicmf_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*);
116:   PetscErrorCode (*lfib)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*);
117:   PetscErrorCode (*adic_lfib)(DALocalInfo*,MatStencil*,void*,void*,void*);
118:   PetscErrorCode (*adicmf_lfib)(DALocalInfo*,MatStencil*,void*,void*,void*);

120:   /* used by DASetBlockFills() */
121:   PetscInt               *ofill,*dfill;

123:   /* used by DASetMatPreallocateOnly() */
124:   PetscTruth             prealloc_only;
125: };

127: /*
128:   Vectors:
129:      Global has on each processor the interior degrees of freedom and
130:          no ghost points. This vector is what the solvers usually see.
131:      Local has on each processor the ghost points as well. This is 
132:           what code to calculate Jacobians, etc. usually sees.
133:   Vector scatters:
134:      gtol - Global representation to local
135:      ltog - Local representation to global (involves no communication)
136:      ltol - Local representation to local representation, updates the
137:             ghostpoint values in the second vector from (correct) interior
138:             values in the first vector.  This is good for explicit
139:             nearest neighbor timestepping.
140: */

143: EXTERN PetscErrorCode PETSCDM_DLLEXPORT VecView_MPI_DA(Vec,PetscViewer);
144: EXTERN PetscErrorCode PETSCDM_DLLEXPORT VecLoadIntoVector_Binary_DA(PetscViewer,Vec);


149: #endif