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,const MatType,Mat*);
23: int (*getinterpolation)(DM,DM,Mat*,Vec*);
24: int (*refine)(DM,MPI_Comm,DM*);
25: int (*getinjection)(DM,DM,VecScatter*);
26: };
28: struct _p_DM {
29: PETSCHEADER(struct _DMOps)
30: };
32: typedef struct _DAOps *DAOps;
33: struct _DAOps {
34: int (*view)(DA,PetscViewer);
35: int (*createglobalvector)(DA,Vec*);
36: int (*getcoloring)(DA,ISColoringType,ISColoring*);
37: int (*getmatrix)(DA,const MatType,Mat*);
38: int (*getinterpolation)(DA,DA,Mat*,Vec*);
39: int (*refine)(DA,MPI_Comm,DA*);
40: int (*getinjection)(DA,DA,VecScatter*);
41: };
43: struct _p_DA {
44: PETSCHEADER(struct _DAOps)
45: int M,N,P; /* array dimensions */
46: int m,n,p; /* processor layout */
47: int w; /* degrees of freedom per node */
48: int s; /* stencil width */
49: int xs,xe,ys,ye,zs,ze; /* range of local values */
50: int Xs,Xe,Ys,Ye,Zs,Ze; /* range including ghost values
51: values above already scaled by w */
52: int *idx,Nl; /* local to global map */
53: int 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: int dim; /* DA dimension (1,2, or 3) */
58: DAInterpolationType interptype;
60: int 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,x) of local nodes, not including ghosts*/
66: char **fieldname; /* names of individual components in vectors */
68: int *lx,*ly,*lz; /* number of nodes in each partition block along 3 axis */
69: Vec natural; /* global vector for storing items in natural order */
70: VecScatter gton; /* vector scatter from global to natural */
72: ISColoring localcoloring; /* set by DAGetColoring() */
73: ISColoring ghostedcoloring;
75: #define DA_MAX_WORK_VECTORS 10 /* work vectors available to users via DAVecGetArray() */
76: Vec localin[DA_MAX_WORK_VECTORS],localout[DA_MAX_WORK_VECTORS];
77: Vec globalin[DA_MAX_WORK_VECTORS],globalout[DA_MAX_WORK_VECTORS];
79: int refine_x,refine_y,refine_z; /* ratio used in refining */
81: #define DA_MAX_AD_ARRAYS 2 /* work arrays for holding derivative type data, via DAGetAdicArray() */
82: void *adarrayin[DA_MAX_AD_ARRAYS],*adarrayout[DA_MAX_AD_ARRAYS];
83: void *adarrayghostedin[DA_MAX_AD_ARRAYS],*adarrayghostedout[DA_MAX_AD_ARRAYS];
84: void *adstartin[DA_MAX_AD_ARRAYS],*adstartout[DA_MAX_AD_ARRAYS];
85: void *adstartghostedin[DA_MAX_AD_ARRAYS],*adstartghostedout[DA_MAX_AD_ARRAYS];
86: int tdof,ghostedtdof;
88: /* work arrays for holding derivative type data, via DAGetAdicMFArray() */
89: void *admfarrayin[DA_MAX_AD_ARRAYS],*admfarrayout[DA_MAX_AD_ARRAYS];
90: void *admfarrayghostedin[DA_MAX_AD_ARRAYS],*admfarrayghostedout[DA_MAX_AD_ARRAYS];
91: void *admfstartin[DA_MAX_AD_ARRAYS],*admfstartout[DA_MAX_AD_ARRAYS];
92: void *admfstartghostedin[DA_MAX_AD_ARRAYS],*admfstartghostedout[DA_MAX_AD_ARRAYS];
94: #define DA_MAX_WORK_ARRAYS 2 /* work arrays for holding work via DAGetArray() */
95: void *arrayin[DA_MAX_WORK_ARRAYS],*arrayout[DA_MAX_WORK_ARRAYS];
96: void *arrayghostedin[DA_MAX_WORK_ARRAYS],*arrayghostedout[DA_MAX_WORK_ARRAYS];
97: void *startin[DA_MAX_WORK_ARRAYS],*startout[DA_MAX_WORK_ARRAYS];
98: void *startghostedin[DA_MAX_WORK_ARRAYS],*startghostedout[DA_MAX_WORK_ARRAYS];
100: DALocalFunction1 lf;
101: DALocalFunction1 lj;
102: DALocalFunction1 adic_lf;
103: DALocalFunction1 adicmf_lf;
104: DALocalFunction1 adifor_lf;
105: DALocalFunction1 adiformf_lf;
107: int (*lfi)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*);
108: int (*adic_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*);
109: int (*adicmf_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*);
111: /* this is used by DASetBlockFills() */
112: int *ofill,*dfill;
113: };
115: /*
116: Vectors:
117: Global has on each processor the interior degrees of freedom and
118: no ghost points. This vector is what the solvers usually see.
119: Local has on each processor the ghost points as well. This is
120: what code to calculate Jacobians, etc. usually sees.
121: Vector scatters:
122: gtol - Global representation to local
123: ltog - Local representation to global (involves no communication)
124: ltol - Local representation to local representation, updates the
125: ghostpoint values in the second vector from (correct) interior
126: values in the first vector. This is good for explicit
127: nearest neighbor timestepping.
128: */
130: EXTERN_C_BEGIN
131: EXTERN int VecView_MPI_DA(Vec,PetscViewer);
132: EXTERN int VecLoadIntoVector_Binary_DA(PetscViewer,Vec);
133: EXTERN_C_END
135: #endif