Actual source code: dagetarray.c
1: /*$Id: dagetarray.c,v 1.13 2001/08/06 21:18:33 bsmith Exp $*/
2:
3: #include petscda.h
5: #undef __FUNCT__
7: /*@C
8: DAVecGetArray - Returns a multiple dimension array that shares data with
9: the underlying vector and is indexed using the global dimensions.
11: Not Collective
13: Input Parameter:
14: + da - the distributed array
15: - vec - the vector, either a vector the same size as one obtained with
16: DACreateGlobalVector() or DACreateLocalVector()
17:
18: Output Parameter:
19: . array - the array
21: Notes:
22: Call DAVecRestoreArray() once you have finished accessing the vector entries.
24: Level: intermediate
26: .keywords: distributed array, get, corners, nodes, local indices, coordinates
28: .seealso: DAGetGhostCorners(), DAGetCorners(), VecGetArray(), VecRestoreArray(), DAVecRestoreArray()
29: @*/
30: int DAVecGetArray(DA da,Vec vec,void **array)
31: {
32: int ierr,xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
35: DAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
36: DAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
37: DAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0);
39: /* Handle case where user passes in global vector as opposed to local */
40: VecGetLocalSize(vec,&N);
41: if (N == xm*ym*zm*dof) {
42: gxm = xm;
43: gym = ym;
44: gzm = zm;
45: gxs = xs;
46: gys = ys;
47: gzs = zs;
48: } else if (N != gxm*gym*gzm*dof) {
49: SETERRQ3(1,"Vector local size %d is not compatible with DA local sizes %d %dn",N,xm*ym*zm*dof,gxm*gym*gzm*dof);
50: }
52: if (dim == 1) {
53: VecGetArray1d(vec,gxm*dof,gxs*dof,(PetscScalar **)array);
54: } else if (dim == 2) {
55: VecGetArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
56: } else if (dim == 3) {
57: VecGetArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
58: } else {
59: SETERRQ1(1,"DA dimension not 1, 2, or 3, it is %dn",dim);
60: }
62: return(0);
63: }
65: #undef __FUNCT__
67: /*@
68: DAVecRestoreArray - Restores a multiple dimension array obtained with DAVecGetArray()
70: Not Collective
72: Input Parameter:
73: + da - the distributed array
74: . vec - the vector, either a vector the same size as one obtained with
75: DACreateGlobalVector() or DACreateLocalVector()
76: - array - the array
78: Level: intermediate
80: .keywords: distributed array, get, corners, nodes, local indices, coordinates
82: .seealso: DAGetGhostCorners(), DAGetCorners(), VecGetArray(), VecRestoreArray(), DAVecGetArray()
83: @*/
84: int DAVecRestoreArray(DA da,Vec vec,void **array)
85: {
86: int ierr,xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
89: DAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
90: DAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
91: DAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0);
93: /* Handle case where user passes in global vector as opposed to local */
94: VecGetLocalSize(vec,&N);
95: if (N == xm*ym*zm*dof) {
96: gxm = xm;
97: gym = ym;
98: gzm = zm;
99: gxs = xs;
100: gys = ys;
101: gzs = zs;
102: }
104: if (dim == 1) {
105: VecRestoreArray1d(vec,gxm*dof,gxs*dof,(PetscScalar **)array);
106: } else if (dim == 2) {
107: VecRestoreArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
108: } else if (dim == 3) {
109: VecRestoreArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
110: } else {
111: SETERRQ1(1,"DA dimension not 1, 2, or 3, it is %dn",dim);
112: }
113: return(0);
114: }