xref: /petsc/src/dm/impls/da/dagetarray.c (revision 4b955ea4ba7f490d0e5a0cddf35875088abbd11f)
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 Parameters:
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   Level: intermediate
18 
19    Notes:
20     Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
21 
22     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
23 
24     If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
25     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
26 
27     appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
28 
29   Fortran Notes:
30     From Fortran use `DMDAVecGetArrayF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
31        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
32        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
33        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
34        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
35 
36   Due to bugs in the compiler `DMDAVecGetArrayF90()` does not work with gfortran versions before 4.5
37 
38 .seealso: `DM`, `DMDA`, `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 {
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   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
51   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
52   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
53 
54   /* Handle case where user passes in global vector as opposed to local */
55   PetscCall(VecGetLocalSize(vec, &N));
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 PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof);
64 
65   if (dim == 1) {
66     PetscCall(VecGetArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
67   } else if (dim == 2) {
68     PetscCall(VecGetArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
69   } else if (dim == 3) {
70     PetscCall(VecGetArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
71   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, 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 Parameters:
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 Note:
89     From Fortran use `DMDAVecRestoreArayF90()`
90 
91 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`,
92           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
93           `DMDStagVecRestoreArray()`
94 @*/
95 PetscErrorCode DMDAVecRestoreArray(DM da, Vec vec, void *array)
96 {
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   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
104   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
105   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
106 
107   /* Handle case where user passes in global vector as opposed to local */
108   PetscCall(VecGetLocalSize(vec, &N));
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 PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof);
117 
118   if (dim == 1) {
119     PetscCall(VecRestoreArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
120   } else if (dim == 2) {
121     PetscCall(VecRestoreArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
122   } else if (dim == 3) {
123     PetscCall(VecRestoreArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
124   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, 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 Parameters:
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   Level: intermediate
142 
143    Notes:
144     Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
145 
146     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
147 
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
150 
151     appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
152 
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.
159 
160   Due to bugs in the compiler `DMDAVecGetArrayF90()` does not work with gfortran versions before 4.5
161 
162   Developer Note:
163   This has code duplication with `DMDAVecGetArray()` and `DMDAVecGetArrayRead()`
164 
165 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecRestoreArrayDOF()`
166           `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
167 @*/
168 PetscErrorCode DMDAVecGetArrayWrite(DM da, Vec vec, void *array)
169 {
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     PetscCall(VecGetArrayWrite(vec, (PetscScalar **)array));
178     PetscFunctionReturn(0);
179   }
180   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
181   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
182   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
183 
184   /* Handle case where user passes in global vector as opposed to local */
185   PetscCall(VecGetLocalSize(vec, &N));
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 PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof);
194 
195   if (dim == 1) {
196     PetscCall(VecGetArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
197   } else if (dim == 2) {
198     PetscCall(VecGetArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
199   } else if (dim == 3) {
200     PetscCall(VecGetArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
201   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, 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 Parameters:
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 Note:
219     From Fortran use `DMDAVecRestoreArayF90()`
220 
221 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayWrite()`,
222           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
223 @*/
224 PetscErrorCode DMDAVecRestoreArrayWrite(DM da, Vec vec, void *array)
225 {
226   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
227 
228   PetscFunctionBegin;
229   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
230   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
231   PetscValidPointer(array, 3);
232   if (da->localSection) {
233     PetscCall(VecRestoreArray(vec, (PetscScalar **)array));
234     PetscFunctionReturn(0);
235   }
236   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
237   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
238   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
239 
240   /* Handle case where user passes in global vector as opposed to local */
241   PetscCall(VecGetLocalSize(vec, &N));
242   if (N == xm * ym * zm * dof) {
243     gxm = xm;
244     gym = ym;
245     gzm = zm;
246     gxs = xs;
247     gys = ys;
248     gzs = zs;
249   } else PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof);
250 
251   if (dim == 1) {
252     PetscCall(VecRestoreArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
253   } else if (dim == 2) {
254     PetscCall(VecRestoreArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
255   } else if (dim == 3) {
256     PetscCall(VecRestoreArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
257   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
258   PetscFunctionReturn(0);
259 }
260 
261 /*@C
262    DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
263       the underlying vector and is indexed using the global dimensions.
264 
265    Logically collective
266 
267    Input Parameters:
268 +  da - the distributed array
269 -  vec - the vector, either a vector the same size as one obtained with
270          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
271 
272    Output Parameter:
273 .  array - the array
274 
275   Level: intermediate
276 
277    Notes:
278     Call `DMDAVecRestoreArrayDOF()` once you have finished accessing the vector entries.
279 
280     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
281 
282     In Fortran 90 you do not need a version of `DMDAVecRestoreArrayDOF()` just use `DMDAVecRestoreArrayF90()` and declare your array with one higher dimension,
283     see src/dm/tutorials/ex11f90.F
284 
285 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecRestoreArrayDOF()`,
286           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, `DMDAVecGetArrayDOFRead()`
287 @*/
288 PetscErrorCode DMDAVecGetArrayDOF(DM da, Vec vec, void *array)
289 {
290   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
291 
292   PetscFunctionBegin;
293   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
294   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
295   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
296 
297   /* Handle case where user passes in global vector as opposed to local */
298   PetscCall(VecGetLocalSize(vec, &N));
299   if (N == xm * ym * zm * dof) {
300     gxm = xm;
301     gym = ym;
302     gzm = zm;
303     gxs = xs;
304     gys = ys;
305     gzs = zs;
306   } else PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof);
307 
308   if (dim == 1) {
309     PetscCall(VecGetArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
310   } else if (dim == 2) {
311     PetscCall(VecGetArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
312   } else if (dim == 3) {
313     PetscCall(VecGetArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
314   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
315   PetscFunctionReturn(0);
316 }
317 
318 /*@
319    DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOF()`
320 
321    Logically collective
322 
323    Input Parameters:
324 +  da - the distributed array
325 .  vec - the vector, either a vector the same size as one obtained with
326          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
327 -  array - the array
328 
329   Level: intermediate
330 
331 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
332           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
333 @*/
334 PetscErrorCode DMDAVecRestoreArrayDOF(DM da, Vec vec, void *array)
335 {
336   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
337 
338   PetscFunctionBegin;
339   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
340   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
341   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
342 
343   /* Handle case where user passes in global vector as opposed to local */
344   PetscCall(VecGetLocalSize(vec, &N));
345   if (N == xm * ym * zm * dof) {
346     gxm = xm;
347     gym = ym;
348     gzm = zm;
349     gxs = xs;
350     gys = ys;
351     gzs = zs;
352   }
353 
354   if (dim == 1) {
355     PetscCall(VecRestoreArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
356   } else if (dim == 2) {
357     PetscCall(VecRestoreArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
358   } else if (dim == 3) {
359     PetscCall(VecRestoreArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
360   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
361   PetscFunctionReturn(0);
362 }
363 
364 /*@C
365    DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
366       the underlying vector and is indexed using the global dimensions.
367 
368    Not collective
369 
370    Input Parameters:
371 +  da - the distributed array
372 -  vec - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
373 
374    Output Parameter:
375 .  array - the array
376 
377   Level: intermediate
378 
379    Notes:
380     Call `DMDAVecRestoreArrayRead()` once you have finished accessing the vector entries.
381 
382     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
383 
384     If vec is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
385     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
386     appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
387 
388   Fortran Notes:
389     From Fortran use `DMDAVecGetArrayReadF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
390        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
391        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
392        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
393        DMDAGetCorners() for a global array or `DMDAGetGhostCorners()` for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
394 
395   Due to bugs in the compiler `DMDAVecGetArrayReadF90()` does not work with gfortran versions before 4.5
396 
397 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayRead()`, `DMDAVecRestoreArrayDOF()`
398           `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
399           `DMStagVecGetArrayRead()`
400 @*/
401 PetscErrorCode DMDAVecGetArrayRead(DM da, Vec vec, void *array)
402 {
403   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
404 
405   PetscFunctionBegin;
406   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
407   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
408   PetscValidPointer(array, 3);
409   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
410   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
411   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
412 
413   /* Handle case where user passes in global vector as opposed to local */
414   PetscCall(VecGetLocalSize(vec, &N));
415   if (N == xm * ym * zm * dof) {
416     gxm = xm;
417     gym = ym;
418     gzm = zm;
419     gxs = xs;
420     gys = ys;
421     gzs = zs;
422   } else PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof);
423 
424   if (dim == 1) {
425     PetscCall(VecGetArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
426   } else if (dim == 2) {
427     PetscCall(VecGetArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
428   } else if (dim == 3) {
429     PetscCall(VecGetArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
430   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
431   PetscFunctionReturn(0);
432 }
433 
434 /*@
435    DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayRead()`
436 
437    Not collective
438 
439    Input Parameters:
440 +  da - the distributed array
441 .  vec - the vector, either a vector the same size as one obtained with
442          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
443 -  array - the array, non-NULL pointer is zeroed
444 
445   Level: intermediate
446 
447   Fortran Note:
448     From Fortran use `DMDAVecRestoreArrayReadF90()`
449 
450 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayRead()`,
451           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`,
452           `DMStagVecRestoreArrayRead()`
453 @*/
454 PetscErrorCode DMDAVecRestoreArrayRead(DM da, Vec vec, void *array)
455 {
456   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
457 
458   PetscFunctionBegin;
459   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
460   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
461   PetscValidPointer(array, 3);
462   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
463   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
464   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
465 
466   /* Handle case where user passes in global vector as opposed to local */
467   PetscCall(VecGetLocalSize(vec, &N));
468   if (N == xm * ym * zm * dof) {
469     gxm = xm;
470     gym = ym;
471     gzm = zm;
472     gxs = xs;
473     gys = ys;
474     gzs = zs;
475   } else PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof);
476 
477   if (dim == 1) {
478     PetscCall(VecRestoreArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
479   } else if (dim == 2) {
480     PetscCall(VecRestoreArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
481   } else if (dim == 3) {
482     PetscCall(VecRestoreArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
483   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
484   PetscFunctionReturn(0);
485 }
486 
487 /*@C
488    DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
489       the underlying vector and is indexed using the global dimensions.
490 
491    Not Collective
492 
493    Input Parameters:
494 +  da - the distributed array
495 -  vec - the vector, either a vector the same size as one obtained with
496          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
497 
498    Output Parameter:
499 .  array - the array
500 
501   Level: intermediate
502 
503    Notes:
504     Call `DMDAVecRestoreArrayDOFRead()` once you have finished accessing the vector entries.
505 
506     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
507 
508    Fortran Note:
509     In Fortran you do not need a version of `DMDAVecRestoreArrayDOF()` just use  `DMDAVecRestoreArrayReadF90()` and declare your array with one higher dimension,
510     see src/dm/tutorials/ex11f90.F
511 
512 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
513           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
514 @*/
515 PetscErrorCode DMDAVecGetArrayDOFRead(DM da, Vec vec, void *array)
516 {
517   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
518 
519   PetscFunctionBegin;
520   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
521   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
522   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
523 
524   /* Handle case where user passes in global vector as opposed to local */
525   PetscCall(VecGetLocalSize(vec, &N));
526   if (N == xm * ym * zm * dof) {
527     gxm = xm;
528     gym = ym;
529     gzm = zm;
530     gxs = xs;
531     gys = ys;
532     gzs = zs;
533   } else PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof);
534 
535   if (dim == 1) {
536     PetscCall(VecGetArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
537   } else if (dim == 2) {
538     PetscCall(VecGetArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
539   } else if (dim == 3) {
540     PetscCall(VecGetArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
541   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
542   PetscFunctionReturn(0);
543 }
544 
545 /*@
546    DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFRead()`
547 
548    Not Collective
549 
550    Input Parameters:
551 +  da - the distributed array
552 .  vec - the vector, either a vector the same size as one obtained with
553          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
554 -  array - the array
555 
556   Level: intermediate
557 
558 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
559           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
560 @*/
561 PetscErrorCode DMDAVecRestoreArrayDOFRead(DM da, Vec vec, void *array)
562 {
563   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
564 
565   PetscFunctionBegin;
566   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
567   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
568   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
569 
570   /* Handle case where user passes in global vector as opposed to local */
571   PetscCall(VecGetLocalSize(vec, &N));
572   if (N == xm * ym * zm * dof) {
573     gxm = xm;
574     gym = ym;
575     gzm = zm;
576     gxs = xs;
577     gys = ys;
578     gzs = zs;
579   }
580 
581   if (dim == 1) {
582     PetscCall(VecRestoreArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
583   } else if (dim == 2) {
584     PetscCall(VecRestoreArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
585   } else if (dim == 3) {
586     PetscCall(VecRestoreArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
587   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
588   PetscFunctionReturn(0);
589 }
590 
591 /*@C
592    DMDAVecGetArrayDOFWrite - Returns a multiple dimension array that shares data with
593       the underlying vector and is indexed using the global dimensions.
594 
595    Not Collective
596 
597    Input Parameters:
598 +  da - the distributed array
599 -  vec - the vector, either a vector the same size as one obtained with
600          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
601 
602    Output Parameter:
603 .  array - the array
604 
605    Notes:
606     Call `DMDAVecRestoreArrayDOFWrite()` once you have finished accessing the vector entries.
607 
608     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
609 
610    Fortran Noe:
611     In Fortran you do not need a version of `DMDAVecRestoreArrayDOF()` just use  `DMDAVecRestoreArrayWriteF90()` and declare your array with one higher dimension,
612     see src/dm/tutorials/ex11f90.F
613 
614   Level: intermediate
615 
616 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
617           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
618 @*/
619 PetscErrorCode DMDAVecGetArrayDOFWrite(DM da, Vec vec, void *array)
620 {
621   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
622 
623   PetscFunctionBegin;
624   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
625   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
626   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
627 
628   /* Handle case where user passes in global vector as opposed to local */
629   PetscCall(VecGetLocalSize(vec, &N));
630   if (N == xm * ym * zm * dof) {
631     gxm = xm;
632     gym = ym;
633     gzm = zm;
634     gxs = xs;
635     gys = ys;
636     gzs = zs;
637   } else PetscCheck(N == gxm * gym * gzm * dof, PETSC_COMM_SELF, PETSC_ERR_ARG_INCOMP, "Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT, N, xm * ym * zm * dof, gxm * gym * gzm * dof);
638 
639   if (dim == 1) {
640     PetscCall(VecGetArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
641   } else if (dim == 2) {
642     PetscCall(VecGetArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
643   } else if (dim == 3) {
644     PetscCall(VecGetArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
645   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
646   PetscFunctionReturn(0);
647 }
648 
649 /*@
650    DMDAVecRestoreArrayDOFWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFWrite()`
651 
652    Not Collective
653 
654    Input Parameters:
655 +  da - the distributed array
656 .  vec - the vector, either a vector the same size as one obtained with
657          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
658 -  array - the array
659 
660   Level: intermediate
661 
662 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
663           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
664 @*/
665 PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM da, Vec vec, void *array)
666 {
667   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
668 
669   PetscFunctionBegin;
670   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
671   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
672   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
673 
674   /* Handle case where user passes in global vector as opposed to local */
675   PetscCall(VecGetLocalSize(vec, &N));
676   if (N == xm * ym * zm * dof) {
677     gxm = xm;
678     gym = ym;
679     gzm = zm;
680     gxs = xs;
681     gys = ys;
682     gzs = zs;
683   }
684 
685   if (dim == 1) {
686     PetscCall(VecRestoreArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
687   } else if (dim == 2) {
688     PetscCall(VecRestoreArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
689   } else if (dim == 3) {
690     PetscCall(VecRestoreArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
691   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
692   PetscFunctionReturn(0);
693 }
694