Actual source code: dagetarray.c

petsc-3.9.4 2018-09-11
Report Typos and Errors

  2:  #include <petsc/private/dmdaimpl.h>

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

  8:    Logically collective on Vec

 10:    Input Parameter:
 11: +  da - the distributed array
 12: -  vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()

 14:    Output Parameter:
 15: .  array - the array

 17:    Notes:
 18:     Call DMDAVecRestoreArray() once you have finished accessing the vector entries.

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

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

 25:     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.

 27:   Fortran Notes: From Fortran use DMDAVecGetArrayF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
 28:        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
 29:        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
 30:        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
 31:        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.

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

 35:   Level: intermediate

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

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

 51:   if (da->defaultSection) {
 52:     VecGetArray(vec,(PetscScalar**)array);
 53:     return(0);
 54:   }
 55:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
 56:   DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
 57:   DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);

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

 70:   if (dim == 1) {
 71:     VecGetArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);
 72:   } else if (dim == 2) {
 73:     VecGetArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
 74:   } else if (dim == 3) {
 75:     VecGetArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
 76:   } 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: }

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

 83:    Logically collective on Vec

 85:    Input Parameter:
 86: +  da - the distributed array
 87: .  vec - the vector, either a vector the same size as one obtained with
 88:          DMCreateGlobalVector() or DMCreateLocalVector()
 89: -  array - the array, non-NULL pointer is zeroed

 91:   Level: intermediate

 93:   Fortran Notes: From Fortran use DMDAVecRestoreArrayF90()

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

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

108:   if (da->defaultSection) {
109:     VecRestoreArray(vec,(PetscScalar**)array);
110:     return(0);
111:   }
112:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
113:   DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
114:   DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);

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

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

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

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:     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayF90() and declare your array with one higher dimension,
157:     see src/dm/examples/tutorials/ex11f90.F

159:   Level: intermediate

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

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

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

175:   /* Handle case where user passes in global vector as opposed to local */
176:   VecGetLocalSize(vec,&N);
177:   if (N == xm*ym*zm*dof) {
178:     gxm = xm;
179:     gym = ym;
180:     gzm = zm;
181:     gxs = xs;
182:     gys = ys;
183:     gzs = zs;
184:   } 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);

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

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

199:    Logically collective

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

207:   Level: intermediate

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

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

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

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

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

244: /*@C
245:    DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
246:       the underlying vector and is indexed using the global dimensions.

248:    Not collective

250:    Input Parameter:
251: +  da - the distributed array
252: -  vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()

254:    Output Parameter:
255: .  array - the array

257:    Notes:
258:     Call DMDAVecRestoreArrayRead() once you have finished accessing the vector entries.

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

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

265:     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.

267:   Fortran Notes: From Fortran use DMDAVecGetArrayReadF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
268:        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
269:        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
270:        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
271:        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.

273:   Due to bugs in the compiler DMDAVecGetArrayReadF90() does not work with gfortran versions before 4.5

275:   Level: intermediate

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

279: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF()
280:           DMDAVecGetArrayDOF()
281: @*/
282: PetscErrorCode  DMDAVecGetArrayRead(DM da,Vec vec,void *array)
283: {
285:   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;

291:   if (da->defaultSection) {
292:     VecGetArrayRead(vec,(const PetscScalar**)array);
293:     return(0);
294:   }
295:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
296:   DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
297:   DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);

299:   /* Handle case where user passes in global vector as opposed to local */
300:   VecGetLocalSize(vec,&N);
301:   if (N == xm*ym*zm*dof) {
302:     gxm = xm;
303:     gym = ym;
304:     gzm = zm;
305:     gxs = xs;
306:     gys = ys;
307:     gzs = zs;
308:   } 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);

310:   if (dim == 1) {
311:     VecGetArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);
312:   } else if (dim == 2) {
313:     VecGetArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
314:   } else if (dim == 3) {
315:     VecGetArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
316:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
317:   return(0);
318: }

320: /*@
321:    DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with DMDAVecGetArrayRead()

323:    Not collective

325:    Input Parameter:
326: +  da - the distributed array
327: .  vec - the vector, either a vector the same size as one obtained with
328:          DMCreateGlobalVector() or DMCreateLocalVector()
329: -  array - the array, non-NULL pointer is zeroed

331:   Level: intermediate

333:   Fortran Notes: From Fortran use DMDAVecRestoreArrayReadF90()

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

337: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray()
338: @*/
339: PetscErrorCode  DMDAVecRestoreArrayRead(DM da,Vec vec,void *array)
340: {
342:   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;

348:   if (da->defaultSection) {
349:     VecRestoreArrayRead(vec,(const PetscScalar**)array);
350:     return(0);
351:   }
352:   DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
353:   DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
354:   DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);

356:   /* Handle case where user passes in global vector as opposed to local */
357:   VecGetLocalSize(vec,&N);
358:   if (N == xm*ym*zm*dof) {
359:     gxm = xm;
360:     gym = ym;
361:     gzm = zm;
362:     gxs = xs;
363:     gys = ys;
364:     gzs = zs;
365:   } 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);

367:   if (dim == 1) {
368:     VecRestoreArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);
369:   } else if (dim == 2) {
370:     VecRestoreArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
371:   } else if (dim == 3) {
372:     VecRestoreArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
373:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
374:   return(0);
375: }

377: /*@C
378:    DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
379:       the underlying vector and is indexed using the global dimensions.

381:    Not Collective

383:    Input Parameter:
384: +  da - the distributed array
385: -  vec - the vector, either a vector the same size as one obtained with
386:          DMCreateGlobalVector() or DMCreateLocalVector()

388:    Output Parameter:
389: .  array - the array

391:    Notes:
392:     Call DMDAVecRestoreArrayDOFRead() once you have finished accessing the vector entries.

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

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

399:   Level: intermediate

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

403: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF()
404: @*/
405: PetscErrorCode  DMDAVecGetArrayDOFRead(DM da,Vec vec,void *array)
406: {
408:   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;

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

415:   /* Handle case where user passes in global vector as opposed to local */
416:   VecGetLocalSize(vec,&N);
417:   if (N == xm*ym*zm*dof) {
418:     gxm = xm;
419:     gym = ym;
420:     gzm = zm;
421:     gxs = xs;
422:     gys = ys;
423:     gzs = zs;
424:   } 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);

426:   if (dim == 1) {
427:     VecGetArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);
428:   } else if (dim == 2) {
429:     VecGetArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);
430:   } else if (dim == 3) {
431:     VecGetArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);
432:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
433:   return(0);
434: }

436: /*@
437:    DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with DMDAVecGetArrayDOFRead()

439:    Not Collective

441:    Input Parameter:
442: +  da - the distributed array
443: .  vec - the vector, either a vector the same size as one obtained with
444:          DMCreateGlobalVector() or DMCreateLocalVector()
445: -  array - the array

447:   Level: intermediate

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

451: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF()
452: @*/
453: PetscErrorCode  DMDAVecRestoreArrayDOFRead(DM da,Vec vec,void *array)
454: {
456:   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;

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

463:   /* Handle case where user passes in global vector as opposed to local */
464:   VecGetLocalSize(vec,&N);
465:   if (N == xm*ym*zm*dof) {
466:     gxm = xm;
467:     gym = ym;
468:     gzm = zm;
469:     gxs = xs;
470:     gys = ys;
471:     gzs = zs;
472:   }

474:   if (dim == 1) {
475:     VecRestoreArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);
476:   } else if (dim == 2) {
477:     VecRestoreArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);
478:   } else if (dim == 3) {
479:     VecRestoreArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);
480:   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
481:   return(0);
482: }