xref: /petsc/src/dm/impls/da/dagetarray.c (revision 3ee9839e4674d2e0e9bcea975330f3458ebcb61e)
1 
2 #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
3 
4 /*@C
5    DMDAVecGetArray - Returns a multiple dimension array that shares data with
6       the underlying vector and is indexed using the global dimensions.
7 
8    Logically collective on da
9 
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()
13 
14    Output Parameter:
15 .  array - the array
16 
17    Notes:
18     Call DMDAVecRestoreArray() once you have finished accessing the vector entries.
19 
20     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
21 
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
24 
25     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
26 
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.
33 
34   Due to bugs in the compiler DMDAVecGetArrayF90() does not work with gfortran versions before 4.5
35 
36   Level: intermediate
37 
38 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF()
39           DMDAVecGetArrayDOF(), DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead()
40 @*/
41 PetscErrorCode  DMDAVecGetArray(DM da,Vec vec,void *array)
42 {
43   PetscErrorCode ierr;
44   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
45 
46   PetscFunctionBegin;
47   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
48   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
49   PetscValidPointer(array, 3);
50   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
51   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
52   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
53 
54   /* Handle case where user passes in global vector as opposed to local */
55   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
56   if (N == xm*ym*zm*dof) {
57     gxm = xm;
58     gym = ym;
59     gzm = zm;
60     gxs = xs;
61     gys = ys;
62     gzs = zs;
63   } 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);
64 
65   if (dim == 1) {
66     ierr = VecGetArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
67   } else if (dim == 2) {
68     ierr = VecGetArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
69   } else if (dim == 3) {
70     ierr = VecGetArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
71   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
72   PetscFunctionReturn(0);
73 }
74 
75 /*@
76    DMDAVecRestoreArray - Restores a multiple dimension array obtained with DMDAVecGetArray()
77 
78    Logically collective on da
79 
80    Input Parameter:
81 +  da - the distributed array
82 .  vec - the vector, either a vector the same size as one obtained with
83          DMCreateGlobalVector() or DMCreateLocalVector()
84 -  array - the array, non-NULL pointer is zeroed
85 
86   Level: intermediate
87 
88   Fortran Notes:
89     From Fortran use DMDAVecRestoreArayF90()
90 
91 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(),
92           DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead()
93 @*/
94 PetscErrorCode  DMDAVecRestoreArray(DM da,Vec vec,void *array)
95 {
96   PetscErrorCode ierr;
97   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
98 
99   PetscFunctionBegin;
100   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
101   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
102   PetscValidPointer(array, 3);
103   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
104   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
105   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
106 
107   /* Handle case where user passes in global vector as opposed to local */
108   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
109   if (N == xm*ym*zm*dof) {
110     gxm = xm;
111     gym = ym;
112     gzm = zm;
113     gxs = xs;
114     gys = ys;
115     gzs = zs;
116   } 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);
117 
118   if (dim == 1) {
119     ierr = VecRestoreArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
120   } else if (dim == 2) {
121     ierr = VecRestoreArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
122   } else if (dim == 3) {
123     ierr = VecRestoreArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
124   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
125   PetscFunctionReturn(0);
126 }
127 
128 /*@C
129    DMDAVecGetArrayWrite - Returns a multiple dimension array that shares data with
130       the underlying vector and is indexed using the global dimensions.
131 
132    Logically collective on Vec
133 
134    Input Parameter:
135 +  da - the distributed array
136 -  vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
137 
138    Output Parameter:
139 .  array - the array
140 
141    Notes:
142     Call DMDAVecRestoreArray() once you have finished accessing the vector entries.
143 
144     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
145 
146     If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
147     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
148 
149     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
150 
151   Fortran Notes:
152     From Fortran use DMDAVecGetArrayF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
153        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
154        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
155        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
156        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
157 
158   Due to bugs in the compiler DMDAVecGetArrayF90() does not work with gfortran versions before 4.5
159 
160   Level: intermediate
161 
162   Developer Notes: This has code duplication with DMDAVecGetArray() and DMDAVecGetArrayRead()
163 
164 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArrayWrite(), DMDAVecRestoreArrayDOF()
165           DMDAVecGetArrayDOF(), DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead()
166 @*/
167 PetscErrorCode  DMDAVecGetArrayWrite(DM da,Vec vec,void *array)
168 {
169   PetscErrorCode ierr;
170   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
171 
172   PetscFunctionBegin;
173   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
174   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
175   PetscValidPointer(array, 3);
176   if (da->localSection) {
177     ierr = VecGetArrayWrite(vec,(PetscScalar**)array);CHKERRQ(ierr);
178     PetscFunctionReturn(0);
179   }
180   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
181   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
182   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
183 
184   /* Handle case where user passes in global vector as opposed to local */
185   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
186   if (N == xm*ym*zm*dof) {
187     gxm = xm;
188     gym = ym;
189     gzm = zm;
190     gxs = xs;
191     gys = ys;
192     gzs = zs;
193   } 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);
194 
195   if (dim == 1) {
196     ierr = VecGetArray1dWrite(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
197   } else if (dim == 2) {
198     ierr = VecGetArray2dWrite(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
199   } else if (dim == 3) {
200     ierr = VecGetArray3dWrite(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
201   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
202   PetscFunctionReturn(0);
203 }
204 
205 /*@
206    DMDAVecRestoreArrayWrite - Restores a multiple dimension array obtained with DMDAVecGetArrayWrite()
207 
208    Logically collective on Vec
209 
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, non-NULL pointer is zeroed
215 
216   Level: intermediate
217 
218   Fortran Notes:
219     From Fortran use DMDAVecRestoreArayF90()
220 
221 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArrayWrite(),
222           DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead()
223 @*/
224 PetscErrorCode  DMDAVecRestoreArrayWrite(DM da,Vec vec,void *array)
225 {
226   PetscErrorCode ierr;
227   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
228 
229   PetscFunctionBegin;
230   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
231   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
232   PetscValidPointer(array, 3);
233   if (da->localSection) {
234     ierr = VecRestoreArray(vec,(PetscScalar**)array);CHKERRQ(ierr);
235     PetscFunctionReturn(0);
236   }
237   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
238   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
239   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
240 
241   /* Handle case where user passes in global vector as opposed to local */
242   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
243   if (N == xm*ym*zm*dof) {
244     gxm = xm;
245     gym = ym;
246     gzm = zm;
247     gxs = xs;
248     gys = ys;
249     gzs = zs;
250   } 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);
251 
252   if (dim == 1) {
253     ierr = VecRestoreArray1dWrite(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
254   } else if (dim == 2) {
255     ierr = VecRestoreArray2dWrite(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
256   } else if (dim == 3) {
257     ierr = VecRestoreArray3dWrite(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
258   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
259   PetscFunctionReturn(0);
260 }
261 
262 /*@C
263    DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
264       the underlying vector and is indexed using the global dimensions.
265 
266    Logically collective
267 
268    Input Parameter:
269 +  da - the distributed array
270 -  vec - the vector, either a vector the same size as one obtained with
271          DMCreateGlobalVector() or DMCreateLocalVector()
272 
273    Output Parameter:
274 .  array - the array
275 
276    Notes:
277     Call DMDAVecRestoreArrayDOF() once you have finished accessing the vector entries.
278 
279     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
280 
281     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayF90() and declare your array with one higher dimension,
282     see src/dm/examples/tutorials/ex11f90.F
283 
284   Level: intermediate
285 
286 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF(),
287           DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead()
288 @*/
289 PetscErrorCode  DMDAVecGetArrayDOF(DM da,Vec vec,void *array)
290 {
291   PetscErrorCode ierr;
292   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
293 
294   PetscFunctionBegin;
295   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
296   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
297   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
298 
299   /* Handle case where user passes in global vector as opposed to local */
300   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
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);
309 
310   if (dim == 1) {
311     ierr = VecGetArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
312   } else if (dim == 2) {
313     ierr = VecGetArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
314   } else if (dim == 3) {
315     ierr = VecGetArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
316   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
317   PetscFunctionReturn(0);
318 }
319 
320 /*@
321    DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with DMDAVecGetArrayDOF()
322 
323    Logically collective
324 
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
330 
331   Level: intermediate
332 
333 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF(),
334           DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead()
335 @*/
336 PetscErrorCode  DMDAVecRestoreArrayDOF(DM da,Vec vec,void *array)
337 {
338   PetscErrorCode ierr;
339   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
340 
341   PetscFunctionBegin;
342   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
343   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
344   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
345 
346   /* Handle case where user passes in global vector as opposed to local */
347   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
348   if (N == xm*ym*zm*dof) {
349     gxm = xm;
350     gym = ym;
351     gzm = zm;
352     gxs = xs;
353     gys = ys;
354     gzs = zs;
355   }
356 
357   if (dim == 1) {
358     ierr = VecRestoreArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
359   } else if (dim == 2) {
360     ierr = VecRestoreArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
361   } else if (dim == 3) {
362     ierr = VecRestoreArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
363   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
364   PetscFunctionReturn(0);
365 }
366 
367 /*@C
368    DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
369       the underlying vector and is indexed using the global dimensions.
370 
371    Not collective
372 
373    Input Parameter:
374 +  da - the distributed array
375 -  vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
376 
377    Output Parameter:
378 .  array - the array
379 
380    Notes:
381     Call DMDAVecRestoreArrayRead() once you have finished accessing the vector entries.
382 
383     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
384 
385     If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
386     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
387 
388     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
389 
390   Fortran Notes:
391     From Fortran use DMDAVecGetArrayReadF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
392        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
393        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
394        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
395        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
396 
397   Due to bugs in the compiler DMDAVecGetArrayReadF90() does not work with gfortran versions before 4.5
398 
399   Level: intermediate
400 
401 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArrayRead(), DMDAVecRestoreArrayDOF()
402           DMDAVecGetArrayDOF(), DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead()
403 @*/
404 PetscErrorCode  DMDAVecGetArrayRead(DM da,Vec vec,void *array)
405 {
406   PetscErrorCode ierr;
407   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
408 
409   PetscFunctionBegin;
410   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
411   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
412   PetscValidPointer(array, 3);
413   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
414   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
415   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
416 
417   /* Handle case where user passes in global vector as opposed to local */
418   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
419   if (N == xm*ym*zm*dof) {
420     gxm = xm;
421     gym = ym;
422     gzm = zm;
423     gxs = xs;
424     gys = ys;
425     gzs = zs;
426   } 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);
427 
428   if (dim == 1) {
429     ierr = VecGetArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
430   } else if (dim == 2) {
431     ierr = VecGetArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
432   } else if (dim == 3) {
433     ierr = VecGetArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
434   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
435   PetscFunctionReturn(0);
436 }
437 
438 /*@
439    DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with DMDAVecGetArrayRead()
440 
441    Not collective
442 
443    Input Parameter:
444 +  da - the distributed array
445 .  vec - the vector, either a vector the same size as one obtained with
446          DMCreateGlobalVector() or DMCreateLocalVector()
447 -  array - the array, non-NULL pointer is zeroed
448 
449   Level: intermediate
450 
451   Fortran Notes:
452     From Fortran use DMDAVecRestoreArrayReadF90()
453 
454 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArrayRead(),
455           DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite()
456 @*/
457 PetscErrorCode  DMDAVecRestoreArrayRead(DM da,Vec vec,void *array)
458 {
459   PetscErrorCode ierr;
460   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
461 
462   PetscFunctionBegin;
463   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
464   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
465   PetscValidPointer(array, 3);
466   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
467   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
468   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
469 
470   /* Handle case where user passes in global vector as opposed to local */
471   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
472   if (N == xm*ym*zm*dof) {
473     gxm = xm;
474     gym = ym;
475     gzm = zm;
476     gxs = xs;
477     gys = ys;
478     gzs = zs;
479   } 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);
480 
481   if (dim == 1) {
482     ierr = VecRestoreArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
483   } else if (dim == 2) {
484     ierr = VecRestoreArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
485   } else if (dim == 3) {
486     ierr = VecRestoreArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
487   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
488   PetscFunctionReturn(0);
489 }
490 
491 /*@C
492    DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
493       the underlying vector and is indexed using the global dimensions.
494 
495    Not Collective
496 
497    Input Parameter:
498 +  da - the distributed array
499 -  vec - the vector, either a vector the same size as one obtained with
500          DMCreateGlobalVector() or DMCreateLocalVector()
501 
502    Output Parameter:
503 .  array - the array
504 
505    Notes:
506     Call DMDAVecRestoreArrayDOFRead() once you have finished accessing the vector entries.
507 
508     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
509 
510     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayReadF90() and declare your array with one higher dimension,
511     see src/dm/examples/tutorials/ex11f90.F
512 
513   Level: intermediate
514 
515 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF(),
516           DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead()
517 @*/
518 PetscErrorCode  DMDAVecGetArrayDOFRead(DM da,Vec vec,void *array)
519 {
520   PetscErrorCode ierr;
521   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
522 
523   PetscFunctionBegin;
524   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
525   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
526   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
527 
528   /* Handle case where user passes in global vector as opposed to local */
529   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
530   if (N == xm*ym*zm*dof) {
531     gxm = xm;
532     gym = ym;
533     gzm = zm;
534     gxs = xs;
535     gys = ys;
536     gzs = zs;
537   } 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);
538 
539   if (dim == 1) {
540     ierr = VecGetArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
541   } else if (dim == 2) {
542     ierr = VecGetArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
543   } else if (dim == 3) {
544     ierr = VecGetArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
545   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
546   PetscFunctionReturn(0);
547 }
548 
549 /*@
550    DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with DMDAVecGetArrayDOFRead()
551 
552    Not Collective
553 
554    Input Parameter:
555 +  da - the distributed array
556 .  vec - the vector, either a vector the same size as one obtained with
557          DMCreateGlobalVector() or DMCreateLocalVector()
558 -  array - the array
559 
560   Level: intermediate
561 
562 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF(),
563           DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead()
564 @*/
565 PetscErrorCode  DMDAVecRestoreArrayDOFRead(DM da,Vec vec,void *array)
566 {
567   PetscErrorCode ierr;
568   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
569 
570   PetscFunctionBegin;
571   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
572   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
573   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
574 
575   /* Handle case where user passes in global vector as opposed to local */
576   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
577   if (N == xm*ym*zm*dof) {
578     gxm = xm;
579     gym = ym;
580     gzm = zm;
581     gxs = xs;
582     gys = ys;
583     gzs = zs;
584   }
585 
586   if (dim == 1) {
587     ierr = VecRestoreArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
588   } else if (dim == 2) {
589     ierr = VecRestoreArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
590   } else if (dim == 3) {
591     ierr = VecRestoreArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
592   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
593   PetscFunctionReturn(0);
594 }
595 
596 
597 
598 
599 
600 
601 
602 
603 
604 
605 
606 
607 
608