Actual source code: dagetarray.c

  1: #define PETSCDM_DLL
  2: 
 3:  #include petscda.h

  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:     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!

 26:   Level: intermediate

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

 30: .seealso: DAGetGhostCorners(), DAGetCorners(), VecGetArray(), VecRestoreArray(), DAVecRestoreArray(), DAVecRestoreArrayDOF()
 31:           DAVecGetarrayDOF()
 32: @*/
 33: PetscErrorCode  DAVecGetArray(DA da,Vec vec,void *array)
 34: {
 36:   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;

 39:   DAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
 40:   DAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
 41:   DAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0);

 43:   /* Handle case where user passes in global vector as opposed to local */
 44:   VecGetLocalSize(vec,&N);
 45:   if (N == xm*ym*zm*dof) {
 46:     gxm = xm;
 47:     gym = ym;
 48:     gzm = zm;
 49:     gxs = xs;
 50:     gys = ys;
 51:     gzs = zs;
 52:   } else if (N != gxm*gym*gzm*dof) {
 53:     SETERRQ3(PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DA local sizes %D %D\n",N,xm*ym*zm*dof,gxm*gym*gzm*dof);
 54:   }

 56:   if (dim == 1) {
 57:     VecGetArray1d(vec,gxm*dof,gxs*dof,(PetscScalar **)array);
 58:   } else if (dim == 2) {
 59:     VecGetArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
 60:   } else if (dim == 3) {
 61:     VecGetArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
 62:   } else {
 63:     SETERRQ1(PETSC_ERR_ARG_CORRUPT,"DA dimension not 1, 2, or 3, it is %D\n",dim);
 64:   }

 66:   return(0);
 67: }

 71: /*@
 72:    DAVecRestoreArray - Restores a multiple dimension array obtained with DAVecGetArray()

 74:    Not Collective

 76:    Input Parameter:
 77: +  da - the distributed array
 78: .  vec - the vector, either a vector the same size as one obtained with 
 79:          DACreateGlobalVector() or DACreateLocalVector()
 80: -  array - the array

 82:   Level: intermediate

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

 86: .seealso: DAGetGhostCorners(), DAGetCorners(), VecGetArray(), VecRestoreArray(), DAVecGetArray()
 87: @*/
 88: PetscErrorCode  DAVecRestoreArray(DA da,Vec vec,void *array)
 89: {
 91:   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;

 94:   DAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
 95:   DAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
 96:   DAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0);

 98:   /* Handle case where user passes in global vector as opposed to local */
 99:   VecGetLocalSize(vec,&N);
100:   if (N == xm*ym*zm*dof) {
101:     gxm = xm;
102:     gym = ym;
103:     gzm = zm;
104:     gxs = xs;
105:     gys = ys;
106:     gzs = zs;
107:   } else if (N != gxm*gym*gzm*dof) {
108:     SETERRQ3(PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DA local sizes %D %D\n",N,xm*ym*zm*dof,gxm*gym*gzm*dof);
109:   }

111:   if (dim == 1) {
112:     VecRestoreArray1d(vec,gxm*dof,gxs*dof,(PetscScalar **)array);
113:   } else if (dim == 2) {
114:     VecRestoreArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
115:   } else if (dim == 3) {
116:     VecRestoreArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
117:   } else {
118:     SETERRQ1(PETSC_ERR_ARG_CORRUPT,"DA dimension not 1, 2, or 3, it is %D\n",dim);
119:   }
120:   return(0);
121: }

125: /*@C
126:    DAVecGetArrayDOF - Returns a multiple dimension array that shares data with 
127:       the underlying vector and is indexed using the global dimensions.

129:    Not Collective

131:    Input Parameter:
132: +  da - the distributed array
133: -  vec - the vector, either a vector the same size as one obtained with 
134:          DACreateGlobalVector() or DACreateLocalVector()
135:    
136:    Output Parameter:
137: .  array - the array

139:    Notes:
140:     Call DAVecRestoreArray() once you have finished accessing the vector entries.

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

144:   Level: intermediate

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

148: .seealso: DAGetGhostCorners(), DAGetCorners(), VecGetArray(), VecRestoreArray(), DAVecRestoreArray(), DAVecGetArray(), DAVecRestoreArrayDOF()
149: @*/
150: PetscErrorCode  DAVecGetArrayDOF(DA da,Vec vec,void *array)
151: {
153:   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;

156:   DAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
157:   DAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
158:   DAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0);

160:   /* Handle case where user passes in global vector as opposed to local */
161:   VecGetLocalSize(vec,&N);
162:   if (N == xm*ym*zm*dof) {
163:     gxm = xm;
164:     gym = ym;
165:     gzm = zm;
166:     gxs = xs;
167:     gys = ys;
168:     gzs = zs;
169:   } else if (N != gxm*gym*gzm*dof) {
170:     SETERRQ3(PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DA local sizes %D %D\n",N,xm*ym*zm*dof,gxm*gym*gzm*dof);
171:   }

173:   if (dim == 1) {
174:     VecGetArray2d(vec,gxm,dof,gxs,0,(PetscScalar ***)array);
175:   } else if (dim == 2) {
176:     VecGetArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);
177:   } else if (dim == 3) {
178:     VecGetArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);
179:   } else {
180:     SETERRQ1(PETSC_ERR_ARG_CORRUPT,"DA dimension not 1, 2, or 3, it is %D\n",dim);
181:   }
182:   return(0);
183: }

187: /*@
188:    DAVecRestoreArrayDOF - Restores a multiple dimension array obtained with DAVecGetArrayDOF()

190:    Not Collective

192:    Input Parameter:
193: +  da - the distributed array
194: .  vec - the vector, either a vector the same size as one obtained with 
195:          DACreateGlobalVector() or DACreateLocalVector()
196: -  array - the array

198:   Level: intermediate

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

202: .seealso: DAGetGhostCorners(), DAGetCorners(), VecGetArray(), VecRestoreArray(), DAVecGetArray(), DAVecGetArrayDOF(), DAVecRestoreArrayDOF()
203: @*/
204: PetscErrorCode  DAVecRestoreArrayDOF(DA da,Vec vec,void *array)
205: {
207:   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;

210:   DAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
211:   DAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
212:   DAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0);

214:   /* Handle case where user passes in global vector as opposed to local */
215:   VecGetLocalSize(vec,&N);
216:   if (N == xm*ym*zm*dof) {
217:     gxm = xm;
218:     gym = ym;
219:     gzm = zm;
220:     gxs = xs;
221:     gys = ys;
222:     gzs = zs;
223:   }

225:   if (dim == 1) {
226:     VecRestoreArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);
227:   } else if (dim == 2) {
228:     VecRestoreArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);
229:   } else if (dim == 3) {
230:     VecRestoreArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);
231:   } else {
232:     SETERRQ1(PETSC_ERR_ARG_CORRUPT,"DA dimension not 1, 2, or 3, it is %D\n",dim);
233:   }
234:   return(0);
235: }