Actual source code: daimpl.h

  1: /* $Id: daimpl.h,v 1.44 2001/08/06 21:18:31 bsmith Exp $ */

  3: /*
  4:    Distributed arrays - communication tools for parallel, rectangular grids.
  5: */

  7: #if !defined(_DAIMPL_H)
  8: #define _DAIMPL_H

 10:  #include petscda.h

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

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

 31: typedef struct _DAOps *DAOps;
 32: struct _DAOps {
 33:   int  (*view)(DA,PetscViewer);
 34:   int  (*createglobalvector)(DA,Vec*);
 35:   int  (*getcoloring)(DA,ISColoringType,ISColoring*);
 36:   int  (*getmatrix)(DA,MatType,Mat*);
 37:   int  (*getinterpolation)(DA,DA,Mat*,Vec*);
 38:   int  (*refine)(DA,MPI_Comm,DA*);
 39: };

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

 58:   int                 nlocal,Nlocal;         /* local size of local vector and global vector */

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

 62:   ISLocalToGlobalMapping ltogmap,ltogmapb;   /* local to global mapping for associated vectors */
 63:   Vec                    coordinates;        /* coordinates (x,y,x) of local nodes, not including ghosts*/
 64:   char                   **fieldname;        /* names of individual components in vectors */

 66:   int                    *lx,*ly,*lz;        /* number of nodes in each partition block along 3 axis */
 67:   Vec                    natural;            /* global vector for storing items in natural order */
 68:   VecScatter             gton;               /* vector scatter from global to natural */

 70:   ISColoring            localcoloring;       /* set by DAGetColoring() */
 71:   ISColoring            ghostedcoloring;

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

 77:   int                    refine_x,refine_y,refine_z; /* ratio used in refining */

 79: #define DA_MAX_AD_ARRAYS 2 /* work arrays for holding derivative type data, via DAGetAdicArray() */
 80:   void                   *adarrayin[DA_MAX_AD_ARRAYS],*adarrayout[DA_MAX_AD_ARRAYS];
 81:   void                   *adarrayghostedin[DA_MAX_AD_ARRAYS],*adarrayghostedout[DA_MAX_AD_ARRAYS];
 82:   void                   *adstartin[DA_MAX_AD_ARRAYS],*adstartout[DA_MAX_AD_ARRAYS];
 83:   void                   *adstartghostedin[DA_MAX_AD_ARRAYS],*adstartghostedout[DA_MAX_AD_ARRAYS];
 84:   int                    tdof,ghostedtdof;

 86:                             /* work arrays for holding derivative type data, via DAGetAdicMFArray() */
 87:   void                   *admfarrayin[DA_MAX_AD_ARRAYS],*admfarrayout[DA_MAX_AD_ARRAYS];
 88:   void                   *admfarrayghostedin[DA_MAX_AD_ARRAYS],*admfarrayghostedout[DA_MAX_AD_ARRAYS];
 89:   void                   *admfstartin[DA_MAX_AD_ARRAYS],*admfstartout[DA_MAX_AD_ARRAYS];
 90:   void                   *admfstartghostedin[DA_MAX_AD_ARRAYS],*admfstartghostedout[DA_MAX_AD_ARRAYS];

 92: #define DA_MAX_WORK_ARRAYS 2 /* work arrays for holding work via DAGetArray() */
 93:   void                   *arrayin[DA_MAX_WORK_ARRAYS],*arrayout[DA_MAX_WORK_ARRAYS];
 94:   void                   *arrayghostedin[DA_MAX_WORK_ARRAYS],*arrayghostedout[DA_MAX_WORK_ARRAYS];
 95:   void                   *startin[DA_MAX_WORK_ARRAYS],*startout[DA_MAX_WORK_ARRAYS];
 96:   void                   *startghostedin[DA_MAX_WORK_ARRAYS],*startghostedout[DA_MAX_WORK_ARRAYS];

 98:   DALocalFunction1       lf;
 99:   DALocalFunction1       lj;
100:   DALocalFunction1       adic_lf;
101:   DALocalFunction1       adicmf_lf;
102:   DALocalFunction1       adifor_lf;
103:   DALocalFunction1       adiformf_lf;

105:   int                    (*lfi)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*);
106:   int                    (*adic_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*);
107:   int                    (*adicmf_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*);
108: };

110: /*
111:   Vectors:
112:      Global has on each processor the interior degrees of freedom and
113:          no ghost points. This vector is what the solvers usually see.
114:      Local has on each processor the ghost points as well. This is 
115:           what code to calculate Jacobians, etc. usually sees.
116:   Vector scatters:
117:      gtol - Global representation to local
118:      ltog - Local representation to global (involves no communication)
119:      ltol - Local representation to local representation, updates the
120:             ghostpoint values in the second vector from (correct) interior
121:             values in the first vector.  This is good for explicit
122:             nearest neighbor timestepping.
123: */

125: EXTERN int DAView_Binary(DA,PetscViewer);
126: EXTERN_C_BEGIN
127: EXTERN int VecView_MPI_DA(Vec,PetscViewer);
128: EXTERN int VecLoadIntoVector_Binary_DA(PetscViewer,Vec);
129: EXTERN_C_END

131: #endif