Actual source code: dagetarray.c
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 da
10: Input Parameters:
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:
28: From Fortran use DMDAVecGetArrayF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
29: 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
30: 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
31: array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
32: DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
34: Due to bugs in the compiler DMDAVecGetArrayF90() does not work with gfortran versions before 4.5
36: Level: intermediate
38: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF()
39: DMDAVecGetArrayDOF(), DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(),
40: DMStagVecGetArray()
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: DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
52: DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
53: DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
55: /* Handle case where user passes in global vector as opposed to local */
56: VecGetLocalSize(vec,&N);
57: if (N == xm*ym*zm*dof) {
58: gxm = xm;
59: gym = ym;
60: gzm = zm;
61: gxs = xs;
62: gys = ys;
63: gzs = zs;
64: } 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);
66: if (dim == 1) {
67: VecGetArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);
68: } else if (dim == 2) {
69: VecGetArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
70: } else if (dim == 3) {
71: VecGetArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
72: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
73: return(0);
74: }
76: /*@
77: DMDAVecRestoreArray - Restores a multiple dimension array obtained with DMDAVecGetArray()
79: Logically collective on da
81: Input Parameters:
82: + da - the distributed array
83: . vec - the vector, either a vector the same size as one obtained with
84: DMCreateGlobalVector() or DMCreateLocalVector()
85: - array - the array, non-NULL pointer is zeroed
87: Level: intermediate
89: Fortran Notes:
90: From Fortran use DMDAVecRestoreArayF90()
92: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(),
93: DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(),
94: DMDStagVecRestoreArray()
95: @*/
96: PetscErrorCode DMDAVecRestoreArray(DM da,Vec vec,void *array)
97: {
99: PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
105: DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
106: DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
107: DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
109: /* Handle case where user passes in global vector as opposed to local */
110: VecGetLocalSize(vec,&N);
111: if (N == xm*ym*zm*dof) {
112: gxm = xm;
113: gym = ym;
114: gzm = zm;
115: gxs = xs;
116: gys = ys;
117: gzs = zs;
118: } 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);
120: if (dim == 1) {
121: VecRestoreArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);
122: } else if (dim == 2) {
123: VecRestoreArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
124: } else if (dim == 3) {
125: VecRestoreArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
126: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
127: return(0);
128: }
130: /*@C
131: DMDAVecGetArrayWrite - Returns a multiple dimension array that shares data with
132: the underlying vector and is indexed using the global dimensions.
134: Logically collective on Vec
136: Input Parameters:
137: + da - the distributed array
138: - vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
140: Output Parameter:
141: . array - the array
143: Notes:
144: Call DMDAVecRestoreArray() once you have finished accessing the vector entries.
146: In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
148: If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
149: a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
151: appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
153: Fortran Notes:
154: From Fortran use DMDAVecGetArrayF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
155: 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
156: 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
157: array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
158: DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
160: Due to bugs in the compiler DMDAVecGetArrayF90() does not work with gfortran versions before 4.5
162: Level: intermediate
164: Developer Notes: This has code duplication with DMDAVecGetArray() and DMDAVecGetArrayRead()
166: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArrayWrite(), DMDAVecRestoreArrayDOF()
167: DMDAVecGetArrayDOF(), DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead()
168: @*/
169: PetscErrorCode DMDAVecGetArrayWrite(DM da,Vec vec,void *array)
170: {
172: PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
178: if (da->localSection) {
179: VecGetArrayWrite(vec,(PetscScalar**)array);
180: return(0);
181: }
182: DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
183: DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
184: DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
186: /* Handle case where user passes in global vector as opposed to local */
187: VecGetLocalSize(vec,&N);
188: if (N == xm*ym*zm*dof) {
189: gxm = xm;
190: gym = ym;
191: gzm = zm;
192: gxs = xs;
193: gys = ys;
194: gzs = zs;
195: } 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);
197: if (dim == 1) {
198: VecGetArray1dWrite(vec,gxm*dof,gxs*dof,(PetscScalar**)array);
199: } else if (dim == 2) {
200: VecGetArray2dWrite(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
201: } else if (dim == 3) {
202: VecGetArray3dWrite(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
203: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
204: return(0);
205: }
207: /*@
208: DMDAVecRestoreArrayWrite - Restores a multiple dimension array obtained with DMDAVecGetArrayWrite()
210: Logically collective on Vec
212: Input Parameters:
213: + da - the distributed array
214: . vec - the vector, either a vector the same size as one obtained with
215: DMCreateGlobalVector() or DMCreateLocalVector()
216: - array - the array, non-NULL pointer is zeroed
218: Level: intermediate
220: Fortran Notes:
221: From Fortran use DMDAVecRestoreArayF90()
223: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArrayWrite(),
224: DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead()
225: @*/
226: PetscErrorCode DMDAVecRestoreArrayWrite(DM da,Vec vec,void *array)
227: {
229: PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
235: if (da->localSection) {
236: VecRestoreArray(vec,(PetscScalar**)array);
237: return(0);
238: }
239: DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
240: DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
241: DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
243: /* Handle case where user passes in global vector as opposed to local */
244: VecGetLocalSize(vec,&N);
245: if (N == xm*ym*zm*dof) {
246: gxm = xm;
247: gym = ym;
248: gzm = zm;
249: gxs = xs;
250: gys = ys;
251: gzs = zs;
252: } 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);
254: if (dim == 1) {
255: VecRestoreArray1dWrite(vec,gxm*dof,gxs*dof,(PetscScalar**)array);
256: } else if (dim == 2) {
257: VecRestoreArray2dWrite(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
258: } else if (dim == 3) {
259: VecRestoreArray3dWrite(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
260: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
261: return(0);
262: }
264: /*@C
265: DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
266: the underlying vector and is indexed using the global dimensions.
268: Logically collective
270: Input Parameters:
271: + da - the distributed array
272: - vec - the vector, either a vector the same size as one obtained with
273: DMCreateGlobalVector() or DMCreateLocalVector()
275: Output Parameter:
276: . array - the array
278: Notes:
279: Call DMDAVecRestoreArrayDOF() once you have finished accessing the vector entries.
281: In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
283: In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use DMDAVecRestoreArrayF90() and declare your array with one higher dimension,
284: see src/dm/tutorials/ex11f90.F
286: Level: intermediate
288: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF(),
289: DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), DMStagVecGetArrayDOF(), DMDAVecGetArrayDOFRead()
290: @*/
291: PetscErrorCode DMDAVecGetArrayDOF(DM da,Vec vec,void *array)
292: {
294: PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
297: DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
298: DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
299: DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
301: /* Handle case where user passes in global vector as opposed to local */
302: VecGetLocalSize(vec,&N);
303: if (N == xm*ym*zm*dof) {
304: gxm = xm;
305: gym = ym;
306: gzm = zm;
307: gxs = xs;
308: gys = ys;
309: gzs = zs;
310: } 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);
312: if (dim == 1) {
313: VecGetArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);
314: } else if (dim == 2) {
315: VecGetArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);
316: } else if (dim == 3) {
317: VecGetArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);
318: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
319: return(0);
320: }
322: /*@
323: DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with DMDAVecGetArrayDOF()
325: Logically collective
327: Input Parameters:
328: + da - the distributed array
329: . vec - the vector, either a vector the same size as one obtained with
330: DMCreateGlobalVector() or DMCreateLocalVector()
331: - array - the array
333: Level: intermediate
335: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF(),
336: DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), DMStagVecRestoreArrayDOF(), DMDAVecRestoreArrayDOFRead()
337: @*/
338: PetscErrorCode DMDAVecRestoreArrayDOF(DM da,Vec vec,void *array)
339: {
341: PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
344: DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
345: DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
346: DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
348: /* Handle case where user passes in global vector as opposed to local */
349: VecGetLocalSize(vec,&N);
350: if (N == xm*ym*zm*dof) {
351: gxm = xm;
352: gym = ym;
353: gzm = zm;
354: gxs = xs;
355: gys = ys;
356: gzs = zs;
357: }
359: if (dim == 1) {
360: VecRestoreArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);
361: } else if (dim == 2) {
362: VecRestoreArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);
363: } else if (dim == 3) {
364: VecRestoreArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);
365: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
366: return(0);
367: }
369: /*@C
370: DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
371: the underlying vector and is indexed using the global dimensions.
373: Not collective
375: Input Parameters:
376: + da - the distributed array
377: - vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
379: Output Parameter:
380: . array - the array
382: Notes:
383: Call DMDAVecRestoreArrayRead() once you have finished accessing the vector entries.
385: In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
387: If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
388: a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
390: appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
392: Fortran Notes:
393: From Fortran use DMDAVecGetArrayReadF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
394: 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
395: 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
396: array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
397: DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
399: Due to bugs in the compiler DMDAVecGetArrayReadF90() does not work with gfortran versions before 4.5
401: Level: intermediate
403: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArrayRead(), DMDAVecRestoreArrayDOF()
404: DMDAVecGetArrayDOF(), DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(),
405: DMStagVecGetArrayRead()
406: @*/
407: PetscErrorCode DMDAVecGetArrayRead(DM da,Vec vec,void *array)
408: {
410: PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
416: DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
417: DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
418: DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
420: /* Handle case where user passes in global vector as opposed to local */
421: VecGetLocalSize(vec,&N);
422: if (N == xm*ym*zm*dof) {
423: gxm = xm;
424: gym = ym;
425: gzm = zm;
426: gxs = xs;
427: gys = ys;
428: gzs = zs;
429: } 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);
431: if (dim == 1) {
432: VecGetArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);
433: } else if (dim == 2) {
434: VecGetArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
435: } else if (dim == 3) {
436: VecGetArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
437: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
438: return(0);
439: }
441: /*@
442: DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with DMDAVecGetArrayRead()
444: Not collective
446: Input Parameters:
447: + da - the distributed array
448: . vec - the vector, either a vector the same size as one obtained with
449: DMCreateGlobalVector() or DMCreateLocalVector()
450: - array - the array, non-NULL pointer is zeroed
452: Level: intermediate
454: Fortran Notes:
455: From Fortran use DMDAVecRestoreArrayReadF90()
457: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArrayRead(),
458: DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(),
459: DMStagVecRestoreArrayRead()
460: @*/
461: PetscErrorCode DMDAVecRestoreArrayRead(DM da,Vec vec,void *array)
462: {
464: PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
470: DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
471: DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
472: DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
474: /* Handle case where user passes in global vector as opposed to local */
475: VecGetLocalSize(vec,&N);
476: if (N == xm*ym*zm*dof) {
477: gxm = xm;
478: gym = ym;
479: gzm = zm;
480: gxs = xs;
481: gys = ys;
482: gzs = zs;
483: } 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);
485: if (dim == 1) {
486: VecRestoreArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);
487: } else if (dim == 2) {
488: VecRestoreArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);
489: } else if (dim == 3) {
490: VecRestoreArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);
491: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
492: return(0);
493: }
495: /*@C
496: DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
497: the underlying vector and is indexed using the global dimensions.
499: Not Collective
501: Input Parameters:
502: + da - the distributed array
503: - vec - the vector, either a vector the same size as one obtained with
504: DMCreateGlobalVector() or DMCreateLocalVector()
506: Output Parameter:
507: . array - the array
509: Notes:
510: Call DMDAVecRestoreArrayDOFRead() once you have finished accessing the vector entries.
512: In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
514: In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use DMDAVecRestoreArrayReadF90() and declare your array with one higher dimension,
515: see src/dm/tutorials/ex11f90.F
517: Level: intermediate
519: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(),
520: DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), DMStagVecGetArrayDOFRead()
521: @*/
522: PetscErrorCode DMDAVecGetArrayDOFRead(DM da,Vec vec,void *array)
523: {
525: PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
528: DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
529: DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
530: DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
532: /* Handle case where user passes in global vector as opposed to local */
533: VecGetLocalSize(vec,&N);
534: if (N == xm*ym*zm*dof) {
535: gxm = xm;
536: gym = ym;
537: gzm = zm;
538: gxs = xs;
539: gys = ys;
540: gzs = zs;
541: } 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);
543: if (dim == 1) {
544: VecGetArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);
545: } else if (dim == 2) {
546: VecGetArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);
547: } else if (dim == 3) {
548: VecGetArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);
549: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
550: return(0);
551: }
553: /*@
554: DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with DMDAVecGetArrayDOFRead()
556: Not Collective
558: Input Parameters:
559: + da - the distributed array
560: . vec - the vector, either a vector the same size as one obtained with
561: DMCreateGlobalVector() or DMCreateLocalVector()
562: - array - the array
564: Level: intermediate
566: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF(),
567: DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), DMStagVecRestoreArrayDOFRead()
568: @*/
569: PetscErrorCode DMDAVecRestoreArrayDOFRead(DM da,Vec vec,void *array)
570: {
572: PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
575: DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
576: DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
577: DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
579: /* Handle case where user passes in global vector as opposed to local */
580: VecGetLocalSize(vec,&N);
581: if (N == xm*ym*zm*dof) {
582: gxm = xm;
583: gym = ym;
584: gzm = zm;
585: gxs = xs;
586: gys = ys;
587: gzs = zs;
588: }
590: if (dim == 1) {
591: VecRestoreArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);
592: } else if (dim == 2) {
593: VecRestoreArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);
594: } else if (dim == 3) {
595: VecRestoreArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);
596: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
597: return(0);
598: }
600: /*@C
601: DMDAVecGetArrayDOFWrite - Returns a multiple dimension array that shares data with
602: the underlying vector and is indexed using the global dimensions.
604: Not Collective
606: Input Parameters:
607: + da - the distributed array
608: - vec - the vector, either a vector the same size as one obtained with
609: DMCreateGlobalVector() or DMCreateLocalVector()
611: Output Parameter:
612: . array - the array
614: Notes:
615: Call DMDAVecRestoreArrayDOFWrite() once you have finished accessing the vector entries.
617: In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
619: In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use DMDAVecRestoreArrayWriteF90() and declare your array with one higher dimension,
620: see src/dm/tutorials/ex11f90.F
622: Level: intermediate
624: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(),
625: DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMStagVecGetArrayDOFWrite(),
626: DMStagVecRestoreArrayDOFRead(), DMStagVecRestoreArrayDOFRead()
627: @*/
628: PetscErrorCode DMDAVecGetArrayDOFWrite(DM da,Vec vec,void *array)
629: {
631: PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
634: DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
635: DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
636: DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
638: /* Handle case where user passes in global vector as opposed to local */
639: VecGetLocalSize(vec,&N);
640: if (N == xm*ym*zm*dof) {
641: gxm = xm;
642: gym = ym;
643: gzm = zm;
644: gxs = xs;
645: gys = ys;
646: gzs = zs;
647: } 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);
649: if (dim == 1) {
650: VecGetArray2dWrite(vec,gxm,dof,gxs,0,(PetscScalar***)array);
651: } else if (dim == 2) {
652: VecGetArray3dWrite(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);
653: } else if (dim == 3) {
654: VecGetArray4dWrite(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);
655: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
656: return(0);
657: }
659: /*@
660: DMDAVecRestoreArrayDOFWrite - Restores a multiple dimension array obtained with DMDAVecGetArrayDOFWrite()
662: Not Collective
664: Input Parameters:
665: + da - the distributed array
666: . vec - the vector, either a vector the same size as one obtained with
667: DMCreateGlobalVector() or DMCreateLocalVector()
668: - array - the array
670: Level: intermediate
672: .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF(),
673: DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMStagVecRestoreArrayDOFWrite(),
674: DMStagVecRestoreArrayDOFRead(), DMStagVecRestoreArrayDOFRead()
675: @*/
676: PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM da,Vec vec,void *array)
677: {
679: PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
682: DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);
683: DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);
684: DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);
686: /* Handle case where user passes in global vector as opposed to local */
687: VecGetLocalSize(vec,&N);
688: if (N == xm*ym*zm*dof) {
689: gxm = xm;
690: gym = ym;
691: gzm = zm;
692: gxs = xs;
693: gys = ys;
694: gzs = zs;
695: }
697: if (dim == 1) {
698: VecRestoreArray2dWrite(vec,gxm,dof,gxs,0,(PetscScalar***)array);
699: } else if (dim == 2) {
700: VecRestoreArray3dWrite(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);
701: } else if (dim == 3) {
702: VecRestoreArray4dWrite(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);
703: } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
704: return(0);
705: }