Actual source code: dagetarray.c

petsc-3.4.5 2014-06-29
  2: #include <petsc-private/dmdaimpl.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:    Collective on Vec

 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()

 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(0:dof-1,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 4.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:   if (da->defaultSection) {
 55:     VecGetArray(vec,(PetscScalar**)array);
 56:     return(0);
 57:   }
 58:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
 59:   DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
 60:   DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);

 62:   /* Handle case where user passes in global vector as opposed to local */
 63:   VecGetLocalSize(vec,&N);
 64:   if (N == xm*ym*zm*dof) {
 65:     gxm = xm;
 66:     gym = ym;
 67:     gzm = zm;
 68:     gxs = xs;
 69:     gys = ys;
 70:     gzs = zs;
 71:   } 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);

 73:   if (dim == 1) {
 74:     VecGetArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);
 75:   } else if (dim == 2) {
 76:     VecGetArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
 77:   } else if (dim == 3) {
 78:     VecGetArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
 79:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
 80:   return(0);
 81: }

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

 88:    Collective on Vec

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

 96:   Level: intermediate

 98:   Fortran Notes: From Fortran use DMDAVecRestoreArrayF90()

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

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

113:   if (da->defaultSection) {
114:     VecRestoreArray(vec,(PetscScalar**)array);
115:     return(0);
116:   }
117:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
118:   DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
119:   DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);

121:   /* Handle case where user passes in global vector as opposed to local */
122:   VecGetLocalSize(vec,&N);
123:   if (N == xm*ym*zm*dof) {
124:     gxm = xm;
125:     gym = ym;
126:     gzm = zm;
127:     gxs = xs;
128:     gys = ys;
129:     gzs = zs;
130:   } 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);

132:   if (dim == 1) {
133:     VecRestoreArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);
134:   } else if (dim == 2) {
135:     VecRestoreArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
136:   } else if (dim == 3) {
137:     VecRestoreArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
138:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
139:   return(0);
140: }

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

148:    Not Collective

150:    Input Parameter:
151: +  da - the distributed array
152: -  vec - the vector, either a vector the same size as one obtained with
153:          DMCreateGlobalVector() or DMCreateLocalVector()

155:    Output Parameter:
156: .  array - the array

158:    Notes:
159:     Call DMDAVecRestoreArrayDOF() once you have finished accessing the vector entries.

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

163:     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayF90() and declare your array with one higher dimension,
164:     see src/dm/examples/tutorials/ex11f90.F

166:   Level: intermediate

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

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

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

182:   /* Handle case where user passes in global vector as opposed to local */
183:   VecGetLocalSize(vec,&N);
184:   if (N == xm*ym*zm*dof) {
185:     gxm = xm;
186:     gym = ym;
187:     gzm = zm;
188:     gxs = xs;
189:     gys = ys;
190:     gzs = zs;
191:   } 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);

193:   if (dim == 1) {
194:     VecGetArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);
195:   } else if (dim == 2) {
196:     VecGetArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);
197:   } else if (dim == 3) {
198:     VecGetArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);
199:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
200:   return(0);
201: }

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

208:    Not Collective

210:    Input Parameter:
211: +  da - the distributed array
212: .  vec - the vector, either a vector the same size as one obtained with
213:          DMCreateGlobalVector() or DMCreateLocalVector()
214: -  array - the array

216:   Level: intermediate

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

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

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

232:   /* Handle case where user passes in global vector as opposed to local */
233:   VecGetLocalSize(vec,&N);
234:   if (N == xm*ym*zm*dof) {
235:     gxm = xm;
236:     gym = ym;
237:     gzm = zm;
238:     gxs = xs;
239:     gys = ys;
240:     gzs = zs;
241:   }

243:   if (dim == 1) {
244:     VecRestoreArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);
245:   } else if (dim == 2) {
246:     VecRestoreArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);
247:   } else if (dim == 3) {
248:     VecRestoreArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);
249:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
250:   return(0);
251: }