xref: /petsc/src/dm/impls/da/dagetarray.c (revision 60225df5d8469840be2bf9c1f64795a92b19f3c2)
147c6ae99SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I   "petscdmda.h"   I*/
347c6ae99SBarry Smith 
40f99b6f4SBarry Smith /*MC
50f99b6f4SBarry Smith   DMDAVecGetArrayF90 - check Fortran Notes at `DMDAVecGetArray()`
671b416bcSMatthew G. Knepley 
771b416bcSMatthew G. Knepley   Level: intermediate
80f99b6f4SBarry Smith M*/
90f99b6f4SBarry Smith 
1047c6ae99SBarry Smith /*@C
11aa219208SBarry Smith   DMDAVecGetArray - Returns a multiple dimension array that shares data with
1247c6ae99SBarry Smith   the underlying vector and is indexed using the global dimensions.
1347c6ae99SBarry Smith 
1420f4b53cSBarry Smith   Logically Collective
1547c6ae99SBarry Smith 
16d8d19677SJose E. Roman   Input Parameters:
1747c6ae99SBarry Smith + da  - the distributed array
180af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
1947c6ae99SBarry Smith 
2047c6ae99SBarry Smith   Output Parameter:
2147c6ae99SBarry Smith . array - the array
2247c6ae99SBarry Smith 
23dce8aebaSBarry Smith   Level: intermediate
24dce8aebaSBarry Smith 
2547c6ae99SBarry Smith   Notes:
26dce8aebaSBarry Smith   Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
2747c6ae99SBarry Smith 
2847c6ae99SBarry Smith   In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
2947c6ae99SBarry Smith 
300af9b551SBarry Smith   If `vec` is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
310af9b551SBarry Smith   a global vector then the ghost points are not accessible. Of course, with a local vector you will have had to do the
32dce8aebaSBarry Smith   appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
3347c6ae99SBarry Smith 
340af9b551SBarry Smith   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
350af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
360af9b551SBarry Smith 
3795452b02SPatrick Sanan   Fortran Notes:
380f99b6f4SBarry Smith   Use `DMDAVecGetArrayF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
39dce8aebaSBarry Smith   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
400f99b6f4SBarry Smith   dimension of the `DMDA`.
4147c6ae99SBarry Smith 
420af9b551SBarry Smith   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise
430af9b551SBarry Smith   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
440af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
4547c6ae99SBarry Smith 
46dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
47db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
48db781477SPatrick Sanan           `DMStagVecGetArray()`
4947c6ae99SBarry Smith @*/
50d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArray(DM da, Vec vec, void *array)
51d71ae5a4SJacob Faibussowitsch {
5247c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
5347c6ae99SBarry Smith 
5447c6ae99SBarry Smith   PetscFunctionBegin;
55a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
566db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
576db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
589566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
599566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
609566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
6147c6ae99SBarry Smith 
6247c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
639566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
6447c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
6547c6ae99SBarry Smith     gxm = xm;
6647c6ae99SBarry Smith     gym = ym;
6747c6ae99SBarry Smith     gzm = zm;
6847c6ae99SBarry Smith     gxs = xs;
6947c6ae99SBarry Smith     gys = ys;
7047c6ae99SBarry Smith     gzs = zs;
7163a3b9bcSJacob Faibussowitsch   } 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);
7247c6ae99SBarry Smith 
7347c6ae99SBarry Smith   if (dim == 1) {
749566063dSJacob Faibussowitsch     PetscCall(VecGetArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
7547c6ae99SBarry Smith   } else if (dim == 2) {
769566063dSJacob Faibussowitsch     PetscCall(VecGetArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
7747c6ae99SBarry Smith   } else if (dim == 3) {
789566063dSJacob Faibussowitsch     PetscCall(VecGetArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
7963a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8147c6ae99SBarry Smith }
8247c6ae99SBarry Smith 
830f99b6f4SBarry Smith /*MC
840f99b6f4SBarry Smith   DMDAVecRestoreArrayF90 - check Fortran Notes at `DMDAVecRestoreArray()`
8553c0d4aeSBarry Smith 
8653c0d4aeSBarry Smith   Level: intermediate
870f99b6f4SBarry Smith M*/
880f99b6f4SBarry Smith 
8947c6ae99SBarry Smith /*@
90dce8aebaSBarry Smith   DMDAVecRestoreArray - Restores a multiple dimension array obtained with `DMDAVecGetArray()`
9147c6ae99SBarry Smith 
9220f4b53cSBarry Smith   Logically Collective
9347c6ae99SBarry Smith 
94d8d19677SJose E. Roman   Input Parameters:
9547c6ae99SBarry Smith + da    - the distributed array
960af9b551SBarry Smith . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
970af9b551SBarry Smith - array - the `array` pointer is zeroed
9847c6ae99SBarry Smith 
9947c6ae99SBarry Smith   Level: intermediate
10047c6ae99SBarry Smith 
101*60225df5SJacob Faibussowitsch   Fortran Notes:
1020af9b551SBarry Smith   Use `DMDAVecRestoreArayF90()`
10347c6ae99SBarry Smith 
1040f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `DMDAVecRestoreArayF90()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`,
105db781477SPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
106db781477SPatrick Sanan           `DMDStagVecRestoreArray()`
10747c6ae99SBarry Smith @*/
108d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArray(DM da, Vec vec, void *array)
109d71ae5a4SJacob Faibussowitsch {
11047c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
11147c6ae99SBarry Smith 
11247c6ae99SBarry Smith   PetscFunctionBegin;
113a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
1146db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
1156db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
1169566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
1179566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
1189566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
11947c6ae99SBarry Smith 
12047c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
1219566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
12247c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
12347c6ae99SBarry Smith     gxm = xm;
12447c6ae99SBarry Smith     gym = ym;
12547c6ae99SBarry Smith     gzm = zm;
12647c6ae99SBarry Smith     gxs = xs;
12747c6ae99SBarry Smith     gys = ys;
12847c6ae99SBarry Smith     gzs = zs;
12963a3b9bcSJacob Faibussowitsch   } 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);
13047c6ae99SBarry Smith 
13147c6ae99SBarry Smith   if (dim == 1) {
1329566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
13347c6ae99SBarry Smith   } else if (dim == 2) {
1349566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
13547c6ae99SBarry Smith   } else if (dim == 3) {
1369566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
13763a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
1383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13947c6ae99SBarry Smith }
14047c6ae99SBarry Smith 
1410f99b6f4SBarry Smith /*MC
1420f99b6f4SBarry Smith   DMDAVecGetArrayWriteF90 - check Fortran Notes at `DMDAVecGetArrayWrite()`
14353c0d4aeSBarry Smith 
14453c0d4aeSBarry Smith   Level: intermediate
1450f99b6f4SBarry Smith M*/
1460f99b6f4SBarry Smith 
14747c6ae99SBarry Smith /*@C
148fdc842d1SBarry Smith   DMDAVecGetArrayWrite - Returns a multiple dimension array that shares data with
149fdc842d1SBarry Smith   the underlying vector and is indexed using the global dimensions.
150fdc842d1SBarry Smith 
15120f4b53cSBarry Smith   Logically Collective
152fdc842d1SBarry Smith 
153d8d19677SJose E. Roman   Input Parameters:
154fdc842d1SBarry Smith + da  - the distributed array
1550af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
156fdc842d1SBarry Smith 
157fdc842d1SBarry Smith   Output Parameter:
158fdc842d1SBarry Smith . array - the array
159fdc842d1SBarry Smith 
160dce8aebaSBarry Smith   Level: intermediate
161dce8aebaSBarry Smith 
162fdc842d1SBarry Smith   Notes:
163dce8aebaSBarry Smith   Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
164fdc842d1SBarry Smith 
165fdc842d1SBarry Smith   In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
166fdc842d1SBarry Smith 
1670af9b551SBarry Smith   if `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
168fdc842d1SBarry Smith   a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
169dce8aebaSBarry Smith   appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
170fdc842d1SBarry Smith 
1710af9b551SBarry Smith   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
1720af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
1730af9b551SBarry Smith 
174fdc842d1SBarry Smith   Fortran Notes:
175*60225df5SJacob Faibussowitsch 
1760f99b6f4SBarry Smith   From Fortran use `DMDAVecGetArrayWriteF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
177dce8aebaSBarry Smith   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
1780f99b6f4SBarry Smith   dimension of the `DMDA`.
179fdc842d1SBarry Smith 
1800af9b551SBarry Smith   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise
1810af9b551SBarry Smith   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
1820af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
183fdc842d1SBarry Smith 
184*60225df5SJacob Faibussowitsch   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*60225df5SJacob Faibussowitsch   array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
186*60225df5SJacob Faibussowitsch   `DMDAGetCorners()` for a global array or `DMDAGetGhostCorners()` for a local array.
187*60225df5SJacob Faibussowitsch 
188*60225df5SJacob Faibussowitsch   Developer Notes:
189dce8aebaSBarry Smith   This has code duplication with `DMDAVecGetArray()` and `DMDAVecGetArrayRead()`
190fdc842d1SBarry Smith 
191dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecRestoreArrayDOF()`
192db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
193fdc842d1SBarry Smith @*/
194d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayWrite(DM da, Vec vec, void *array)
195d71ae5a4SJacob Faibussowitsch {
196fdc842d1SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
197fdc842d1SBarry Smith 
198fdc842d1SBarry Smith   PetscFunctionBegin;
199fdc842d1SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
200fdc842d1SBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
201fdc842d1SBarry Smith   PetscValidPointer(array, 3);
2021bb6d2a8SBarry Smith   if (da->localSection) {
2039566063dSJacob Faibussowitsch     PetscCall(VecGetArrayWrite(vec, (PetscScalar **)array));
2043ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
205fdc842d1SBarry Smith   }
2069566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
2079566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
2089566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
209fdc842d1SBarry Smith 
210fdc842d1SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
2119566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
212fdc842d1SBarry Smith   if (N == xm * ym * zm * dof) {
213fdc842d1SBarry Smith     gxm = xm;
214fdc842d1SBarry Smith     gym = ym;
215fdc842d1SBarry Smith     gzm = zm;
216fdc842d1SBarry Smith     gxs = xs;
217fdc842d1SBarry Smith     gys = ys;
218fdc842d1SBarry Smith     gzs = zs;
21963a3b9bcSJacob Faibussowitsch   } 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);
220fdc842d1SBarry Smith 
221fdc842d1SBarry Smith   if (dim == 1) {
2229566063dSJacob Faibussowitsch     PetscCall(VecGetArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
223fdc842d1SBarry Smith   } else if (dim == 2) {
2249566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
225fdc842d1SBarry Smith   } else if (dim == 3) {
2269566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
22763a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
2283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
229fdc842d1SBarry Smith }
230fdc842d1SBarry Smith 
2310f99b6f4SBarry Smith /*MC
2320f99b6f4SBarry Smith   DMDAVecRestoreArrayWriteF90 - check Fortran Notes at `DMDAVecRestoreArrayWrite()`
23353c0d4aeSBarry Smith 
23453c0d4aeSBarry Smith   Level: intermediate
2350f99b6f4SBarry Smith M*/
2360f99b6f4SBarry Smith 
237fdc842d1SBarry Smith /*@
238dce8aebaSBarry Smith   DMDAVecRestoreArrayWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayWrite()`
239fdc842d1SBarry Smith 
24020f4b53cSBarry Smith   Logically Collective
241fdc842d1SBarry Smith 
242d8d19677SJose E. Roman   Input Parameters:
243fdc842d1SBarry Smith + da    - the distributed array
2440af9b551SBarry Smith . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
2450af9b551SBarry Smith - array - the `array` pointer is zeroed
246fdc842d1SBarry Smith 
247fdc842d1SBarry Smith   Level: intermediate
248fdc842d1SBarry Smith 
249*60225df5SJacob Faibussowitsch   Fortran Notes:
2500f99b6f4SBarry Smith   Use `DMDAVecRestoreArayWriteF90()`
251fdc842d1SBarry Smith 
252dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayWrite()`,
253db781477SPatrick Sanan           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
254fdc842d1SBarry Smith @*/
255d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayWrite(DM da, Vec vec, void *array)
256d71ae5a4SJacob Faibussowitsch {
257fdc842d1SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
258fdc842d1SBarry Smith 
259fdc842d1SBarry Smith   PetscFunctionBegin;
260fdc842d1SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
261fdc842d1SBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
262fdc842d1SBarry Smith   PetscValidPointer(array, 3);
2631bb6d2a8SBarry Smith   if (da->localSection) {
2649566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(vec, (PetscScalar **)array));
2653ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
266fdc842d1SBarry Smith   }
2679566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
2689566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
2699566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
270fdc842d1SBarry Smith 
271fdc842d1SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
2729566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
273fdc842d1SBarry Smith   if (N == xm * ym * zm * dof) {
274fdc842d1SBarry Smith     gxm = xm;
275fdc842d1SBarry Smith     gym = ym;
276fdc842d1SBarry Smith     gzm = zm;
277fdc842d1SBarry Smith     gxs = xs;
278fdc842d1SBarry Smith     gys = ys;
279fdc842d1SBarry Smith     gzs = zs;
28063a3b9bcSJacob Faibussowitsch   } 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);
281fdc842d1SBarry Smith 
282fdc842d1SBarry Smith   if (dim == 1) {
2839566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
284fdc842d1SBarry Smith   } else if (dim == 2) {
2859566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
286fdc842d1SBarry Smith   } else if (dim == 3) {
2879566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
28863a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
2893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
290fdc842d1SBarry Smith }
291fdc842d1SBarry Smith 
292fdc842d1SBarry Smith /*@C
293aa219208SBarry Smith   DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
29447c6ae99SBarry Smith   the underlying vector and is indexed using the global dimensions.
29547c6ae99SBarry Smith 
29620f4b53cSBarry Smith   Logically Collective
29747c6ae99SBarry Smith 
298d8d19677SJose E. Roman   Input Parameters:
29947c6ae99SBarry Smith + da  - the distributed array
3000af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
30147c6ae99SBarry Smith 
30247c6ae99SBarry Smith   Output Parameter:
3030af9b551SBarry Smith . array - the `array` pointer is zeroed
30447c6ae99SBarry Smith 
305dce8aebaSBarry Smith   Level: intermediate
306dce8aebaSBarry Smith 
30747c6ae99SBarry Smith   Notes:
308dce8aebaSBarry Smith   Call `DMDAVecRestoreArrayDOF()` once you have finished accessing the vector entries.
30947c6ae99SBarry Smith 
3100af9b551SBarry Smith   In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]
3110af9b551SBarry Smith 
3120af9b551SBarry Smith   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
3130af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
31447c6ae99SBarry Smith 
3150f99b6f4SBarry Smith   Fortran Notes:
3160f99b6f4SBarry Smith   Use `DMDAVecGetArrayF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
3170f99b6f4SBarry Smith   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
3180f99b6f4SBarry Smith   dimension of the `DMDA`.
3190f99b6f4SBarry Smith 
3200af9b551SBarry Smith   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when ndof is 1) otherwise
3210af9b551SBarry Smith   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
3220af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
3231b82215eSBarry Smith 
324dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecRestoreArrayDOF()`,
3254ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, `DMDAVecGetArrayDOFRead()`
32647c6ae99SBarry Smith @*/
327d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOF(DM da, Vec vec, void *array)
328d71ae5a4SJacob Faibussowitsch {
32947c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
33047c6ae99SBarry Smith 
33147c6ae99SBarry Smith   PetscFunctionBegin;
3329566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
3339566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
3349566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
33547c6ae99SBarry Smith 
33647c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
3379566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
33847c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
33947c6ae99SBarry Smith     gxm = xm;
34047c6ae99SBarry Smith     gym = ym;
34147c6ae99SBarry Smith     gzm = zm;
34247c6ae99SBarry Smith     gxs = xs;
34347c6ae99SBarry Smith     gys = ys;
34447c6ae99SBarry Smith     gzs = zs;
34563a3b9bcSJacob Faibussowitsch   } 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);
34647c6ae99SBarry Smith 
34747c6ae99SBarry Smith   if (dim == 1) {
3489566063dSJacob Faibussowitsch     PetscCall(VecGetArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
34947c6ae99SBarry Smith   } else if (dim == 2) {
3509566063dSJacob Faibussowitsch     PetscCall(VecGetArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
35147c6ae99SBarry Smith   } else if (dim == 3) {
3529566063dSJacob Faibussowitsch     PetscCall(VecGetArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
35363a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
3543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
35547c6ae99SBarry Smith }
35647c6ae99SBarry Smith 
35747c6ae99SBarry Smith /*@
358dce8aebaSBarry Smith   DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOF()`
35947c6ae99SBarry Smith 
36020f4b53cSBarry Smith   Logically Collective
36147c6ae99SBarry Smith 
362d8d19677SJose E. Roman   Input Parameters:
36347c6ae99SBarry Smith + da    - the distributed array
3640af9b551SBarry Smith . vec   - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
3650af9b551SBarry Smith - array - the `array` point is zeroed
36647c6ae99SBarry Smith 
36747c6ae99SBarry Smith   Level: intermediate
36847c6ae99SBarry Smith 
369*60225df5SJacob Faibussowitsch   Fortran Notes:
3700f99b6f4SBarry Smith   Use `DMDAVecRestoreArayF90()`
3710f99b6f4SBarry Smith 
372dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
3734ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
37447c6ae99SBarry Smith @*/
375d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOF(DM da, Vec vec, void *array)
376d71ae5a4SJacob Faibussowitsch {
37747c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
37847c6ae99SBarry Smith 
37947c6ae99SBarry Smith   PetscFunctionBegin;
3809566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
3819566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
3829566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
38347c6ae99SBarry Smith 
38447c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
3859566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
38647c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
38747c6ae99SBarry Smith     gxm = xm;
38847c6ae99SBarry Smith     gym = ym;
38947c6ae99SBarry Smith     gzm = zm;
39047c6ae99SBarry Smith     gxs = xs;
39147c6ae99SBarry Smith     gys = ys;
39247c6ae99SBarry Smith     gzs = zs;
39347c6ae99SBarry Smith   }
39447c6ae99SBarry Smith 
39547c6ae99SBarry Smith   if (dim == 1) {
3969566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
39747c6ae99SBarry Smith   } else if (dim == 2) {
3989566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
39947c6ae99SBarry Smith   } else if (dim == 3) {
4009566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
40163a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
4023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40347c6ae99SBarry Smith }
40447c6ae99SBarry Smith 
4050f99b6f4SBarry Smith /*MC
4060f99b6f4SBarry Smith   DMDAVecGetArrayReadF90 - check Fortran Notes at `DMDAVecGetArrayRead()`
40753c0d4aeSBarry Smith 
40853c0d4aeSBarry Smith   Level: intermediate
4090f99b6f4SBarry Smith M*/
4100f99b6f4SBarry Smith 
4115edff71fSBarry Smith /*@C
4125edff71fSBarry Smith   DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
4135edff71fSBarry Smith   the underlying vector and is indexed using the global dimensions.
4145edff71fSBarry Smith 
41520f4b53cSBarry Smith   Not Collective
4165edff71fSBarry Smith 
417d8d19677SJose E. Roman   Input Parameters:
4185edff71fSBarry Smith + da  - the distributed array
4190af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
4205edff71fSBarry Smith 
4215edff71fSBarry Smith   Output Parameter:
4225edff71fSBarry Smith . array - the array
4235edff71fSBarry Smith 
424dce8aebaSBarry Smith   Level: intermediate
425dce8aebaSBarry Smith 
4265edff71fSBarry Smith   Notes:
427dce8aebaSBarry Smith   Call `DMDAVecRestoreArrayRead()` once you have finished accessing the vector entries.
4285edff71fSBarry Smith 
4295edff71fSBarry Smith   In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
4305edff71fSBarry Smith 
4310af9b551SBarry Smith   If `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
4325edff71fSBarry Smith   a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
433dce8aebaSBarry Smith   appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
4345edff71fSBarry Smith 
4350af9b551SBarry Smith   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
4360af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
4370af9b551SBarry Smith 
43895452b02SPatrick Sanan   Fortran Notes:
4390f99b6f4SBarry Smith   Use `DMDAVecGetArrayReadF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
440dce8aebaSBarry Smith   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
4410f99b6f4SBarry Smith   dimension of the `DMDA`.
4420f99b6f4SBarry Smith 
4430af9b551SBarry Smith   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise
4440af9b551SBarry Smith   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
4450af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
4465edff71fSBarry Smith 
4470f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAVecGetArrayReadF90()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayRead()`, `DMDAVecRestoreArrayDOF()`
448*60225df5SJacob Faibussowitsch           `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`,
449db781477SPatrick Sanan           `DMStagVecGetArrayRead()`
4505edff71fSBarry Smith @*/
451d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayRead(DM da, Vec vec, void *array)
452d71ae5a4SJacob Faibussowitsch {
4535edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
4545edff71fSBarry Smith 
4555edff71fSBarry Smith   PetscFunctionBegin;
456a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
4575edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
4585edff71fSBarry Smith   PetscValidPointer(array, 3);
4599566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
4609566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
4619566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
4625edff71fSBarry Smith 
4635edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
4649566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
4655edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
4665edff71fSBarry Smith     gxm = xm;
4675edff71fSBarry Smith     gym = ym;
4685edff71fSBarry Smith     gzm = zm;
4695edff71fSBarry Smith     gxs = xs;
4705edff71fSBarry Smith     gys = ys;
4715edff71fSBarry Smith     gzs = zs;
47263a3b9bcSJacob Faibussowitsch   } 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);
4735edff71fSBarry Smith 
4745edff71fSBarry Smith   if (dim == 1) {
4759566063dSJacob Faibussowitsch     PetscCall(VecGetArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
4765edff71fSBarry Smith   } else if (dim == 2) {
4779566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
4785edff71fSBarry Smith   } else if (dim == 3) {
4799566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
48063a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
4813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4825edff71fSBarry Smith }
4835edff71fSBarry Smith 
4840f99b6f4SBarry Smith /*MC
4850f99b6f4SBarry Smith   DMDAVecRestoreArrayReadF90 - check Fortran Notes at `DMDAVecRestoreArrayRead()`
48653c0d4aeSBarry Smith 
48753c0d4aeSBarry Smith   Level: intermediate
4880f99b6f4SBarry Smith M*/
4890f99b6f4SBarry Smith 
4905edff71fSBarry Smith /*@
491dce8aebaSBarry Smith   DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayRead()`
4925edff71fSBarry Smith 
49320f4b53cSBarry Smith   Not Collective
4945edff71fSBarry Smith 
495d8d19677SJose E. Roman   Input Parameters:
4965edff71fSBarry Smith + da    - the distributed array
4970af9b551SBarry Smith . vec   - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
4980af9b551SBarry Smith - array - the `array` pointer is zeroed
4995edff71fSBarry Smith 
5005edff71fSBarry Smith   Level: intermediate
5015edff71fSBarry Smith 
502*60225df5SJacob Faibussowitsch   Fortran Notes:
5030f99b6f4SBarry Smith   Use `DMDAVecRestoreArrayReadF90()`
5045edff71fSBarry Smith 
5050f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAVecRestoreArrayReadF90()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayRead()`,
506db781477SPatrick Sanan           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`,
507db781477SPatrick Sanan           `DMStagVecRestoreArrayRead()`
5085edff71fSBarry Smith @*/
509d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayRead(DM da, Vec vec, void *array)
510d71ae5a4SJacob Faibussowitsch {
5115edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
5125edff71fSBarry Smith 
5135edff71fSBarry Smith   PetscFunctionBegin;
514a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
5155edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
5165edff71fSBarry Smith   PetscValidPointer(array, 3);
5179566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
5189566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
5199566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5205edff71fSBarry Smith 
5215edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
5229566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
5235edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
5245edff71fSBarry Smith     gxm = xm;
5255edff71fSBarry Smith     gym = ym;
5265edff71fSBarry Smith     gzm = zm;
5275edff71fSBarry Smith     gxs = xs;
5285edff71fSBarry Smith     gys = ys;
5295edff71fSBarry Smith     gzs = zs;
53063a3b9bcSJacob Faibussowitsch   } 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);
5315edff71fSBarry Smith 
5325edff71fSBarry Smith   if (dim == 1) {
5339566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
5345edff71fSBarry Smith   } else if (dim == 2) {
5359566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
5365edff71fSBarry Smith   } else if (dim == 3) {
5379566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
53863a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
5393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5405edff71fSBarry Smith }
5415edff71fSBarry Smith 
5425edff71fSBarry Smith /*@C
5435edff71fSBarry Smith   DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
5445edff71fSBarry Smith   the underlying vector and is indexed using the global dimensions.
5455edff71fSBarry Smith 
5465edff71fSBarry Smith   Not Collective
5475edff71fSBarry Smith 
548d8d19677SJose E. Roman   Input Parameters:
5495edff71fSBarry Smith + da  - the distributed array
5500af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
5515edff71fSBarry Smith 
5525edff71fSBarry Smith   Output Parameter:
5535edff71fSBarry Smith . array - the array
5545edff71fSBarry Smith 
555dce8aebaSBarry Smith   Level: intermediate
556dce8aebaSBarry Smith 
5575edff71fSBarry Smith   Notes:
558dce8aebaSBarry Smith   Call `DMDAVecRestoreArrayDOFRead()` once you have finished accessing the vector entries.
5595edff71fSBarry Smith 
5605edff71fSBarry Smith   In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
5615edff71fSBarry Smith 
5620af9b551SBarry Smith   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
5630af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
5640af9b551SBarry Smith 
565*60225df5SJacob Faibussowitsch   Fortran Notes:
5660f99b6f4SBarry Smith   Use  `DMDAVecGetArrayReadF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
5670f99b6f4SBarry Smith   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
5680f99b6f4SBarry Smith   dimension of the `DMDA`.
5690f99b6f4SBarry Smith 
5700af9b551SBarry Smith   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise
5710af9b551SBarry Smith   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
5720af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
5735edff71fSBarry Smith 
574dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
5754ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
5765edff71fSBarry Smith @*/
577d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOFRead(DM da, Vec vec, void *array)
578d71ae5a4SJacob Faibussowitsch {
5795edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
5805edff71fSBarry Smith 
5815edff71fSBarry Smith   PetscFunctionBegin;
5829566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
5839566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
5849566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5855edff71fSBarry Smith 
5865edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
5879566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
5885edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
5895edff71fSBarry Smith     gxm = xm;
5905edff71fSBarry Smith     gym = ym;
5915edff71fSBarry Smith     gzm = zm;
5925edff71fSBarry Smith     gxs = xs;
5935edff71fSBarry Smith     gys = ys;
5945edff71fSBarry Smith     gzs = zs;
59563a3b9bcSJacob Faibussowitsch   } 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);
5965edff71fSBarry Smith 
5975edff71fSBarry Smith   if (dim == 1) {
5989566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
5995edff71fSBarry Smith   } else if (dim == 2) {
6009566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
6015edff71fSBarry Smith   } else if (dim == 3) {
6029566063dSJacob Faibussowitsch     PetscCall(VecGetArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
60363a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
6043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6055edff71fSBarry Smith }
6065edff71fSBarry Smith 
6075edff71fSBarry Smith /*@
608dce8aebaSBarry Smith   DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFRead()`
6095edff71fSBarry Smith 
6105edff71fSBarry Smith   Not Collective
6115edff71fSBarry Smith 
612d8d19677SJose E. Roman   Input Parameters:
6135edff71fSBarry Smith + da    - the distributed array
6140af9b551SBarry Smith . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
6150af9b551SBarry Smith - array - the `array` pointer is zeroed
6165edff71fSBarry Smith 
6175edff71fSBarry Smith   Level: intermediate
6185edff71fSBarry Smith 
619*60225df5SJacob Faibussowitsch   Fortran Notes:
6200f99b6f4SBarry Smith   Use `DMDAVecRestoreArrayReadF90()`
6210f99b6f4SBarry Smith 
622dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
6234ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
6245edff71fSBarry Smith @*/
625d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOFRead(DM da, Vec vec, void *array)
626d71ae5a4SJacob Faibussowitsch {
6275edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
6285edff71fSBarry Smith 
6295edff71fSBarry Smith   PetscFunctionBegin;
6309566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
6319566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
6329566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
6335edff71fSBarry Smith 
6345edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
6359566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
6365edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
6375edff71fSBarry Smith     gxm = xm;
6385edff71fSBarry Smith     gym = ym;
6395edff71fSBarry Smith     gzm = zm;
6405edff71fSBarry Smith     gxs = xs;
6415edff71fSBarry Smith     gys = ys;
6425edff71fSBarry Smith     gzs = zs;
6435edff71fSBarry Smith   }
6445edff71fSBarry Smith 
6455edff71fSBarry Smith   if (dim == 1) {
6469566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
6475edff71fSBarry Smith   } else if (dim == 2) {
6489566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
6495edff71fSBarry Smith   } else if (dim == 3) {
6509566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
65163a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
6523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6535edff71fSBarry Smith }
6545edff71fSBarry Smith 
6551e5d2365SBarry Smith /*@C
6561e5d2365SBarry Smith   DMDAVecGetArrayDOFWrite - Returns a multiple dimension array that shares data with
6571e5d2365SBarry Smith   the underlying vector and is indexed using the global dimensions.
6581e5d2365SBarry Smith 
6591e5d2365SBarry Smith   Not Collective
6601e5d2365SBarry Smith 
661d8d19677SJose E. Roman   Input Parameters:
6621e5d2365SBarry Smith + da  - the distributed array
6630af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
6641e5d2365SBarry Smith 
6651e5d2365SBarry Smith   Output Parameter:
6661e5d2365SBarry Smith . array - the array
6671e5d2365SBarry Smith 
6680f99b6f4SBarry Smith   Level: intermediate
6690f99b6f4SBarry Smith 
6701e5d2365SBarry Smith   Notes:
671dce8aebaSBarry Smith   Call `DMDAVecRestoreArrayDOFWrite()` once you have finished accessing the vector entries.
6721e5d2365SBarry Smith 
6731e5d2365SBarry Smith   In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
6741e5d2365SBarry Smith 
6750af9b551SBarry Smith   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
6760af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
6770af9b551SBarry Smith 
678*60225df5SJacob Faibussowitsch   Fortran Notes:
6790f99b6f4SBarry Smith   Use  `DMDAVecGetArrayWriteF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
6800f99b6f4SBarry Smith   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
6810f99b6f4SBarry Smith   dimension of the `DMDA`.
6821e5d2365SBarry Smith 
6830af9b551SBarry Smith   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise
6840af9b551SBarry Smith   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
6850af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
6861e5d2365SBarry Smith 
6870f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `DMDAVecRestoreArrayWriteF90()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
688*60225df5SJacob Faibussowitsch           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
6891e5d2365SBarry Smith @*/
690d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOFWrite(DM da, Vec vec, void *array)
691d71ae5a4SJacob Faibussowitsch {
6921e5d2365SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
6931e5d2365SBarry Smith 
6941e5d2365SBarry Smith   PetscFunctionBegin;
6959566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
6969566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
6979566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
6981e5d2365SBarry Smith 
6991e5d2365SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
7009566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
7011e5d2365SBarry Smith   if (N == xm * ym * zm * dof) {
7021e5d2365SBarry Smith     gxm = xm;
7031e5d2365SBarry Smith     gym = ym;
7041e5d2365SBarry Smith     gzm = zm;
7051e5d2365SBarry Smith     gxs = xs;
7061e5d2365SBarry Smith     gys = ys;
7071e5d2365SBarry Smith     gzs = zs;
70863a3b9bcSJacob Faibussowitsch   } 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);
7091e5d2365SBarry Smith 
7101e5d2365SBarry Smith   if (dim == 1) {
7119566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
7121e5d2365SBarry Smith   } else if (dim == 2) {
7139566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
7141e5d2365SBarry Smith   } else if (dim == 3) {
7159566063dSJacob Faibussowitsch     PetscCall(VecGetArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
71663a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
7173ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7181e5d2365SBarry Smith }
7191e5d2365SBarry Smith 
7201e5d2365SBarry Smith /*@
721dce8aebaSBarry Smith   DMDAVecRestoreArrayDOFWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFWrite()`
7221e5d2365SBarry Smith 
7231e5d2365SBarry Smith   Not Collective
7241e5d2365SBarry Smith 
725d8d19677SJose E. Roman   Input Parameters:
7261e5d2365SBarry Smith + da    - the distributed array
7270af9b551SBarry Smith . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
7280af9b551SBarry Smith - array - the `array` pointer is zeroed
7291e5d2365SBarry Smith 
7301e5d2365SBarry Smith   Level: intermediate
7311e5d2365SBarry Smith 
732*60225df5SJacob Faibussowitsch   Fortran Notes:
7330f99b6f4SBarry Smith   Use `DMDAVecRestoreArrayWriteF90()`
7340f99b6f4SBarry Smith 
735dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
736*60225df5SJacob Faibussowitsch           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
7371e5d2365SBarry Smith @*/
738d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM da, Vec vec, void *array)
739d71ae5a4SJacob Faibussowitsch {
7401e5d2365SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
7411e5d2365SBarry Smith 
7421e5d2365SBarry Smith   PetscFunctionBegin;
7439566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
7449566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
7459566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
7461e5d2365SBarry Smith 
7471e5d2365SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
7489566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
7491e5d2365SBarry Smith   if (N == xm * ym * zm * dof) {
7501e5d2365SBarry Smith     gxm = xm;
7511e5d2365SBarry Smith     gym = ym;
7521e5d2365SBarry Smith     gzm = zm;
7531e5d2365SBarry Smith     gxs = xs;
7541e5d2365SBarry Smith     gys = ys;
7551e5d2365SBarry Smith     gzs = zs;
7561e5d2365SBarry Smith   }
7571e5d2365SBarry Smith 
7581e5d2365SBarry Smith   if (dim == 1) {
7599566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
7601e5d2365SBarry Smith   } else if (dim == 2) {
7619566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
7621e5d2365SBarry Smith   } else if (dim == 3) {
7639566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
76463a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
7653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7661e5d2365SBarry Smith }
767