xref: /petsc/src/dm/impls/da/dagetarray.c (revision 0619917b5a674bb687c64e7daba2ab22be99af31)
1 
2 #include <petsc/private/dmdaimpl.h> /*I   "petscdmda.h"   I*/
3 
4 /*MC
5   DMDAVecGetArrayF90 - check Fortran Notes at `DMDAVecGetArray()`
6 
7   Level: intermediate
8 M*/
9 
10 /*@C
11   DMDAVecGetArray - Returns a multiple dimension array that shares data with
12   the underlying vector and is indexed using the global dimensions.
13 
14   Logically Collective
15 
16   Input Parameters:
17 + da  - the distributed array
18 - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
19 
20   Output Parameter:
21 . array - the array
22 
23   Level: intermediate
24 
25   Notes:
26   Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
27 
28   In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
29 
30   If `vec` is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
31   a global vector then the ghost points are not accessible. Of course, with a local vector you will have had to do the
32   appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
33 
34   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
35   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
36 
37   Fortran Notes:
38   Use `DMDAVecGetArrayF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
39   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
40   dimension of the `DMDA`.
41 
42   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise
43   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
44   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
45 
46 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
47           `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
48           `DMStagVecGetArray()`
49 @*/
50 PetscErrorCode DMDAVecGetArray(DM da, Vec vec, void *array)
51 {
52   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
53 
54   PetscFunctionBegin;
55   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
56   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
57   PetscAssertPointer(array, 3);
58   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
59   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
60   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
61 
62   /* Handle case where user passes in global vector as opposed to local */
63   PetscCall(VecGetLocalSize(vec, &N));
64   if (N == xm * ym * zm * dof) {
65     gxm = xm;
66     gym = ym;
67     gzm = zm;
68     gxs = xs;
69     gys = ys;
70     gzs = zs;
71   } 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);
72 
73   if (dim == 1) {
74     PetscCall(VecGetArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
75   } else if (dim == 2) {
76     PetscCall(VecGetArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
77   } else if (dim == 3) {
78     PetscCall(VecGetArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
79   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
80   PetscFunctionReturn(PETSC_SUCCESS);
81 }
82 
83 /*MC
84   DMDAVecRestoreArrayF90 - check Fortran Notes at `DMDAVecRestoreArray()`
85 
86   Level: intermediate
87 M*/
88 
89 /*@
90   DMDAVecRestoreArray - Restores a multiple dimension array obtained with `DMDAVecGetArray()`
91 
92   Logically Collective
93 
94   Input Parameters:
95 + da    - the distributed array
96 . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
97 - array - the `array` pointer is zeroed
98 
99   Level: intermediate
100 
101   Fortran Notes:
102   Use `DMDAVecRestoreArayF90()`
103 
104 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `DMDAVecRestoreArayF90()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`,
105           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
106           `DMDStagVecRestoreArray()`
107 @*/
108 PetscErrorCode DMDAVecRestoreArray(DM da, Vec vec, void *array)
109 {
110   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
111 
112   PetscFunctionBegin;
113   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
114   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
115   PetscAssertPointer(array, 3);
116   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
117   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
118   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
119 
120   /* Handle case where user passes in global vector as opposed to local */
121   PetscCall(VecGetLocalSize(vec, &N));
122   if (N == xm * ym * zm * dof) {
123     gxm = xm;
124     gym = ym;
125     gzm = zm;
126     gxs = xs;
127     gys = ys;
128     gzs = zs;
129   } 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);
130 
131   if (dim == 1) {
132     PetscCall(VecRestoreArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
133   } else if (dim == 2) {
134     PetscCall(VecRestoreArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
135   } else if (dim == 3) {
136     PetscCall(VecRestoreArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
137   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
138   PetscFunctionReturn(PETSC_SUCCESS);
139 }
140 
141 /*MC
142   DMDAVecGetArrayWriteF90 - check Fortran Notes at `DMDAVecGetArrayWrite()`
143 
144   Level: intermediate
145 M*/
146 
147 /*@C
148   DMDAVecGetArrayWrite - Returns a multiple dimension array that shares data with
149   the underlying vector and is indexed using the global dimensions.
150 
151   Logically Collective
152 
153   Input Parameters:
154 + da  - the distributed array
155 - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
156 
157   Output Parameter:
158 . array - the array
159 
160   Level: intermediate
161 
162   Notes:
163   Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
164 
165   In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
166 
167   if `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
168   a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
169   appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
170 
171   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
172   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
173 
174   Fortran Notes:
175 
176   From Fortran use `DMDAVecGetArrayWriteF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
177   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
178   dimension of the `DMDA`.
179 
180   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise
181   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
182   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
183 
184   The order of the indices is array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) (when dof is 1) otherwise
185   array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
186   `DMDAGetCorners()` for a global array or `DMDAGetGhostCorners()` for a local array.
187 
188   Developer Notes:
189   This has code duplication with `DMDAVecGetArray()` and `DMDAVecGetArrayRead()`
190 
191 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecRestoreArrayDOF()`
192           `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
193 @*/
194 PetscErrorCode DMDAVecGetArrayWrite(DM da, Vec vec, void *array)
195 {
196   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
197 
198   PetscFunctionBegin;
199   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
200   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
201   PetscAssertPointer(array, 3);
202   if (da->localSection) {
203     PetscCall(VecGetArrayWrite(vec, (PetscScalar **)array));
204     PetscFunctionReturn(PETSC_SUCCESS);
205   }
206   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
207   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
208   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
209 
210   /* Handle case where user passes in global vector as opposed to local */
211   PetscCall(VecGetLocalSize(vec, &N));
212   if (N == xm * ym * zm * dof) {
213     gxm = xm;
214     gym = ym;
215     gzm = zm;
216     gxs = xs;
217     gys = ys;
218     gzs = zs;
219   } 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);
220 
221   if (dim == 1) {
222     PetscCall(VecGetArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
223   } else if (dim == 2) {
224     PetscCall(VecGetArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
225   } else if (dim == 3) {
226     PetscCall(VecGetArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
227   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
228   PetscFunctionReturn(PETSC_SUCCESS);
229 }
230 
231 /*MC
232   DMDAVecRestoreArrayWriteF90 - check Fortran Notes at `DMDAVecRestoreArrayWrite()`
233 
234   Level: intermediate
235 M*/
236 
237 /*@
238   DMDAVecRestoreArrayWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayWrite()`
239 
240   Logically Collective
241 
242   Input Parameters:
243 + da    - the distributed array
244 . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
245 - array - the `array` pointer is zeroed
246 
247   Level: intermediate
248 
249   Fortran Notes:
250   Use `DMDAVecRestoreArayWriteF90()`
251 
252 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayWrite()`,
253           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
254 @*/
255 PetscErrorCode DMDAVecRestoreArrayWrite(DM da, Vec vec, void *array)
256 {
257   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
258 
259   PetscFunctionBegin;
260   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
261   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
262   PetscAssertPointer(array, 3);
263   if (da->localSection) {
264     PetscCall(VecRestoreArray(vec, (PetscScalar **)array));
265     PetscFunctionReturn(PETSC_SUCCESS);
266   }
267   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
268   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
269   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
270 
271   /* Handle case where user passes in global vector as opposed to local */
272   PetscCall(VecGetLocalSize(vec, &N));
273   if (N == xm * ym * zm * dof) {
274     gxm = xm;
275     gym = ym;
276     gzm = zm;
277     gxs = xs;
278     gys = ys;
279     gzs = zs;
280   } 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);
281 
282   if (dim == 1) {
283     PetscCall(VecRestoreArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
284   } else if (dim == 2) {
285     PetscCall(VecRestoreArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
286   } else if (dim == 3) {
287     PetscCall(VecRestoreArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
288   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
289   PetscFunctionReturn(PETSC_SUCCESS);
290 }
291 
292 /*@C
293   DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
294   the underlying vector and is indexed using the global dimensions.
295 
296   Logically Collective
297 
298   Input Parameters:
299 + da  - the distributed array
300 - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
301 
302   Output Parameter:
303 . array - the `array` pointer is zeroed
304 
305   Level: intermediate
306 
307   Notes:
308   Call `DMDAVecRestoreArrayDOF()` once you have finished accessing the vector entries.
309 
310   In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]
311 
312   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1][0:ndof-1]` where the values are obtained from
313   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
314 
315   Fortran Notes:
316   Use `DMDAVecGetArrayF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
317   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
318   dimension of the `DMDA`.
319 
320   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when ndof is 1) otherwise
321   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
322   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
323 
324 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecRestoreArrayDOF()`,
325           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, `DMDAVecGetArrayDOFRead()`
326 @*/
327 PetscErrorCode DMDAVecGetArrayDOF(DM da, Vec vec, void *array)
328 {
329   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
330 
331   PetscFunctionBegin;
332   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
333   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
334   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
335 
336   /* Handle case where user passes in global vector as opposed to local */
337   PetscCall(VecGetLocalSize(vec, &N));
338   if (N == xm * ym * zm * dof) {
339     gxm = xm;
340     gym = ym;
341     gzm = zm;
342     gxs = xs;
343     gys = ys;
344     gzs = zs;
345   } 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);
346 
347   if (dim == 1) {
348     PetscCall(VecGetArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
349   } else if (dim == 2) {
350     PetscCall(VecGetArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
351   } else if (dim == 3) {
352     PetscCall(VecGetArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
353   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
354   PetscFunctionReturn(PETSC_SUCCESS);
355 }
356 
357 /*@
358   DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOF()`
359 
360   Logically Collective
361 
362   Input Parameters:
363 + da    - the distributed array
364 . vec   - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
365 - array - the `array` point is zeroed
366 
367   Level: intermediate
368 
369   Fortran Notes:
370   Use `DMDAVecRestoreArayF90()`
371 
372 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
373           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
374 @*/
375 PetscErrorCode DMDAVecRestoreArrayDOF(DM da, Vec vec, void *array)
376 {
377   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
378 
379   PetscFunctionBegin;
380   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
381   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
382   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
383 
384   /* Handle case where user passes in global vector as opposed to local */
385   PetscCall(VecGetLocalSize(vec, &N));
386   if (N == xm * ym * zm * dof) {
387     gxm = xm;
388     gym = ym;
389     gzm = zm;
390     gxs = xs;
391     gys = ys;
392     gzs = zs;
393   }
394 
395   if (dim == 1) {
396     PetscCall(VecRestoreArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
397   } else if (dim == 2) {
398     PetscCall(VecRestoreArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
399   } else if (dim == 3) {
400     PetscCall(VecRestoreArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
401   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
402   PetscFunctionReturn(PETSC_SUCCESS);
403 }
404 
405 /*MC
406   DMDAVecGetArrayReadF90 - check Fortran Notes at `DMDAVecGetArrayRead()`
407 
408   Level: intermediate
409 M*/
410 
411 /*@C
412   DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
413   the underlying vector and is indexed using the global dimensions.
414 
415   Not Collective
416 
417   Input Parameters:
418 + da  - the distributed array
419 - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
420 
421   Output Parameter:
422 . array - the array
423 
424   Level: intermediate
425 
426   Notes:
427   Call `DMDAVecRestoreArrayRead()` once you have finished accessing the vector entries.
428 
429   In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
430 
431   If `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
432   a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
433   appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
434 
435   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
436   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
437 
438   Fortran Notes:
439   Use `DMDAVecGetArrayReadF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
440   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
441   dimension of the `DMDA`.
442 
443   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise
444   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
445   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
446 
447 .seealso: `DM`, `DMDA`, `DMDAVecGetArrayReadF90()`, `DMDAGetGhostCorners()`,
448 `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayRead()`,
449 `DMDAVecRestoreArrayDOF()`, `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`,
450 `DMDAVecRestoreArray()`, `DMStagVecGetArrayRead()`
451 @*/
452 PetscErrorCode DMDAVecGetArrayRead(DM da, Vec vec, void *array)
453 {
454   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
455 
456   PetscFunctionBegin;
457   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
458   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
459   PetscAssertPointer(array, 3);
460   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
461   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
462   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
463 
464   /* Handle case where user passes in global vector as opposed to local */
465   PetscCall(VecGetLocalSize(vec, &N));
466   if (N == xm * ym * zm * dof) {
467     gxm = xm;
468     gym = ym;
469     gzm = zm;
470     gxs = xs;
471     gys = ys;
472     gzs = zs;
473   } 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);
474 
475   if (dim == 1) {
476     PetscCall(VecGetArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
477   } else if (dim == 2) {
478     PetscCall(VecGetArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
479   } else if (dim == 3) {
480     PetscCall(VecGetArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
481   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
482   PetscFunctionReturn(PETSC_SUCCESS);
483 }
484 
485 /*MC
486   DMDAVecRestoreArrayReadF90 - check Fortran Notes at `DMDAVecRestoreArrayRead()`
487 
488   Level: intermediate
489 M*/
490 
491 /*@
492   DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayRead()`
493 
494   Not Collective
495 
496   Input Parameters:
497 + da    - the distributed array
498 . vec   - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
499 - array - the `array` pointer is zeroed
500 
501   Level: intermediate
502 
503   Fortran Notes:
504   Use `DMDAVecRestoreArrayReadF90()`
505 
506 .seealso: `DM`, `DMDA`, `DMDAVecRestoreArrayReadF90()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayRead()`,
507           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`,
508           `DMStagVecRestoreArrayRead()`
509 @*/
510 PetscErrorCode DMDAVecRestoreArrayRead(DM da, Vec vec, void *array)
511 {
512   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
513 
514   PetscFunctionBegin;
515   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
516   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
517   PetscAssertPointer(array, 3);
518   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
519   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
520   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
521 
522   /* Handle case where user passes in global vector as opposed to local */
523   PetscCall(VecGetLocalSize(vec, &N));
524   if (N == xm * ym * zm * dof) {
525     gxm = xm;
526     gym = ym;
527     gzm = zm;
528     gxs = xs;
529     gys = ys;
530     gzs = zs;
531   } 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);
532 
533   if (dim == 1) {
534     PetscCall(VecRestoreArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
535   } else if (dim == 2) {
536     PetscCall(VecRestoreArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
537   } else if (dim == 3) {
538     PetscCall(VecRestoreArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
539   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
540   PetscFunctionReturn(PETSC_SUCCESS);
541 }
542 
543 /*@C
544   DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
545   the underlying vector and is indexed using the global dimensions.
546 
547   Not Collective
548 
549   Input Parameters:
550 + da  - the distributed array
551 - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
552 
553   Output Parameter:
554 . array - the array
555 
556   Level: intermediate
557 
558   Notes:
559   Call `DMDAVecRestoreArrayDOFRead()` once you have finished accessing the vector entries.
560 
561   In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
562 
563   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
564   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
565 
566   Fortran Notes:
567   Use  `DMDAVecGetArrayReadF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
568   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
569   dimension of the `DMDA`.
570 
571   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise
572   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
573   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
574 
575 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
576           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
577 @*/
578 PetscErrorCode DMDAVecGetArrayDOFRead(DM da, Vec vec, void *array)
579 {
580   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
581 
582   PetscFunctionBegin;
583   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
584   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
585   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
586 
587   /* Handle case where user passes in global vector as opposed to local */
588   PetscCall(VecGetLocalSize(vec, &N));
589   if (N == xm * ym * zm * dof) {
590     gxm = xm;
591     gym = ym;
592     gzm = zm;
593     gxs = xs;
594     gys = ys;
595     gzs = zs;
596   } 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);
597 
598   if (dim == 1) {
599     PetscCall(VecGetArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
600   } else if (dim == 2) {
601     PetscCall(VecGetArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
602   } else if (dim == 3) {
603     PetscCall(VecGetArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
604   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
605   PetscFunctionReturn(PETSC_SUCCESS);
606 }
607 
608 /*@
609   DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFRead()`
610 
611   Not Collective
612 
613   Input Parameters:
614 + da    - the distributed array
615 . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
616 - array - the `array` pointer is zeroed
617 
618   Level: intermediate
619 
620   Fortran Notes:
621   Use `DMDAVecRestoreArrayReadF90()`
622 
623 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
624           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
625 @*/
626 PetscErrorCode DMDAVecRestoreArrayDOFRead(DM da, Vec vec, void *array)
627 {
628   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
629 
630   PetscFunctionBegin;
631   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
632   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
633   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
634 
635   /* Handle case where user passes in global vector as opposed to local */
636   PetscCall(VecGetLocalSize(vec, &N));
637   if (N == xm * ym * zm * dof) {
638     gxm = xm;
639     gym = ym;
640     gzm = zm;
641     gxs = xs;
642     gys = ys;
643     gzs = zs;
644   }
645 
646   if (dim == 1) {
647     PetscCall(VecRestoreArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
648   } else if (dim == 2) {
649     PetscCall(VecRestoreArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
650   } else if (dim == 3) {
651     PetscCall(VecRestoreArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
652   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
653   PetscFunctionReturn(PETSC_SUCCESS);
654 }
655 
656 /*@C
657   DMDAVecGetArrayDOFWrite - Returns a multiple dimension array that shares data with
658   the underlying vector and is indexed using the global dimensions.
659 
660   Not Collective
661 
662   Input Parameters:
663 + da  - the distributed array
664 - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
665 
666   Output Parameter:
667 . array - the array
668 
669   Level: intermediate
670 
671   Notes:
672   Call `DMDAVecRestoreArrayDOFWrite()` once you have finished accessing the vector entries.
673 
674   In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
675 
676   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1][0:dof-1]` where the values are obtained from
677   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
678 
679   Fortran Notes:
680   Use  `DMDAVecGetArrayWriteF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
681   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
682   dimension of the `DMDA`.
683 
684   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise
685   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
686   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
687 
688 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `DMDAVecRestoreArrayWriteF90()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
689           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
690 @*/
691 PetscErrorCode DMDAVecGetArrayDOFWrite(DM da, Vec vec, void *array)
692 {
693   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
694 
695   PetscFunctionBegin;
696   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
697   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
698   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
699 
700   /* Handle case where user passes in global vector as opposed to local */
701   PetscCall(VecGetLocalSize(vec, &N));
702   if (N == xm * ym * zm * dof) {
703     gxm = xm;
704     gym = ym;
705     gzm = zm;
706     gxs = xs;
707     gys = ys;
708     gzs = zs;
709   } 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);
710 
711   if (dim == 1) {
712     PetscCall(VecGetArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
713   } else if (dim == 2) {
714     PetscCall(VecGetArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
715   } else if (dim == 3) {
716     PetscCall(VecGetArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
717   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
718   PetscFunctionReturn(PETSC_SUCCESS);
719 }
720 
721 /*@
722   DMDAVecRestoreArrayDOFWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFWrite()`
723 
724   Not Collective
725 
726   Input Parameters:
727 + da    - the distributed array
728 . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
729 - array - the `array` pointer is zeroed
730 
731   Level: intermediate
732 
733   Fortran Notes:
734   Use `DMDAVecRestoreArrayWriteF90()`
735 
736 .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
737           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
738 @*/
739 PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM da, Vec vec, void *array)
740 {
741   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
742 
743   PetscFunctionBegin;
744   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
745   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
746   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
747 
748   /* Handle case where user passes in global vector as opposed to local */
749   PetscCall(VecGetLocalSize(vec, &N));
750   if (N == xm * ym * zm * dof) {
751     gxm = xm;
752     gym = ym;
753     gzm = zm;
754     gxs = xs;
755     gys = ys;
756     gzs = zs;
757   }
758 
759   if (dim == 1) {
760     PetscCall(VecRestoreArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
761   } else if (dim == 2) {
762     PetscCall(VecRestoreArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
763   } else if (dim == 3) {
764     PetscCall(VecRestoreArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
765   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
766   PetscFunctionReturn(PETSC_SUCCESS);
767 }
768