xref: /petsc/src/dm/impls/da/dagetarray.c (revision 862e4a309d45a165aaa4da0d704ba733429d833a)
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    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           `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 Notes:
89     From Fortran use DMDAVecRestoreArayF90()
90 
91 .seealso: `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    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   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
170 
171   PetscFunctionBegin;
172   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
173   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
174   PetscValidPointer(array, 3);
175   if (da->localSection) {
176     PetscCall(VecGetArrayWrite(vec, (PetscScalar **)array));
177     PetscFunctionReturn(0);
178   }
179   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
180   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
181   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
182 
183   /* Handle case where user passes in global vector as opposed to local */
184   PetscCall(VecGetLocalSize(vec, &N));
185   if (N == xm * ym * zm * dof) {
186     gxm = xm;
187     gym = ym;
188     gzm = zm;
189     gxs = xs;
190     gys = ys;
191     gzs = zs;
192   } 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);
193 
194   if (dim == 1) {
195     PetscCall(VecGetArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
196   } else if (dim == 2) {
197     PetscCall(VecGetArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
198   } else if (dim == 3) {
199     PetscCall(VecGetArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
200   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
201   PetscFunctionReturn(0);
202 }
203 
204 /*@
205    DMDAVecRestoreArrayWrite - Restores a multiple dimension array obtained with DMDAVecGetArrayWrite()
206 
207    Logically collective on Vec
208 
209    Input Parameters:
210 +  da - the distributed array
211 .  vec - the vector, either a vector the same size as one obtained with
212          DMCreateGlobalVector() or DMCreateLocalVector()
213 -  array - the array, non-NULL pointer is zeroed
214 
215   Level: intermediate
216 
217   Fortran Notes:
218     From Fortran use DMDAVecRestoreArayF90()
219 
220 .seealso: `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayWrite()`,
221           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
222 @*/
223 PetscErrorCode DMDAVecRestoreArrayWrite(DM da, Vec vec, void *array)
224 {
225   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
226 
227   PetscFunctionBegin;
228   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
229   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
230   PetscValidPointer(array, 3);
231   if (da->localSection) {
232     PetscCall(VecRestoreArray(vec, (PetscScalar **)array));
233     PetscFunctionReturn(0);
234   }
235   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
236   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
237   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
238 
239   /* Handle case where user passes in global vector as opposed to local */
240   PetscCall(VecGetLocalSize(vec, &N));
241   if (N == xm * ym * zm * dof) {
242     gxm = xm;
243     gym = ym;
244     gzm = zm;
245     gxs = xs;
246     gys = ys;
247     gzs = zs;
248   } 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);
249 
250   if (dim == 1) {
251     PetscCall(VecRestoreArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
252   } else if (dim == 2) {
253     PetscCall(VecRestoreArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
254   } else if (dim == 3) {
255     PetscCall(VecRestoreArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
256   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
257   PetscFunctionReturn(0);
258 }
259 
260 /*@C
261    DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
262       the underlying vector and is indexed using the global dimensions.
263 
264    Logically collective
265 
266    Input Parameters:
267 +  da - the distributed array
268 -  vec - the vector, either a vector the same size as one obtained with
269          DMCreateGlobalVector() or DMCreateLocalVector()
270 
271    Output Parameter:
272 .  array - the array
273 
274    Notes:
275     Call DMDAVecRestoreArrayDOF() once you have finished accessing the vector entries.
276 
277     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
278 
279     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayF90() and declare your array with one higher dimension,
280     see src/dm/tutorials/ex11f90.F
281 
282   Level: intermediate
283 
284 .seealso: `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecRestoreArrayDOF()`,
285           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, `DMDAVecGetArrayDOFRead()`
286 @*/
287 PetscErrorCode DMDAVecGetArrayDOF(DM da, Vec vec, void *array)
288 {
289   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
290 
291   PetscFunctionBegin;
292   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
293   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
294   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
295 
296   /* Handle case where user passes in global vector as opposed to local */
297   PetscCall(VecGetLocalSize(vec, &N));
298   if (N == xm * ym * zm * dof) {
299     gxm = xm;
300     gym = ym;
301     gzm = zm;
302     gxs = xs;
303     gys = ys;
304     gzs = zs;
305   } 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);
306 
307   if (dim == 1) {
308     PetscCall(VecGetArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
309   } else if (dim == 2) {
310     PetscCall(VecGetArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
311   } else if (dim == 3) {
312     PetscCall(VecGetArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
313   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
314   PetscFunctionReturn(0);
315 }
316 
317 /*@
318    DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with DMDAVecGetArrayDOF()
319 
320    Logically collective
321 
322    Input Parameters:
323 +  da - the distributed array
324 .  vec - the vector, either a vector the same size as one obtained with
325          DMCreateGlobalVector() or DMCreateLocalVector()
326 -  array - the array
327 
328   Level: intermediate
329 
330 .seealso: `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
331           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
332 @*/
333 PetscErrorCode DMDAVecRestoreArrayDOF(DM da, Vec vec, void *array)
334 {
335   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
336 
337   PetscFunctionBegin;
338   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
339   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
340   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
341 
342   /* Handle case where user passes in global vector as opposed to local */
343   PetscCall(VecGetLocalSize(vec, &N));
344   if (N == xm * ym * zm * dof) {
345     gxm = xm;
346     gym = ym;
347     gzm = zm;
348     gxs = xs;
349     gys = ys;
350     gzs = zs;
351   }
352 
353   if (dim == 1) {
354     PetscCall(VecRestoreArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
355   } else if (dim == 2) {
356     PetscCall(VecRestoreArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
357   } else if (dim == 3) {
358     PetscCall(VecRestoreArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
359   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
360   PetscFunctionReturn(0);
361 }
362 
363 /*@C
364    DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
365       the underlying vector and is indexed using the global dimensions.
366 
367    Not collective
368 
369    Input Parameters:
370 +  da - the distributed array
371 -  vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
372 
373    Output Parameter:
374 .  array - the array
375 
376    Notes:
377     Call DMDAVecRestoreArrayRead() once you have finished accessing the vector entries.
378 
379     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
380 
381     If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
382     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
383 
384     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
385 
386   Fortran Notes:
387     From Fortran use DMDAVecGetArrayReadF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
388        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
389        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
390        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
391        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
392 
393   Due to bugs in the compiler DMDAVecGetArrayReadF90() does not work with gfortran versions before 4.5
394 
395   Level: intermediate
396 
397 .seealso: `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 Notes:
448     From Fortran use DMDAVecRestoreArrayReadF90()
449 
450 .seealso: `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    Notes:
502     Call DMDAVecRestoreArrayDOFRead() once you have finished accessing the vector entries.
503 
504     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
505 
506     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayReadF90() and declare your array with one higher dimension,
507     see src/dm/tutorials/ex11f90.F
508 
509   Level: intermediate
510 
511 .seealso: `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
512           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
513 @*/
514 PetscErrorCode DMDAVecGetArrayDOFRead(DM da, Vec vec, void *array)
515 {
516   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
517 
518   PetscFunctionBegin;
519   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
520   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
521   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
522 
523   /* Handle case where user passes in global vector as opposed to local */
524   PetscCall(VecGetLocalSize(vec, &N));
525   if (N == xm * ym * zm * dof) {
526     gxm = xm;
527     gym = ym;
528     gzm = zm;
529     gxs = xs;
530     gys = ys;
531     gzs = zs;
532   } 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);
533 
534   if (dim == 1) {
535     PetscCall(VecGetArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
536   } else if (dim == 2) {
537     PetscCall(VecGetArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
538   } else if (dim == 3) {
539     PetscCall(VecGetArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
540   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
541   PetscFunctionReturn(0);
542 }
543 
544 /*@
545    DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with DMDAVecGetArrayDOFRead()
546 
547    Not Collective
548 
549    Input Parameters:
550 +  da - the distributed array
551 .  vec - the vector, either a vector the same size as one obtained with
552          DMCreateGlobalVector() or DMCreateLocalVector()
553 -  array - the array
554 
555   Level: intermediate
556 
557 .seealso: `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
558           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
559 @*/
560 PetscErrorCode DMDAVecRestoreArrayDOFRead(DM da, Vec vec, void *array)
561 {
562   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
563 
564   PetscFunctionBegin;
565   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
566   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
567   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
568 
569   /* Handle case where user passes in global vector as opposed to local */
570   PetscCall(VecGetLocalSize(vec, &N));
571   if (N == xm * ym * zm * dof) {
572     gxm = xm;
573     gym = ym;
574     gzm = zm;
575     gxs = xs;
576     gys = ys;
577     gzs = zs;
578   }
579 
580   if (dim == 1) {
581     PetscCall(VecRestoreArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
582   } else if (dim == 2) {
583     PetscCall(VecRestoreArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
584   } else if (dim == 3) {
585     PetscCall(VecRestoreArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
586   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
587   PetscFunctionReturn(0);
588 }
589 
590 /*@C
591    DMDAVecGetArrayDOFWrite - Returns a multiple dimension array that shares data with
592       the underlying vector and is indexed using the global dimensions.
593 
594    Not Collective
595 
596    Input Parameters:
597 +  da - the distributed array
598 -  vec - the vector, either a vector the same size as one obtained with
599          DMCreateGlobalVector() or DMCreateLocalVector()
600 
601    Output Parameter:
602 .  array - the array
603 
604    Notes:
605     Call DMDAVecRestoreArrayDOFWrite() once you have finished accessing the vector entries.
606 
607     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
608 
609     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayWriteF90() and declare your array with one higher dimension,
610     see src/dm/tutorials/ex11f90.F
611 
612   Level: intermediate
613 
614 .seealso: `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
615           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
616 @*/
617 PetscErrorCode DMDAVecGetArrayDOFWrite(DM da, Vec vec, void *array)
618 {
619   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
620 
621   PetscFunctionBegin;
622   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
623   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
624   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
625 
626   /* Handle case where user passes in global vector as opposed to local */
627   PetscCall(VecGetLocalSize(vec, &N));
628   if (N == xm * ym * zm * dof) {
629     gxm = xm;
630     gym = ym;
631     gzm = zm;
632     gxs = xs;
633     gys = ys;
634     gzs = zs;
635   } 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);
636 
637   if (dim == 1) {
638     PetscCall(VecGetArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
639   } else if (dim == 2) {
640     PetscCall(VecGetArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
641   } else if (dim == 3) {
642     PetscCall(VecGetArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
643   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
644   PetscFunctionReturn(0);
645 }
646 
647 /*@
648    DMDAVecRestoreArrayDOFWrite - Restores a multiple dimension array obtained with DMDAVecGetArrayDOFWrite()
649 
650    Not Collective
651 
652    Input Parameters:
653 +  da - the distributed array
654 .  vec - the vector, either a vector the same size as one obtained with
655          DMCreateGlobalVector() or DMCreateLocalVector()
656 -  array - the array
657 
658   Level: intermediate
659 
660 .seealso: `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
661           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
662 @*/
663 PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM da, Vec vec, void *array)
664 {
665   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
666 
667   PetscFunctionBegin;
668   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
669   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
670   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
671 
672   /* Handle case where user passes in global vector as opposed to local */
673   PetscCall(VecGetLocalSize(vec, &N));
674   if (N == xm * ym * zm * dof) {
675     gxm = xm;
676     gym = ym;
677     gzm = zm;
678     gxs = xs;
679     gys = ys;
680     gzs = zs;
681   }
682 
683   if (dim == 1) {
684     PetscCall(VecRestoreArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
685   } else if (dim == 2) {
686     PetscCall(VecRestoreArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
687   } else if (dim == 3) {
688     PetscCall(VecRestoreArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
689   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
690   PetscFunctionReturn(0);
691 }
692