xref: /petsc/src/dm/impls/da/dagetarray.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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()`
60f99b6f4SBarry Smith M*/
70f99b6f4SBarry Smith 
847c6ae99SBarry Smith /*@C
9aa219208SBarry Smith    DMDAVecGetArray - Returns a multiple dimension array that shares data with
1047c6ae99SBarry Smith       the underlying vector and is indexed using the global dimensions.
1147c6ae99SBarry Smith 
12d083f849SBarry Smith    Logically collective on da
1347c6ae99SBarry Smith 
14d8d19677SJose E. Roman    Input Parameters:
1547c6ae99SBarry Smith +  da - the distributed array
16dce8aebaSBarry Smith -  vec - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
1747c6ae99SBarry Smith 
1847c6ae99SBarry Smith    Output Parameter:
1947c6ae99SBarry Smith .  array - the array
2047c6ae99SBarry Smith 
21dce8aebaSBarry Smith   Level: intermediate
22dce8aebaSBarry Smith 
2347c6ae99SBarry Smith    Notes:
24dce8aebaSBarry Smith     Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
2547c6ae99SBarry Smith 
2647c6ae99SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
2747c6ae99SBarry Smith 
287e57d48aSBarry Smith     If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
297e57d48aSBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
3047c6ae99SBarry Smith 
31dce8aebaSBarry Smith     appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
3247c6ae99SBarry Smith 
3395452b02SPatrick Sanan   Fortran Notes:
340f99b6f4SBarry Smith   Use `DMDAVecGetArrayF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
35dce8aebaSBarry 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
360f99b6f4SBarry Smith   dimension of the `DMDA`.
3747c6ae99SBarry Smith 
380f99b6f4SBarry 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
390f99b6f4SBarry Smith   array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
400f99b6f4SBarry Smith   `DMDAGetCorners()` for a global array or `DMDAGetGhostCorners()` for a local array.
4147c6ae99SBarry Smith 
42dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
43db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
44db781477SPatrick Sanan           `DMStagVecGetArray()`
4547c6ae99SBarry Smith @*/
46d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArray(DM da, Vec vec, void *array)
47d71ae5a4SJacob Faibussowitsch {
4847c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
4947c6ae99SBarry Smith 
5047c6ae99SBarry Smith   PetscFunctionBegin;
51a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
526db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
536db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
549566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
559566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
569566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5747c6ae99SBarry Smith 
5847c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
599566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
6047c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
6147c6ae99SBarry Smith     gxm = xm;
6247c6ae99SBarry Smith     gym = ym;
6347c6ae99SBarry Smith     gzm = zm;
6447c6ae99SBarry Smith     gxs = xs;
6547c6ae99SBarry Smith     gys = ys;
6647c6ae99SBarry Smith     gzs = zs;
6763a3b9bcSJacob 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);
6847c6ae99SBarry Smith 
6947c6ae99SBarry Smith   if (dim == 1) {
709566063dSJacob Faibussowitsch     PetscCall(VecGetArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
7147c6ae99SBarry Smith   } else if (dim == 2) {
729566063dSJacob Faibussowitsch     PetscCall(VecGetArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
7347c6ae99SBarry Smith   } else if (dim == 3) {
749566063dSJacob Faibussowitsch     PetscCall(VecGetArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
7563a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
76*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7747c6ae99SBarry Smith }
7847c6ae99SBarry Smith 
790f99b6f4SBarry Smith /*MC
800f99b6f4SBarry Smith     DMDAVecRestoreArrayF90 - check Fortran Notes at `DMDAVecRestoreArray()`
810f99b6f4SBarry Smith M*/
820f99b6f4SBarry Smith 
8347c6ae99SBarry Smith /*@
84dce8aebaSBarry Smith    DMDAVecRestoreArray - Restores a multiple dimension array obtained with `DMDAVecGetArray()`
8547c6ae99SBarry Smith 
86d083f849SBarry Smith    Logically collective on da
8747c6ae99SBarry Smith 
88d8d19677SJose E. Roman    Input Parameters:
8947c6ae99SBarry Smith +  da - the distributed array
9047c6ae99SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
910f99b6f4SBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
92e5c84f05SJed Brown -  array - the array, non-NULL pointer is zeroed
9347c6ae99SBarry Smith 
9447c6ae99SBarry Smith   Level: intermediate
9547c6ae99SBarry Smith 
96dce8aebaSBarry Smith   Fortran Note:
97dce8aebaSBarry Smith   From Fortran use `DMDAVecRestoreArayF90()`
9847c6ae99SBarry Smith 
990f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `DMDAVecRestoreArayF90()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`,
100db781477SPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
101db781477SPatrick Sanan           `DMDStagVecRestoreArray()`
10247c6ae99SBarry Smith @*/
103d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArray(DM da, Vec vec, void *array)
104d71ae5a4SJacob Faibussowitsch {
10547c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
10647c6ae99SBarry Smith 
10747c6ae99SBarry Smith   PetscFunctionBegin;
108a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
1096db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
1106db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
1119566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
1129566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
1139566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
11447c6ae99SBarry Smith 
11547c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
1169566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
11747c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
11847c6ae99SBarry Smith     gxm = xm;
11947c6ae99SBarry Smith     gym = ym;
12047c6ae99SBarry Smith     gzm = zm;
12147c6ae99SBarry Smith     gxs = xs;
12247c6ae99SBarry Smith     gys = ys;
12347c6ae99SBarry Smith     gzs = zs;
12463a3b9bcSJacob 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);
12547c6ae99SBarry Smith 
12647c6ae99SBarry Smith   if (dim == 1) {
1279566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
12847c6ae99SBarry Smith   } else if (dim == 2) {
1299566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
13047c6ae99SBarry Smith   } else if (dim == 3) {
1319566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
13263a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
133*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13447c6ae99SBarry Smith }
13547c6ae99SBarry Smith 
1360f99b6f4SBarry Smith /*MC
1370f99b6f4SBarry Smith     DMDAVecGetArrayWriteF90 - check Fortran Notes at `DMDAVecGetArrayWrite()`
1380f99b6f4SBarry Smith M*/
1390f99b6f4SBarry Smith 
14047c6ae99SBarry Smith /*@C
141fdc842d1SBarry Smith    DMDAVecGetArrayWrite - Returns a multiple dimension array that shares data with
142fdc842d1SBarry Smith       the underlying vector and is indexed using the global dimensions.
143fdc842d1SBarry Smith 
144fdc842d1SBarry Smith    Logically collective on Vec
145fdc842d1SBarry Smith 
146d8d19677SJose E. Roman    Input Parameters:
147fdc842d1SBarry Smith +  da - the distributed array
148dce8aebaSBarry Smith -  vec - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
149fdc842d1SBarry Smith 
150fdc842d1SBarry Smith    Output Parameter:
151fdc842d1SBarry Smith .  array - the array
152fdc842d1SBarry Smith 
153dce8aebaSBarry Smith   Level: intermediate
154dce8aebaSBarry Smith 
155fdc842d1SBarry Smith    Notes:
156dce8aebaSBarry Smith     Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
157fdc842d1SBarry Smith 
158fdc842d1SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
159fdc842d1SBarry Smith 
160dce8aebaSBarry Smith     If vec is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
161fdc842d1SBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
162dce8aebaSBarry Smith     appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
163fdc842d1SBarry Smith 
164fdc842d1SBarry Smith    Fortran Notes:
1650f99b6f4SBarry Smith    From Fortran use `DMDAVecGetArrayWriteF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
166dce8aebaSBarry 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
1670f99b6f4SBarry Smith    dimension of the `DMDA`.
168fdc842d1SBarry Smith 
1690f99b6f4SBarry 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
1700f99b6f4SBarry Smith    array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
1710f99b6f4SBarry Smith    `DMDAGetCorners()` for a global array or `DMDAGetGhostCorners()` for a local array.
172fdc842d1SBarry Smith 
173dce8aebaSBarry Smith   Developer Note:
174dce8aebaSBarry Smith   This has code duplication with `DMDAVecGetArray()` and `DMDAVecGetArrayRead()`
175fdc842d1SBarry Smith 
176dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecRestoreArrayDOF()`
177db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
178fdc842d1SBarry Smith @*/
179d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayWrite(DM da, Vec vec, void *array)
180d71ae5a4SJacob Faibussowitsch {
181fdc842d1SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
182fdc842d1SBarry Smith 
183fdc842d1SBarry Smith   PetscFunctionBegin;
184fdc842d1SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
185fdc842d1SBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
186fdc842d1SBarry Smith   PetscValidPointer(array, 3);
1871bb6d2a8SBarry Smith   if (da->localSection) {
1889566063dSJacob Faibussowitsch     PetscCall(VecGetArrayWrite(vec, (PetscScalar **)array));
189*3ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
190fdc842d1SBarry Smith   }
1919566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
1929566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
1939566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
194fdc842d1SBarry Smith 
195fdc842d1SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
1969566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
197fdc842d1SBarry Smith   if (N == xm * ym * zm * dof) {
198fdc842d1SBarry Smith     gxm = xm;
199fdc842d1SBarry Smith     gym = ym;
200fdc842d1SBarry Smith     gzm = zm;
201fdc842d1SBarry Smith     gxs = xs;
202fdc842d1SBarry Smith     gys = ys;
203fdc842d1SBarry Smith     gzs = zs;
20463a3b9bcSJacob 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);
205fdc842d1SBarry Smith 
206fdc842d1SBarry Smith   if (dim == 1) {
2079566063dSJacob Faibussowitsch     PetscCall(VecGetArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
208fdc842d1SBarry Smith   } else if (dim == 2) {
2099566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
210fdc842d1SBarry Smith   } else if (dim == 3) {
2119566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
21263a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
213*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
214fdc842d1SBarry Smith }
215fdc842d1SBarry Smith 
2160f99b6f4SBarry Smith /*MC
2170f99b6f4SBarry Smith     DMDAVecRestoreArrayWriteF90 - check Fortran Notes at `DMDAVecRestoreArrayWrite()`
2180f99b6f4SBarry Smith M*/
2190f99b6f4SBarry Smith 
220fdc842d1SBarry Smith /*@
221dce8aebaSBarry Smith    DMDAVecRestoreArrayWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayWrite()`
222fdc842d1SBarry Smith 
223dce8aebaSBarry Smith    Logically collective on vec
224fdc842d1SBarry Smith 
225d8d19677SJose E. Roman    Input Parameters:
226fdc842d1SBarry Smith +  da - the distributed array
227fdc842d1SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
228dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
229fdc842d1SBarry Smith -  array - the array, non-NULL pointer is zeroed
230fdc842d1SBarry Smith 
231fdc842d1SBarry Smith   Level: intermediate
232fdc842d1SBarry Smith 
233dce8aebaSBarry Smith   Fortran Note:
2340f99b6f4SBarry Smith   Use `DMDAVecRestoreArayWriteF90()`
235fdc842d1SBarry Smith 
236dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayWrite()`,
237db781477SPatrick Sanan           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
238fdc842d1SBarry Smith @*/
239d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayWrite(DM da, Vec vec, void *array)
240d71ae5a4SJacob Faibussowitsch {
241fdc842d1SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
242fdc842d1SBarry Smith 
243fdc842d1SBarry Smith   PetscFunctionBegin;
244fdc842d1SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
245fdc842d1SBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
246fdc842d1SBarry Smith   PetscValidPointer(array, 3);
2471bb6d2a8SBarry Smith   if (da->localSection) {
2489566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(vec, (PetscScalar **)array));
249*3ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
250fdc842d1SBarry Smith   }
2519566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
2529566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
2539566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
254fdc842d1SBarry Smith 
255fdc842d1SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
2569566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
257fdc842d1SBarry Smith   if (N == xm * ym * zm * dof) {
258fdc842d1SBarry Smith     gxm = xm;
259fdc842d1SBarry Smith     gym = ym;
260fdc842d1SBarry Smith     gzm = zm;
261fdc842d1SBarry Smith     gxs = xs;
262fdc842d1SBarry Smith     gys = ys;
263fdc842d1SBarry Smith     gzs = zs;
26463a3b9bcSJacob 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);
265fdc842d1SBarry Smith 
266fdc842d1SBarry Smith   if (dim == 1) {
2679566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
268fdc842d1SBarry Smith   } else if (dim == 2) {
2699566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
270fdc842d1SBarry Smith   } else if (dim == 3) {
2719566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
27263a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
273*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
274fdc842d1SBarry Smith }
275fdc842d1SBarry Smith 
276fdc842d1SBarry Smith /*@C
277aa219208SBarry Smith    DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
27847c6ae99SBarry Smith       the underlying vector and is indexed using the global dimensions.
27947c6ae99SBarry Smith 
2802ddbf755SMatthew G. Knepley    Logically collective
28147c6ae99SBarry Smith 
282d8d19677SJose E. Roman    Input Parameters:
28347c6ae99SBarry Smith +  da - the distributed array
28447c6ae99SBarry Smith -  vec - the vector, either a vector the same size as one obtained with
285dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
28647c6ae99SBarry Smith 
28747c6ae99SBarry Smith    Output Parameter:
28847c6ae99SBarry Smith .  array - the array
28947c6ae99SBarry Smith 
290dce8aebaSBarry Smith   Level: intermediate
291dce8aebaSBarry Smith 
29247c6ae99SBarry Smith    Notes:
293dce8aebaSBarry Smith     Call `DMDAVecRestoreArrayDOF()` once you have finished accessing the vector entries.
29447c6ae99SBarry Smith 
29547c6ae99SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
29647c6ae99SBarry Smith 
2970f99b6f4SBarry Smith    Fortran Notes:
2980f99b6f4SBarry Smith     Use `DMDAVecGetArrayF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
2990f99b6f4SBarry 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
3000f99b6f4SBarry Smith    dimension of the `DMDA`.
3010f99b6f4SBarry Smith 
3020f99b6f4SBarry 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
3030f99b6f4SBarry Smith    array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
3040f99b6f4SBarry Smith    `DMDAGetCorners()` for a global array or `DMDAGetGhostCorners()` for a local array.
3051b82215eSBarry Smith 
306dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecRestoreArrayDOF()`,
3074ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, `DMDAVecGetArrayDOFRead()`
30847c6ae99SBarry Smith @*/
309d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOF(DM da, Vec vec, void *array)
310d71ae5a4SJacob Faibussowitsch {
31147c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
31247c6ae99SBarry Smith 
31347c6ae99SBarry Smith   PetscFunctionBegin;
3149566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
3159566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
3169566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
31747c6ae99SBarry Smith 
31847c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
3199566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
32047c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
32147c6ae99SBarry Smith     gxm = xm;
32247c6ae99SBarry Smith     gym = ym;
32347c6ae99SBarry Smith     gzm = zm;
32447c6ae99SBarry Smith     gxs = xs;
32547c6ae99SBarry Smith     gys = ys;
32647c6ae99SBarry Smith     gzs = zs;
32763a3b9bcSJacob 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);
32847c6ae99SBarry Smith 
32947c6ae99SBarry Smith   if (dim == 1) {
3309566063dSJacob Faibussowitsch     PetscCall(VecGetArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
33147c6ae99SBarry Smith   } else if (dim == 2) {
3329566063dSJacob Faibussowitsch     PetscCall(VecGetArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
33347c6ae99SBarry Smith   } else if (dim == 3) {
3349566063dSJacob Faibussowitsch     PetscCall(VecGetArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
33563a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
336*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
33747c6ae99SBarry Smith }
33847c6ae99SBarry Smith 
33947c6ae99SBarry Smith /*@
340dce8aebaSBarry Smith    DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOF()`
34147c6ae99SBarry Smith 
3422ddbf755SMatthew G. Knepley    Logically collective
34347c6ae99SBarry Smith 
344d8d19677SJose E. Roman    Input Parameters:
34547c6ae99SBarry Smith +  da - the distributed array
34647c6ae99SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
347dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
34847c6ae99SBarry Smith -  array - the array
34947c6ae99SBarry Smith 
35047c6ae99SBarry Smith   Level: intermediate
35147c6ae99SBarry Smith 
3520f99b6f4SBarry Smith   Fortran Note:
3530f99b6f4SBarry Smith   Use `DMDAVecRestoreArayF90()`
3540f99b6f4SBarry Smith 
355dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
3564ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
35747c6ae99SBarry Smith @*/
358d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOF(DM da, Vec vec, void *array)
359d71ae5a4SJacob Faibussowitsch {
36047c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
36147c6ae99SBarry Smith 
36247c6ae99SBarry Smith   PetscFunctionBegin;
3639566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
3649566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
3659566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
36647c6ae99SBarry Smith 
36747c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
3689566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
36947c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
37047c6ae99SBarry Smith     gxm = xm;
37147c6ae99SBarry Smith     gym = ym;
37247c6ae99SBarry Smith     gzm = zm;
37347c6ae99SBarry Smith     gxs = xs;
37447c6ae99SBarry Smith     gys = ys;
37547c6ae99SBarry Smith     gzs = zs;
37647c6ae99SBarry Smith   }
37747c6ae99SBarry Smith 
37847c6ae99SBarry Smith   if (dim == 1) {
3799566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
38047c6ae99SBarry Smith   } else if (dim == 2) {
3819566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
38247c6ae99SBarry Smith   } else if (dim == 3) {
3839566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
38463a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
385*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38647c6ae99SBarry Smith }
38747c6ae99SBarry Smith 
3880f99b6f4SBarry Smith /*MC
3890f99b6f4SBarry Smith     DMDAVecGetArrayReadF90 - check Fortran Notes at `DMDAVecGetArrayRead()`
3900f99b6f4SBarry Smith M*/
3910f99b6f4SBarry Smith 
3925edff71fSBarry Smith /*@C
3935edff71fSBarry Smith    DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
3945edff71fSBarry Smith       the underlying vector and is indexed using the global dimensions.
3955edff71fSBarry Smith 
3962ddbf755SMatthew G. Knepley    Not collective
3975edff71fSBarry Smith 
398d8d19677SJose E. Roman    Input Parameters:
3995edff71fSBarry Smith +  da - the distributed array
400dce8aebaSBarry Smith -  vec - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
4015edff71fSBarry Smith 
4025edff71fSBarry Smith    Output Parameter:
4035edff71fSBarry Smith .  array - the array
4045edff71fSBarry Smith 
405dce8aebaSBarry Smith   Level: intermediate
406dce8aebaSBarry Smith 
4075edff71fSBarry Smith    Notes:
408dce8aebaSBarry Smith     Call `DMDAVecRestoreArrayRead()` once you have finished accessing the vector entries.
4095edff71fSBarry Smith 
4105edff71fSBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
4115edff71fSBarry Smith 
412dce8aebaSBarry Smith     If vec is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
4135edff71fSBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
414dce8aebaSBarry Smith     appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
4155edff71fSBarry Smith 
41695452b02SPatrick Sanan   Fortran Notes:
4170f99b6f4SBarry Smith   Use `DMDAVecGetArrayReadF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
418dce8aebaSBarry 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
4190f99b6f4SBarry Smith   dimension of the `DMDA`.
4200f99b6f4SBarry Smith 
4210f99b6f4SBarry 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
4225edff71fSBarry Smith   array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
4230f99b6f4SBarry Smith   `DMDAGetCorners()` for a global array or `DMDAGetGhostCorners()` for a local array.
4245edff71fSBarry Smith 
4250f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAVecGetArrayReadF90()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayRead()`, `DMDAVecRestoreArrayDOF()`
426db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
427db781477SPatrick Sanan           `DMStagVecGetArrayRead()`
4285edff71fSBarry Smith @*/
429d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayRead(DM da, Vec vec, void *array)
430d71ae5a4SJacob Faibussowitsch {
4315edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
4325edff71fSBarry Smith 
4335edff71fSBarry Smith   PetscFunctionBegin;
434a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
4355edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
4365edff71fSBarry Smith   PetscValidPointer(array, 3);
4379566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
4389566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
4399566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
4405edff71fSBarry Smith 
4415edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
4429566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
4435edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
4445edff71fSBarry Smith     gxm = xm;
4455edff71fSBarry Smith     gym = ym;
4465edff71fSBarry Smith     gzm = zm;
4475edff71fSBarry Smith     gxs = xs;
4485edff71fSBarry Smith     gys = ys;
4495edff71fSBarry Smith     gzs = zs;
45063a3b9bcSJacob 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);
4515edff71fSBarry Smith 
4525edff71fSBarry Smith   if (dim == 1) {
4539566063dSJacob Faibussowitsch     PetscCall(VecGetArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
4545edff71fSBarry Smith   } else if (dim == 2) {
4559566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
4565edff71fSBarry Smith   } else if (dim == 3) {
4579566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
45863a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
459*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4605edff71fSBarry Smith }
4615edff71fSBarry Smith 
4620f99b6f4SBarry Smith /*MC
4630f99b6f4SBarry Smith     DMDAVecRestoreArrayReadF90 - check Fortran Notes at `DMDAVecRestoreArrayRead()`
4640f99b6f4SBarry Smith M*/
4650f99b6f4SBarry Smith 
4665edff71fSBarry Smith /*@
467dce8aebaSBarry Smith    DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayRead()`
4685edff71fSBarry Smith 
4692ddbf755SMatthew G. Knepley    Not collective
4705edff71fSBarry Smith 
471d8d19677SJose E. Roman    Input Parameters:
4725edff71fSBarry Smith +  da - the distributed array
4735edff71fSBarry Smith .  vec - the vector, either a vector the same size as one obtained with
474dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
4755edff71fSBarry Smith -  array - the array, non-NULL pointer is zeroed
4765edff71fSBarry Smith 
4775edff71fSBarry Smith   Level: intermediate
4785edff71fSBarry Smith 
479dce8aebaSBarry Smith   Fortran Note:
4800f99b6f4SBarry Smith   Use `DMDAVecRestoreArrayReadF90()`
4815edff71fSBarry Smith 
4820f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAVecRestoreArrayReadF90()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayRead()`,
483db781477SPatrick Sanan           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`,
484db781477SPatrick Sanan           `DMStagVecRestoreArrayRead()`
4855edff71fSBarry Smith @*/
486d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayRead(DM da, Vec vec, void *array)
487d71ae5a4SJacob Faibussowitsch {
4885edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
4895edff71fSBarry Smith 
4905edff71fSBarry Smith   PetscFunctionBegin;
491a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
4925edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
4935edff71fSBarry Smith   PetscValidPointer(array, 3);
4949566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
4959566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
4969566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
4975edff71fSBarry Smith 
4985edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
4999566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
5005edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
5015edff71fSBarry Smith     gxm = xm;
5025edff71fSBarry Smith     gym = ym;
5035edff71fSBarry Smith     gzm = zm;
5045edff71fSBarry Smith     gxs = xs;
5055edff71fSBarry Smith     gys = ys;
5065edff71fSBarry Smith     gzs = zs;
50763a3b9bcSJacob 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);
5085edff71fSBarry Smith 
5095edff71fSBarry Smith   if (dim == 1) {
5109566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
5115edff71fSBarry Smith   } else if (dim == 2) {
5129566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
5135edff71fSBarry Smith   } else if (dim == 3) {
5149566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
51563a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
516*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5175edff71fSBarry Smith }
5185edff71fSBarry Smith 
5195edff71fSBarry Smith /*@C
5205edff71fSBarry Smith    DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
5215edff71fSBarry Smith       the underlying vector and is indexed using the global dimensions.
5225edff71fSBarry Smith 
5235edff71fSBarry Smith    Not Collective
5245edff71fSBarry Smith 
525d8d19677SJose E. Roman    Input Parameters:
5265edff71fSBarry Smith +  da - the distributed array
5275edff71fSBarry Smith -  vec - the vector, either a vector the same size as one obtained with
528dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
5295edff71fSBarry Smith 
5305edff71fSBarry Smith    Output Parameter:
5315edff71fSBarry Smith .  array - the array
5325edff71fSBarry Smith 
533dce8aebaSBarry Smith   Level: intermediate
534dce8aebaSBarry Smith 
5355edff71fSBarry Smith    Notes:
536dce8aebaSBarry Smith    Call `DMDAVecRestoreArrayDOFRead()` once you have finished accessing the vector entries.
5375edff71fSBarry Smith 
5385edff71fSBarry Smith    In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
5395edff71fSBarry Smith 
540dce8aebaSBarry Smith    Fortran Note:
5410f99b6f4SBarry Smith    Use  `DMDAVecGetArrayReadF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
5420f99b6f4SBarry 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
5430f99b6f4SBarry Smith    dimension of the `DMDA`.
5440f99b6f4SBarry Smith 
5450f99b6f4SBarry 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
5460f99b6f4SBarry Smith    array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
5470f99b6f4SBarry Smith    `DMDAGetCorners()` for a global array or `DMDAGetGhostCorners()` for a local array.
5485edff71fSBarry Smith 
549dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
5504ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
5515edff71fSBarry Smith @*/
552d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOFRead(DM da, Vec vec, void *array)
553d71ae5a4SJacob Faibussowitsch {
5545edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
5555edff71fSBarry Smith 
5565edff71fSBarry Smith   PetscFunctionBegin;
5579566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
5589566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
5599566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5605edff71fSBarry Smith 
5615edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
5629566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
5635edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
5645edff71fSBarry Smith     gxm = xm;
5655edff71fSBarry Smith     gym = ym;
5665edff71fSBarry Smith     gzm = zm;
5675edff71fSBarry Smith     gxs = xs;
5685edff71fSBarry Smith     gys = ys;
5695edff71fSBarry Smith     gzs = zs;
57063a3b9bcSJacob 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);
5715edff71fSBarry Smith 
5725edff71fSBarry Smith   if (dim == 1) {
5739566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
5745edff71fSBarry Smith   } else if (dim == 2) {
5759566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
5765edff71fSBarry Smith   } else if (dim == 3) {
5779566063dSJacob Faibussowitsch     PetscCall(VecGetArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
57863a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
579*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5805edff71fSBarry Smith }
5815edff71fSBarry Smith 
5825edff71fSBarry Smith /*@
583dce8aebaSBarry Smith    DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFRead()`
5845edff71fSBarry Smith 
5855edff71fSBarry Smith    Not Collective
5865edff71fSBarry Smith 
587d8d19677SJose E. Roman    Input Parameters:
5885edff71fSBarry Smith +  da - the distributed array
5895edff71fSBarry Smith .  vec - the vector, either a vector the same size as one obtained with
590dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
5915edff71fSBarry Smith -  array - the array
5925edff71fSBarry Smith 
5935edff71fSBarry Smith   Level: intermediate
5945edff71fSBarry Smith 
5950f99b6f4SBarry Smith   Fortran Note:
5960f99b6f4SBarry Smith   Use `DMDAVecRestoreArrayReadF90()`
5970f99b6f4SBarry Smith 
598dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
5994ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
6005edff71fSBarry Smith @*/
601d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOFRead(DM da, Vec vec, void *array)
602d71ae5a4SJacob Faibussowitsch {
6035edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
6045edff71fSBarry Smith 
6055edff71fSBarry Smith   PetscFunctionBegin;
6069566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
6079566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
6089566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
6095edff71fSBarry Smith 
6105edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
6119566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
6125edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
6135edff71fSBarry Smith     gxm = xm;
6145edff71fSBarry Smith     gym = ym;
6155edff71fSBarry Smith     gzm = zm;
6165edff71fSBarry Smith     gxs = xs;
6175edff71fSBarry Smith     gys = ys;
6185edff71fSBarry Smith     gzs = zs;
6195edff71fSBarry Smith   }
6205edff71fSBarry Smith 
6215edff71fSBarry Smith   if (dim == 1) {
6229566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
6235edff71fSBarry Smith   } else if (dim == 2) {
6249566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
6255edff71fSBarry Smith   } else if (dim == 3) {
6269566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
62763a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
628*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6295edff71fSBarry Smith }
6305edff71fSBarry Smith 
6311e5d2365SBarry Smith /*@C
6321e5d2365SBarry Smith    DMDAVecGetArrayDOFWrite - Returns a multiple dimension array that shares data with
6331e5d2365SBarry Smith       the underlying vector and is indexed using the global dimensions.
6341e5d2365SBarry Smith 
6351e5d2365SBarry Smith    Not Collective
6361e5d2365SBarry Smith 
637d8d19677SJose E. Roman    Input Parameters:
6381e5d2365SBarry Smith +  da - the distributed array
6391e5d2365SBarry Smith -  vec - the vector, either a vector the same size as one obtained with
640dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
6411e5d2365SBarry Smith 
6421e5d2365SBarry Smith    Output Parameter:
6431e5d2365SBarry Smith .  array - the array
6441e5d2365SBarry Smith 
6450f99b6f4SBarry Smith    Level: intermediate
6460f99b6f4SBarry Smith 
6471e5d2365SBarry Smith    Notes:
648dce8aebaSBarry Smith     Call `DMDAVecRestoreArrayDOFWrite()` once you have finished accessing the vector entries.
6491e5d2365SBarry Smith 
6501e5d2365SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
6511e5d2365SBarry Smith 
6520f99b6f4SBarry Smith    Fortran Note:
6530f99b6f4SBarry Smith    Use  `DMDAVecGetArrayWriteF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
6540f99b6f4SBarry 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
6550f99b6f4SBarry Smith    dimension of the `DMDA`.
6561e5d2365SBarry Smith 
6570f99b6f4SBarry 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
6580f99b6f4SBarry Smith    array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
6590f99b6f4SBarry Smith    `DMDAGetCorners()` for a global array or `DMDAGetGhostCorners()` for a local array.
6601e5d2365SBarry Smith 
6610f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `DMDAVecRestoreArrayWriteF90()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
6624ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
6631e5d2365SBarry Smith @*/
664d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOFWrite(DM da, Vec vec, void *array)
665d71ae5a4SJacob Faibussowitsch {
6661e5d2365SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
6671e5d2365SBarry Smith 
6681e5d2365SBarry Smith   PetscFunctionBegin;
6699566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
6709566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
6719566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
6721e5d2365SBarry Smith 
6731e5d2365SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
6749566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
6751e5d2365SBarry Smith   if (N == xm * ym * zm * dof) {
6761e5d2365SBarry Smith     gxm = xm;
6771e5d2365SBarry Smith     gym = ym;
6781e5d2365SBarry Smith     gzm = zm;
6791e5d2365SBarry Smith     gxs = xs;
6801e5d2365SBarry Smith     gys = ys;
6811e5d2365SBarry Smith     gzs = zs;
68263a3b9bcSJacob 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);
6831e5d2365SBarry Smith 
6841e5d2365SBarry Smith   if (dim == 1) {
6859566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
6861e5d2365SBarry Smith   } else if (dim == 2) {
6879566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
6881e5d2365SBarry Smith   } else if (dim == 3) {
6899566063dSJacob Faibussowitsch     PetscCall(VecGetArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
69063a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
691*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6921e5d2365SBarry Smith }
6931e5d2365SBarry Smith 
6941e5d2365SBarry Smith /*@
695dce8aebaSBarry Smith    DMDAVecRestoreArrayDOFWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFWrite()`
6961e5d2365SBarry Smith 
6971e5d2365SBarry Smith    Not Collective
6981e5d2365SBarry Smith 
699d8d19677SJose E. Roman    Input Parameters:
7001e5d2365SBarry Smith +  da - the distributed array
7011e5d2365SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
702dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
7031e5d2365SBarry Smith -  array - the array
7041e5d2365SBarry Smith 
7051e5d2365SBarry Smith   Level: intermediate
7061e5d2365SBarry Smith 
7070f99b6f4SBarry Smith   Fortran Note:
7080f99b6f4SBarry Smith   Use `DMDAVecRestoreArrayWriteF90()`
7090f99b6f4SBarry Smith 
710dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
7114ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
7121e5d2365SBarry Smith @*/
713d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM da, Vec vec, void *array)
714d71ae5a4SJacob Faibussowitsch {
7151e5d2365SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
7161e5d2365SBarry Smith 
7171e5d2365SBarry Smith   PetscFunctionBegin;
7189566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
7199566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
7209566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
7211e5d2365SBarry Smith 
7221e5d2365SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
7239566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
7241e5d2365SBarry Smith   if (N == xm * ym * zm * dof) {
7251e5d2365SBarry Smith     gxm = xm;
7261e5d2365SBarry Smith     gym = ym;
7271e5d2365SBarry Smith     gzm = zm;
7281e5d2365SBarry Smith     gxs = xs;
7291e5d2365SBarry Smith     gys = ys;
7301e5d2365SBarry Smith     gzs = zs;
7311e5d2365SBarry Smith   }
7321e5d2365SBarry Smith 
7331e5d2365SBarry Smith   if (dim == 1) {
7349566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
7351e5d2365SBarry Smith   } else if (dim == 2) {
7369566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
7371e5d2365SBarry Smith   } else if (dim == 3) {
7389566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
73963a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
740*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7411e5d2365SBarry Smith }
742