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 private/dmimpl.h
10: typedef struct _DAOps *DAOps;
11: struct _DAOps {
12: DMOPS(DA)
13: };
15: struct _p_DA {
16: PETSCHEADER(struct _DAOps);
17: DMHEADER
18: PetscInt M,N,P; /* array dimensions */
19: PetscInt m,n,p; /* processor layout */
20: PetscInt w; /* degrees of freedom per node */
21: PetscInt s; /* stencil width */
22: PetscInt xs,xe,ys,ye,zs,ze; /* range of local values */
23: PetscInt Xs,Xe,Ys,Ye,Zs,Ze; /* range including ghost values
24: values above already scaled by w */
25: PetscInt *idx,Nl; /* local to global map */
26: PetscInt base; /* global number of 1st local node */
27: DAPeriodicType wrap; /* indicates type of periodic boundaries */
28: VecScatter gtol,ltog,ltol; /* scatters, see below for details */
29: DAStencilType stencil_type; /* stencil, either box or star */
30: PetscInt dim; /* DA dimension (1,2, or 3) */
31: DAInterpolationType interptype;
33: PetscInt nlocal,Nlocal; /* local size of local vector and global vector */
35: AO ao; /* application ordering context */
37: ISLocalToGlobalMapping ltogmap,ltogmapb; /* local to global mapping for associated vectors */
38: Vec coordinates; /* coordinates (x,y,z) of local nodes, not including ghosts*/
39: DA da_coordinates; /* da for getting ghost values of coordinates */
40: Vec ghosted_coordinates;/* coordinates with ghost nodes */
41: char **fieldname; /* names of individual components in vectors */
43: PetscInt *lx,*ly,*lz; /* number of nodes in each partition block along 3 axis */
44: Vec natural; /* global vector for storing items in natural order */
45: VecScatter gton; /* vector scatter from global to natural */
46: PetscMPIInt *neighbors; /* ranks of all neighbors and self */
48: ISColoring localcoloring; /* set by DAGetColoring() */
49: ISColoring ghostedcoloring;
51: DAElementType elementtype;
52: PetscInt ne; /* number of elements */
53: PetscInt *e; /* the elements */
55: PetscInt refine_x,refine_y,refine_z; /* ratio used in refining */
57: #define DA_MAX_AD_ARRAYS 2 /* work arrays for holding derivative type data, via DAGetAdicArray() */
58: void *adarrayin[DA_MAX_AD_ARRAYS],*adarrayout[DA_MAX_AD_ARRAYS];
59: void *adarrayghostedin[DA_MAX_AD_ARRAYS],*adarrayghostedout[DA_MAX_AD_ARRAYS];
60: void *adstartin[DA_MAX_AD_ARRAYS],*adstartout[DA_MAX_AD_ARRAYS];
61: void *adstartghostedin[DA_MAX_AD_ARRAYS],*adstartghostedout[DA_MAX_AD_ARRAYS];
62: PetscInt tdof,ghostedtdof;
64: /* work arrays for holding derivative type data, via DAGetAdicMFArray() */
65: void *admfarrayin[DA_MAX_AD_ARRAYS],*admfarrayout[DA_MAX_AD_ARRAYS];
66: void *admfarrayghostedin[DA_MAX_AD_ARRAYS],*admfarrayghostedout[DA_MAX_AD_ARRAYS];
67: void *admfstartin[DA_MAX_AD_ARRAYS],*admfstartout[DA_MAX_AD_ARRAYS];
68: void *admfstartghostedin[DA_MAX_AD_ARRAYS],*admfstartghostedout[DA_MAX_AD_ARRAYS];
70: #define DA_MAX_WORK_ARRAYS 2 /* work arrays for holding work via DAGetArray() */
71: void *arrayin[DA_MAX_WORK_ARRAYS],*arrayout[DA_MAX_WORK_ARRAYS];
72: void *arrayghostedin[DA_MAX_WORK_ARRAYS],*arrayghostedout[DA_MAX_WORK_ARRAYS];
73: void *startin[DA_MAX_WORK_ARRAYS],*startout[DA_MAX_WORK_ARRAYS];
74: void *startghostedin[DA_MAX_WORK_ARRAYS],*startghostedout[DA_MAX_WORK_ARRAYS];
76: DALocalFunction1 lf;
77: DALocalFunction1 lj;
78: DALocalFunction1 adic_lf;
79: DALocalFunction1 adicmf_lf;
80: DALocalFunction1 adifor_lf;
81: DALocalFunction1 adiformf_lf;
83: PetscErrorCode (*lfi)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*);
84: PetscErrorCode (*adic_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*);
85: PetscErrorCode (*adicmf_lfi)(DALocalInfo*,MatStencil*,void*,void*,void*);
86: PetscErrorCode (*lfib)(DALocalInfo*,MatStencil*,void*,PetscScalar*,void*);
87: PetscErrorCode (*adic_lfib)(DALocalInfo*,MatStencil*,void*,void*,void*);
88: PetscErrorCode (*adicmf_lfib)(DALocalInfo*,MatStencil*,void*,void*,void*);
90: /* used by DASetBlockFills() */
91: PetscInt *ofill,*dfill;
93: /* used by DASetMatPreallocateOnly() */
94: PetscTruth prealloc_only;
95: };
97: /*
98: Vectors:
99: Global has on each processor the interior degrees of freedom and
100: no ghost points. This vector is what the solvers usually see.
101: Local has on each processor the ghost points as well. This is
102: what code to calculate Jacobians, etc. usually sees.
103: Vector scatters:
104: gtol - Global representation to local
105: ltog - Local representation to global (involves no communication)
106: ltol - Local representation to local representation, updates the
107: ghostpoint values in the second vector from (correct) interior
108: values in the first vector. This is good for explicit
109: nearest neighbor timestepping.
110: */
113: EXTERN PetscErrorCode VecView_MPI_DA(Vec,PetscViewer);
114: EXTERN PetscErrorCode VecLoadIntoVector_Binary_DA(PetscViewer,Vec);
116: EXTERN PetscErrorCode DAView_Private(DA);
120: #endif