Actual source code: dagetarray.c

petsc-3.3-p7 2013-05-11
  1: 
  2: #include <petscdmda.h>    /*I   "petscdmda.h"   I*/

  6: /*@C
  7:    DMDAVecGetArray - Returns a multiple dimension array that shares data with 
  8:       the underlying vector and is indexed using the global dimensions.

 10:    Not Collective

 12:    Input Parameter:
 13: +  da - the distributed array
 14: -  vec - the vector, either a vector the same size as one obtained with 
 15:          DMCreateGlobalVector() or DMCreateLocalVector()
 16:    
 17:    Output Parameter:
 18: .  array - the array

 20:    Notes:
 21:     Call DMDAVecRestoreArray() once you have finished accessing the vector entries.

 23:     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!

 25:     If vec is a local vector (obtained with DMCreateLocalVector() etc) then they ghost point locations are accessable. If it is 
 26:     a global vector then the ghost points are not accessable. Of course with the local vector you will have had to do the 

 28:     appropriate DMLocalToGlobalBegin() and DMLocalToGlobalEnd() to have correct values in the ghost locations.

 30:   Fortran Notes: From Fortran use DMDAVecGetArrayF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate 
 31:        dimension. For a DMDA created with a dof of 1 use the dimension of the DMDA, for a DMDA created with a dof greater than 1 use one more than the 
 32:        dimension of the DMDA. The order of the indices is array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) (when dof is 1) otherwise
 33:        array(1:dof,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from 
 34:        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include finclude/petscdmda.h90 to access this routine.

 36:   Due to bugs in the compiler DMDAVecGetArrayF90() does not work with gfortran versions before 2.5

 38:   Level: intermediate

 40: .keywords: distributed array, get, corners, nodes, local indices, coordinates

 42: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF()
 43:           DMDAVecGetArrayDOF()
 44: @*/
 45: PetscErrorCode  DMDAVecGetArray(DM da,Vec vec,void *array)
 46: {
 48:   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;

 54:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
 55:   DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
 56:   DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);

 58:   /* Handle case where user passes in global vector as opposed to local */
 59:   VecGetLocalSize(vec,&N);
 60:   if (N == xm*ym*zm*dof) {
 61:     gxm = xm;
 62:     gym = ym;
 63:     gzm = zm;
 64:     gxs = xs;
 65:     gys = ys;
 66:     gzs = zs;
 67:   } else if (N != gxm*gym*gzm*dof) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMDA local sizes %D %D\n",N,xm*ym*zm*dof,gxm*gym*gzm*dof);

 69:   if (dim == 1) {
 70:     VecGetArray1d(vec,gxm*dof,gxs*dof,(PetscScalar **)array);
 71:   } else if (dim == 2) {
 72:     VecGetArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
 73:   } else if (dim == 3) {
 74:     VecGetArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
 75:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);

 77:   return(0);
 78: }

 82: /*@
 83:    DMDAVecRestoreArray - Restores a multiple dimension array obtained with DMDAVecGetArray()

 85:    Not Collective

 87:    Input Parameter:
 88: +  da - the distributed array
 89: .  vec - the vector, either a vector the same size as one obtained with 
 90:          DMCreateGlobalVector() or DMCreateLocalVector()
 91: -  array - the array

 93:   Level: intermediate

 95:   Fortran Notes: From Fortran use DMDAVecRestoreArrayF90()

 97: .keywords: distributed array, get, corners, nodes, local indices, coordinates

 99: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray()
100: @*/
101: PetscErrorCode  DMDAVecRestoreArray(DM da,Vec vec,void *array)
102: {
104:   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;

110:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
111:   DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
112:   DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);

114:   /* Handle case where user passes in global vector as opposed to local */
115:   VecGetLocalSize(vec,&N);
116:   if (N == xm*ym*zm*dof) {
117:     gxm = xm;
118:     gym = ym;
119:     gzm = zm;
120:     gxs = xs;
121:     gys = ys;
122:     gzs = zs;
123:   } else if (N != gxm*gym*gzm*dof) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMDA local sizes %D %D\n",N,xm*ym*zm*dof,gxm*gym*gzm*dof);

125:   if (dim == 1) {
126:     VecRestoreArray1d(vec,gxm*dof,gxs*dof,(PetscScalar **)array);
127:   } else if (dim == 2) {
128:     VecRestoreArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
129:   } else if (dim == 3) {
130:     VecRestoreArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
131:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
132:   return(0);
133: }

137: /*@C
138:    DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with 
139:       the underlying vector and is indexed using the global dimensions.

141:    Not Collective

143:    Input Parameter:
144: +  da - the distributed array
145: -  vec - the vector, either a vector the same size as one obtained with 
146:          DMCreateGlobalVector() or DMCreateLocalVector()
147:    
148:    Output Parameter:
149: .  array - the array

151:    Notes:
152:     Call DMDAVecRestoreArrayDOF() once you have finished accessing the vector entries.

154:     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!

156:   Level: intermediate

158: .keywords: distributed array, get, corners, nodes, local indices, coordinates

160: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF()
161: @*/
162: PetscErrorCode  DMDAVecGetArrayDOF(DM da,Vec vec,void *array)
163: {
165:   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;

168:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
169:   DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
170:   DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);

172:   /* Handle case where user passes in global vector as opposed to local */
173:   VecGetLocalSize(vec,&N);
174:   if (N == xm*ym*zm*dof) {
175:     gxm = xm;
176:     gym = ym;
177:     gzm = zm;
178:     gxs = xs;
179:     gys = ys;
180:     gzs = zs;
181:   } else if (N != gxm*gym*gzm*dof) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMDA local sizes %D %D\n",N,xm*ym*zm*dof,gxm*gym*gzm*dof);

183:   if (dim == 1) {
184:     VecGetArray2d(vec,gxm,dof,gxs,0,(PetscScalar ***)array);
185:   } else if (dim == 2) {
186:     VecGetArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);
187:   } else if (dim == 3) {
188:     VecGetArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);
189:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
190:   return(0);
191: }

195: /*@
196:    DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with DMDAVecGetArrayDOF()

198:    Not Collective

200:    Input Parameter:
201: +  da - the distributed array
202: .  vec - the vector, either a vector the same size as one obtained with 
203:          DMCreateGlobalVector() or DMCreateLocalVector()
204: -  array - the array

206:   Level: intermediate

208: .keywords: distributed array, get, corners, nodes, local indices, coordinates

210: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF()
211: @*/
212: PetscErrorCode  DMDAVecRestoreArrayDOF(DM da,Vec vec,void *array)
213: {
215:   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;

218:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
219:   DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
220:   DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);

222:   /* Handle case where user passes in global vector as opposed to local */
223:   VecGetLocalSize(vec,&N);
224:   if (N == xm*ym*zm*dof) {
225:     gxm = xm;
226:     gym = ym;
227:     gzm = zm;
228:     gxs = xs;
229:     gys = ys;
230:     gzs = zs;
231:   }

233:   if (dim == 1) {
234:     VecRestoreArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);
235:   } else if (dim == 2) {
236:     VecRestoreArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);
237:   } else if (dim == 3) {
238:     VecRestoreArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);
239:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
240:   return(0);
241: }