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:   Level: intermediate

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

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

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

 40:   /* Handle case where user passes in global vector as opposed to local */
 41:   VecGetLocalSize(vec,&N);
 42:   if (N == xm*ym*zm*dof) {
 43:     gxm = xm;
 44:     gym = ym;
 45:     gzm = zm;
 46:     gxs = xs;
 47:     gys = ys;
 48:     gzs = zs;
 49:   } else if (N != gxm*gym*gzm*dof) {
 50:     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);
 51:   }

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

 63:   return(0);
 64: }

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

 71:    Not Collective

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

 79:   Level: intermediate

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

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

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

 95:   /* Handle case where user passes in global vector as opposed to local */
 96:   VecGetLocalSize(vec,&N);
 97:   if (N == xm*ym*zm*dof) {
 98:     gxm = xm;
 99:     gym = ym;
100:     gzm = zm;
101:     gxs = xs;
102:     gys = ys;
103:     gzs = zs;
104:   } else if (N != gxm*gym*gzm*dof) {
105:     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);
106:   }

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

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

126:    Not Collective

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

136:    Notes:
137:     Call DAVecRestoreArray() once you have finished accessing the vector entries.

139:   Level: intermediate

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

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

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

155:   /* Handle case where user passes in global vector as opposed to local */
156:   VecGetLocalSize(vec,&N);
157:   if (N == xm*ym*zm*dof) {
158:     gxm = xm;
159:     gym = ym;
160:     gzm = zm;
161:     gxs = xs;
162:     gys = ys;
163:     gzs = zs;
164:   } else if (N != gxm*gym*gzm*dof) {
165:     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);
166:   }

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

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

185:    Not Collective

187:    Input Parameter:
188: +  da - the distributed array
189: .  vec - the vector, either a vector the same size as one obtained with 
190:          DACreateGlobalVector() or DACreateLocalVector()
191: -  array - the array

193:   Level: intermediate

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

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

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

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

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