xref: /petsc/src/dm/impls/da/dagetarray.c (revision ce78bad369055609e946c9d2c25ea67a45873e27)
1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I   "petscdmda.h"   I*/
247c6ae99SBarry Smith 
30f99b6f4SBarry Smith /*MC
412b4a537SBarry Smith   DMDAVecGetArrayF90 - see Fortran Notes at `DMDAVecGetArray()`
571b416bcSMatthew G. Knepley 
671b416bcSMatthew G. Knepley   Level: intermediate
70f99b6f4SBarry Smith M*/
80f99b6f4SBarry Smith 
947c6ae99SBarry Smith /*@C
10aa219208SBarry Smith   DMDAVecGetArray - Returns a multiple dimension array that shares data with
1112b4a537SBarry Smith   the underlying vector and is indexed using the global or local dimensions of a `DMDA`.
1247c6ae99SBarry Smith 
1320f4b53cSBarry Smith   Logically Collective
1447c6ae99SBarry Smith 
15d8d19677SJose E. Roman   Input Parameters:
1672fd0fbdSBarry Smith + da  - the `DMDA`
170af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
1847c6ae99SBarry Smith 
1947c6ae99SBarry Smith   Output Parameter:
2047c6ae99SBarry Smith . array - the array
2147c6ae99SBarry Smith 
22dce8aebaSBarry Smith   Level: intermediate
23dce8aebaSBarry Smith 
2447c6ae99SBarry Smith   Notes:
25dce8aebaSBarry Smith   Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
2647c6ae99SBarry Smith 
2747c6ae99SBarry Smith   In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
2847c6ae99SBarry Smith 
290af9b551SBarry Smith   If `vec` is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
300af9b551SBarry Smith   a global vector then the ghost points are not accessible. Of course, with a local vector you will have had to do the
31dce8aebaSBarry Smith   appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
3247c6ae99SBarry Smith 
330af9b551SBarry Smith   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
340af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
350af9b551SBarry Smith 
3695452b02SPatrick Sanan   Fortran Notes:
37*ce78bad3SBarry Smith   Use `DMDAVecGetArray()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
38dce8aebaSBarry 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
390f99b6f4SBarry Smith   dimension of the `DMDA`.
4047c6ae99SBarry Smith 
410af9b551SBarry 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
420af9b551SBarry Smith   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
430af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
4447c6ae99SBarry Smith 
4512b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
46db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
47db781477SPatrick Sanan           `DMStagVecGetArray()`
4847c6ae99SBarry Smith @*/
49d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArray(DM da, Vec vec, void *array)
50d71ae5a4SJacob Faibussowitsch {
5147c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
5247c6ae99SBarry Smith 
5347c6ae99SBarry Smith   PetscFunctionBegin;
54a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
556db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
564f572ea9SToby Isaac   PetscAssertPointer(array, 3);
579566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
589566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
599566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
6047c6ae99SBarry Smith 
6147c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
629566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
6347c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
6447c6ae99SBarry Smith     gxm = xm;
6547c6ae99SBarry Smith     gym = ym;
6647c6ae99SBarry Smith     gzm = zm;
6747c6ae99SBarry Smith     gxs = xs;
6847c6ae99SBarry Smith     gys = ys;
6947c6ae99SBarry Smith     gzs = zs;
7063a3b9bcSJacob 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);
7147c6ae99SBarry Smith 
7247c6ae99SBarry Smith   if (dim == 1) {
739566063dSJacob Faibussowitsch     PetscCall(VecGetArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
7447c6ae99SBarry Smith   } else if (dim == 2) {
759566063dSJacob Faibussowitsch     PetscCall(VecGetArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
7647c6ae99SBarry Smith   } else if (dim == 3) {
779566063dSJacob Faibussowitsch     PetscCall(VecGetArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
7863a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
793ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8047c6ae99SBarry Smith }
8147c6ae99SBarry Smith 
820f99b6f4SBarry Smith /*MC
8312b4a537SBarry Smith   DMDAVecRestoreArrayF90 - see Fortran Notes at `DMDAVecRestoreArray()`
8453c0d4aeSBarry Smith 
8553c0d4aeSBarry Smith   Level: intermediate
860f99b6f4SBarry Smith M*/
870f99b6f4SBarry Smith 
8847c6ae99SBarry Smith /*@
89dce8aebaSBarry Smith   DMDAVecRestoreArray - Restores a multiple dimension array obtained with `DMDAVecGetArray()`
9047c6ae99SBarry Smith 
9120f4b53cSBarry Smith   Logically Collective
9247c6ae99SBarry Smith 
93d8d19677SJose E. Roman   Input Parameters:
9472fd0fbdSBarry Smith + da    - the `DMDA`
950af9b551SBarry Smith . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
9612b4a537SBarry Smith - array - the `array` pointer
9747c6ae99SBarry Smith 
9847c6ae99SBarry Smith   Level: intermediate
9947c6ae99SBarry Smith 
100*ce78bad3SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`,
101db781477SPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
102f8d70eaaSPierre Jolivet           `DMStagVecRestoreArray()`
10347c6ae99SBarry Smith @*/
104d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArray(DM da, Vec vec, void *array)
105d71ae5a4SJacob Faibussowitsch {
10647c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
10747c6ae99SBarry Smith 
10847c6ae99SBarry Smith   PetscFunctionBegin;
109a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
1106db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
1114f572ea9SToby Isaac   PetscAssertPointer(array, 3);
1129566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
1139566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
1149566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
11547c6ae99SBarry Smith 
11647c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
1179566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
11847c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
11947c6ae99SBarry Smith     gxm = xm;
12047c6ae99SBarry Smith     gym = ym;
12147c6ae99SBarry Smith     gzm = zm;
12247c6ae99SBarry Smith     gxs = xs;
12347c6ae99SBarry Smith     gys = ys;
12447c6ae99SBarry Smith     gzs = zs;
12563a3b9bcSJacob 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);
12647c6ae99SBarry Smith 
12747c6ae99SBarry Smith   if (dim == 1) {
1289566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
12947c6ae99SBarry Smith   } else if (dim == 2) {
1309566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
13147c6ae99SBarry Smith   } else if (dim == 3) {
1329566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
13363a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
1343ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13547c6ae99SBarry Smith }
13647c6ae99SBarry Smith 
1370f99b6f4SBarry Smith /*MC
13812b4a537SBarry Smith   DMDAVecGetArrayWriteF90 - see Fortran Notes at `DMDAVecGetArrayWrite()`
13953c0d4aeSBarry Smith 
14053c0d4aeSBarry Smith   Level: intermediate
1410f99b6f4SBarry Smith M*/
1420f99b6f4SBarry Smith 
14347c6ae99SBarry Smith /*@C
144fdc842d1SBarry Smith   DMDAVecGetArrayWrite - Returns a multiple dimension array that shares data with
14512b4a537SBarry Smith   the underlying vector and is indexed using the global or local dimensions of a `DMDA`.
146fdc842d1SBarry Smith 
14720f4b53cSBarry Smith   Logically Collective
148fdc842d1SBarry Smith 
149d8d19677SJose E. Roman   Input Parameters:
15072fd0fbdSBarry Smith + da  - the `DMDA`
1510af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
152fdc842d1SBarry Smith 
153fdc842d1SBarry Smith   Output Parameter:
154fdc842d1SBarry Smith . array - the array
155fdc842d1SBarry Smith 
156dce8aebaSBarry Smith   Level: intermediate
157dce8aebaSBarry Smith 
158fdc842d1SBarry Smith   Notes:
159dce8aebaSBarry Smith   Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
160fdc842d1SBarry Smith 
161fdc842d1SBarry Smith   In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
162fdc842d1SBarry Smith 
1630af9b551SBarry Smith   if `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
164fdc842d1SBarry Smith   a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
165dce8aebaSBarry Smith   appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
166fdc842d1SBarry Smith 
1670af9b551SBarry Smith   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
1680af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
1690af9b551SBarry Smith 
170fdc842d1SBarry Smith   Fortran Notes:
171*ce78bad3SBarry Smith   Use `DMDAVecGetArrayWrite()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
172dce8aebaSBarry 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
1730f99b6f4SBarry Smith   dimension of the `DMDA`.
174fdc842d1SBarry Smith 
1750af9b551SBarry 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
1760af9b551SBarry Smith   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
1770af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
178fdc842d1SBarry Smith 
17912b4a537SBarry Smith   Developer Note:
180dce8aebaSBarry Smith   This has code duplication with `DMDAVecGetArray()` and `DMDAVecGetArrayRead()`
181fdc842d1SBarry Smith 
18212b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecRestoreArrayDOF()`
183db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
184fdc842d1SBarry Smith @*/
185d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayWrite(DM da, Vec vec, void *array)
186d71ae5a4SJacob Faibussowitsch {
187fdc842d1SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
188fdc842d1SBarry Smith 
189fdc842d1SBarry Smith   PetscFunctionBegin;
190fdc842d1SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
191fdc842d1SBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
1924f572ea9SToby Isaac   PetscAssertPointer(array, 3);
1931bb6d2a8SBarry Smith   if (da->localSection) {
1949566063dSJacob Faibussowitsch     PetscCall(VecGetArrayWrite(vec, (PetscScalar **)array));
1953ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
196fdc842d1SBarry Smith   }
1979566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
1989566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
1999566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
200fdc842d1SBarry Smith 
201fdc842d1SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
2029566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
203fdc842d1SBarry Smith   if (N == xm * ym * zm * dof) {
204fdc842d1SBarry Smith     gxm = xm;
205fdc842d1SBarry Smith     gym = ym;
206fdc842d1SBarry Smith     gzm = zm;
207fdc842d1SBarry Smith     gxs = xs;
208fdc842d1SBarry Smith     gys = ys;
209fdc842d1SBarry Smith     gzs = zs;
21063a3b9bcSJacob 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);
211fdc842d1SBarry Smith 
212fdc842d1SBarry Smith   if (dim == 1) {
2139566063dSJacob Faibussowitsch     PetscCall(VecGetArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
214fdc842d1SBarry Smith   } else if (dim == 2) {
2159566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
216fdc842d1SBarry Smith   } else if (dim == 3) {
2179566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
21863a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
2193ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
220fdc842d1SBarry Smith }
221fdc842d1SBarry Smith 
2220f99b6f4SBarry Smith /*MC
22312b4a537SBarry Smith   DMDAVecRestoreArrayWriteF90 - see Fortran Notes at `DMDAVecRestoreArrayWrite()`
22453c0d4aeSBarry Smith 
22553c0d4aeSBarry Smith   Level: intermediate
2260f99b6f4SBarry Smith M*/
2270f99b6f4SBarry Smith 
228fdc842d1SBarry Smith /*@
229dce8aebaSBarry Smith   DMDAVecRestoreArrayWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayWrite()`
230fdc842d1SBarry Smith 
23120f4b53cSBarry Smith   Logically Collective
232fdc842d1SBarry Smith 
233d8d19677SJose E. Roman   Input Parameters:
23472fd0fbdSBarry Smith + da    - the `DMDA`
2350af9b551SBarry Smith . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
23612b4a537SBarry Smith - array - the `array` pointer
237fdc842d1SBarry Smith 
238fdc842d1SBarry Smith   Level: intermediate
239fdc842d1SBarry Smith 
24012b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayWrite()`,
241db781477SPatrick Sanan           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
242fdc842d1SBarry Smith @*/
243d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayWrite(DM da, Vec vec, void *array)
244d71ae5a4SJacob Faibussowitsch {
245fdc842d1SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
246fdc842d1SBarry Smith 
247fdc842d1SBarry Smith   PetscFunctionBegin;
248fdc842d1SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
249fdc842d1SBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
2504f572ea9SToby Isaac   PetscAssertPointer(array, 3);
2511bb6d2a8SBarry Smith   if (da->localSection) {
2529566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(vec, (PetscScalar **)array));
2533ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
254fdc842d1SBarry Smith   }
2559566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
2569566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
2579566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
258fdc842d1SBarry Smith 
259fdc842d1SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
2609566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
261fdc842d1SBarry Smith   if (N == xm * ym * zm * dof) {
262fdc842d1SBarry Smith     gxm = xm;
263fdc842d1SBarry Smith     gym = ym;
264fdc842d1SBarry Smith     gzm = zm;
265fdc842d1SBarry Smith     gxs = xs;
266fdc842d1SBarry Smith     gys = ys;
267fdc842d1SBarry Smith     gzs = zs;
26863a3b9bcSJacob 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);
269fdc842d1SBarry Smith 
270fdc842d1SBarry Smith   if (dim == 1) {
2719566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
272fdc842d1SBarry Smith   } else if (dim == 2) {
2739566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
274fdc842d1SBarry Smith   } else if (dim == 3) {
2759566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
27663a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
2773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
278fdc842d1SBarry Smith }
279fdc842d1SBarry Smith 
280fdc842d1SBarry Smith /*@C
281aa219208SBarry Smith   DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
28212b4a537SBarry Smith   the underlying vector and is indexed using the global or local dimensions of a `DMDA`
28347c6ae99SBarry Smith 
28420f4b53cSBarry Smith   Logically Collective
28547c6ae99SBarry Smith 
286d8d19677SJose E. Roman   Input Parameters:
28772fd0fbdSBarry Smith + da  - the `DMDA`
2880af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
28947c6ae99SBarry Smith 
29047c6ae99SBarry Smith   Output Parameter:
29112b4a537SBarry Smith . array - the `array` pointer
29247c6ae99SBarry Smith 
293dce8aebaSBarry Smith   Level: intermediate
294dce8aebaSBarry Smith 
29547c6ae99SBarry Smith   Notes:
296dce8aebaSBarry Smith   Call `DMDAVecRestoreArrayDOF()` once you have finished accessing the vector entries.
29747c6ae99SBarry Smith 
2980af9b551SBarry Smith   In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]
2990af9b551SBarry Smith 
3000af9b551SBarry 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
3010af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
30247c6ae99SBarry Smith 
3030f99b6f4SBarry Smith   Fortran Notes:
304*ce78bad3SBarry Smith   Use `DMDAVecGetArray()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
3050f99b6f4SBarry 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
3060f99b6f4SBarry Smith   dimension of the `DMDA`.
3070f99b6f4SBarry Smith 
3080af9b551SBarry 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
3090af9b551SBarry Smith   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
3100af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
3111b82215eSBarry Smith 
31212b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecRestoreArrayDOF()`,
3134ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, `DMDAVecGetArrayDOFRead()`
31447c6ae99SBarry Smith @*/
315d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOF(DM da, Vec vec, void *array)
316d71ae5a4SJacob Faibussowitsch {
31747c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
31847c6ae99SBarry Smith 
31947c6ae99SBarry Smith   PetscFunctionBegin;
3209566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
3219566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
3229566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
32347c6ae99SBarry Smith 
32447c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
3259566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
32647c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
32747c6ae99SBarry Smith     gxm = xm;
32847c6ae99SBarry Smith     gym = ym;
32947c6ae99SBarry Smith     gzm = zm;
33047c6ae99SBarry Smith     gxs = xs;
33147c6ae99SBarry Smith     gys = ys;
33247c6ae99SBarry Smith     gzs = zs;
33363a3b9bcSJacob 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);
33447c6ae99SBarry Smith 
33547c6ae99SBarry Smith   if (dim == 1) {
3369566063dSJacob Faibussowitsch     PetscCall(VecGetArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
33747c6ae99SBarry Smith   } else if (dim == 2) {
3389566063dSJacob Faibussowitsch     PetscCall(VecGetArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
33947c6ae99SBarry Smith   } else if (dim == 3) {
3409566063dSJacob Faibussowitsch     PetscCall(VecGetArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
34163a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
3423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34347c6ae99SBarry Smith }
34447c6ae99SBarry Smith 
34547c6ae99SBarry Smith /*@
346dce8aebaSBarry Smith   DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOF()`
34747c6ae99SBarry Smith 
34820f4b53cSBarry Smith   Logically Collective
34947c6ae99SBarry Smith 
350d8d19677SJose E. Roman   Input Parameters:
35172fd0fbdSBarry Smith + da    - the `DMDA`
3520af9b551SBarry Smith . vec   - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
35312b4a537SBarry Smith - array - the `array` point
35447c6ae99SBarry Smith 
35547c6ae99SBarry Smith   Level: intermediate
35647c6ae99SBarry Smith 
35712b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
3584ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
35947c6ae99SBarry Smith @*/
360d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOF(DM da, Vec vec, void *array)
361d71ae5a4SJacob Faibussowitsch {
36247c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
36347c6ae99SBarry Smith 
36447c6ae99SBarry Smith   PetscFunctionBegin;
3659566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
3669566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
3679566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
36847c6ae99SBarry Smith 
36947c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
3709566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
37147c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
37247c6ae99SBarry Smith     gxm = xm;
37347c6ae99SBarry Smith     gym = ym;
37447c6ae99SBarry Smith     gzm = zm;
37547c6ae99SBarry Smith     gxs = xs;
37647c6ae99SBarry Smith     gys = ys;
37747c6ae99SBarry Smith     gzs = zs;
37847c6ae99SBarry Smith   }
37947c6ae99SBarry Smith 
38047c6ae99SBarry Smith   if (dim == 1) {
3819566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
38247c6ae99SBarry Smith   } else if (dim == 2) {
3839566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
38447c6ae99SBarry Smith   } else if (dim == 3) {
3859566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
38663a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
3873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
38847c6ae99SBarry Smith }
38947c6ae99SBarry Smith 
3900f99b6f4SBarry Smith /*MC
39112b4a537SBarry Smith   DMDAVecGetArrayReadF90 - see Fortran Notes at `DMDAVecGetArrayRead()`
39253c0d4aeSBarry Smith 
39353c0d4aeSBarry Smith   Level: intermediate
3940f99b6f4SBarry Smith M*/
3950f99b6f4SBarry Smith 
3965edff71fSBarry Smith /*@C
3975edff71fSBarry Smith   DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
39812b4a537SBarry Smith   the underlying vector and is indexed using the global or local dimensions of a `DMDA`.
3995edff71fSBarry Smith 
40020f4b53cSBarry Smith   Not Collective
4015edff71fSBarry Smith 
402d8d19677SJose E. Roman   Input Parameters:
40372fd0fbdSBarry Smith + da  - the `DMDA`
4040af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
4055edff71fSBarry Smith 
4065edff71fSBarry Smith   Output Parameter:
4075edff71fSBarry Smith . array - the array
4085edff71fSBarry Smith 
409dce8aebaSBarry Smith   Level: intermediate
410dce8aebaSBarry Smith 
4115edff71fSBarry Smith   Notes:
412dce8aebaSBarry Smith   Call `DMDAVecRestoreArrayRead()` once you have finished accessing the vector entries.
4135edff71fSBarry Smith 
4145edff71fSBarry Smith   In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
4155edff71fSBarry Smith 
4160af9b551SBarry Smith   If `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
4175edff71fSBarry Smith   a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
418dce8aebaSBarry Smith   appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
4195edff71fSBarry Smith 
4200af9b551SBarry Smith   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
4210af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
4220af9b551SBarry Smith 
42395452b02SPatrick Sanan   Fortran Notes:
424*ce78bad3SBarry Smith   Use `DMDAVecGetArrayRead()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
425dce8aebaSBarry 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
4260f99b6f4SBarry Smith   dimension of the `DMDA`.
4270f99b6f4SBarry Smith 
4280af9b551SBarry 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
4290af9b551SBarry Smith   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
4300af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
4315edff71fSBarry Smith 
432*ce78bad3SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`,
43342747ad1SJacob Faibussowitsch `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayRead()`,
43442747ad1SJacob Faibussowitsch `DMDAVecRestoreArrayDOF()`, `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`,
43542747ad1SJacob Faibussowitsch `DMDAVecRestoreArray()`, `DMStagVecGetArrayRead()`
4365edff71fSBarry Smith @*/
437d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayRead(DM da, Vec vec, void *array)
438d71ae5a4SJacob Faibussowitsch {
4395edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
4405edff71fSBarry Smith 
4415edff71fSBarry Smith   PetscFunctionBegin;
442a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
4435edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
4444f572ea9SToby Isaac   PetscAssertPointer(array, 3);
4459566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
4469566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
4479566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
4485edff71fSBarry Smith 
4495edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
4509566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
4515edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
4525edff71fSBarry Smith     gxm = xm;
4535edff71fSBarry Smith     gym = ym;
4545edff71fSBarry Smith     gzm = zm;
4555edff71fSBarry Smith     gxs = xs;
4565edff71fSBarry Smith     gys = ys;
4575edff71fSBarry Smith     gzs = zs;
45863a3b9bcSJacob 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);
4595edff71fSBarry Smith 
4605edff71fSBarry Smith   if (dim == 1) {
4619566063dSJacob Faibussowitsch     PetscCall(VecGetArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
4625edff71fSBarry Smith   } else if (dim == 2) {
4639566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
4645edff71fSBarry Smith   } else if (dim == 3) {
4659566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
46663a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
4673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4685edff71fSBarry Smith }
4695edff71fSBarry Smith 
4700f99b6f4SBarry Smith /*MC
47112b4a537SBarry Smith   DMDAVecRestoreArrayReadF90 - see Fortran Notes at `DMDAVecRestoreArrayRead()`
47253c0d4aeSBarry Smith 
47353c0d4aeSBarry Smith   Level: intermediate
4740f99b6f4SBarry Smith M*/
4750f99b6f4SBarry Smith 
4765edff71fSBarry Smith /*@
477dce8aebaSBarry Smith   DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayRead()`
4785edff71fSBarry Smith 
47920f4b53cSBarry Smith   Not Collective
4805edff71fSBarry Smith 
481d8d19677SJose E. Roman   Input Parameters:
48272fd0fbdSBarry Smith + da    - the `DMDA`
4830af9b551SBarry Smith . vec   - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
48412b4a537SBarry Smith - array - the `array` pointer
4855edff71fSBarry Smith 
4865edff71fSBarry Smith   Level: intermediate
4875edff71fSBarry Smith 
488*ce78bad3SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayRead()`,
489db781477SPatrick Sanan           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`,
490db781477SPatrick Sanan           `DMStagVecRestoreArrayRead()`
4915edff71fSBarry Smith @*/
492d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayRead(DM da, Vec vec, void *array)
493d71ae5a4SJacob Faibussowitsch {
4945edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
4955edff71fSBarry Smith 
4965edff71fSBarry Smith   PetscFunctionBegin;
497a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
4985edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
4994f572ea9SToby Isaac   PetscAssertPointer(array, 3);
5009566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
5019566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
5029566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5035edff71fSBarry Smith 
5045edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
5059566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
5065edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
5075edff71fSBarry Smith     gxm = xm;
5085edff71fSBarry Smith     gym = ym;
5095edff71fSBarry Smith     gzm = zm;
5105edff71fSBarry Smith     gxs = xs;
5115edff71fSBarry Smith     gys = ys;
5125edff71fSBarry Smith     gzs = zs;
51363a3b9bcSJacob 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);
5145edff71fSBarry Smith 
5155edff71fSBarry Smith   if (dim == 1) {
5169566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
5175edff71fSBarry Smith   } else if (dim == 2) {
5189566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
5195edff71fSBarry Smith   } else if (dim == 3) {
5209566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
52163a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
5223ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5235edff71fSBarry Smith }
5245edff71fSBarry Smith 
5255edff71fSBarry Smith /*@C
5265edff71fSBarry Smith   DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
52712b4a537SBarry Smith   the underlying vector and is indexed using the global or local dimensions of a `DMDA`
5285edff71fSBarry Smith 
5295edff71fSBarry Smith   Not Collective
5305edff71fSBarry Smith 
531d8d19677SJose E. Roman   Input Parameters:
53272fd0fbdSBarry Smith + da  - the `DMDA`
5330af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
5345edff71fSBarry Smith 
5355edff71fSBarry Smith   Output Parameter:
5365edff71fSBarry Smith . array - the array
5375edff71fSBarry Smith 
538dce8aebaSBarry Smith   Level: intermediate
539dce8aebaSBarry Smith 
5405edff71fSBarry Smith   Notes:
541dce8aebaSBarry Smith   Call `DMDAVecRestoreArrayDOFRead()` once you have finished accessing the vector entries.
5425edff71fSBarry Smith 
5435edff71fSBarry Smith   In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
5445edff71fSBarry Smith 
5450af9b551SBarry Smith   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
5460af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
5470af9b551SBarry Smith 
54860225df5SJacob Faibussowitsch   Fortran Notes:
549*ce78bad3SBarry Smith   Use  `DMDAVecGetArrayRead()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
5500f99b6f4SBarry 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
5510f99b6f4SBarry Smith   dimension of the `DMDA`.
5520f99b6f4SBarry Smith 
5530af9b551SBarry 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
5540af9b551SBarry Smith   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
5550af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
5565edff71fSBarry Smith 
55712b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
5584ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
5595edff71fSBarry Smith @*/
560d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOFRead(DM da, Vec vec, void *array)
561d71ae5a4SJacob Faibussowitsch {
5625edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
5635edff71fSBarry Smith 
5645edff71fSBarry Smith   PetscFunctionBegin;
5659566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
5669566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
5679566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5685edff71fSBarry Smith 
5695edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
5709566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
5715edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
5725edff71fSBarry Smith     gxm = xm;
5735edff71fSBarry Smith     gym = ym;
5745edff71fSBarry Smith     gzm = zm;
5755edff71fSBarry Smith     gxs = xs;
5765edff71fSBarry Smith     gys = ys;
5775edff71fSBarry Smith     gzs = zs;
57863a3b9bcSJacob 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);
5795edff71fSBarry Smith 
5805edff71fSBarry Smith   if (dim == 1) {
5819566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
5825edff71fSBarry Smith   } else if (dim == 2) {
5839566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
5845edff71fSBarry Smith   } else if (dim == 3) {
5859566063dSJacob Faibussowitsch     PetscCall(VecGetArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
58663a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
5873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5885edff71fSBarry Smith }
5895edff71fSBarry Smith 
5905edff71fSBarry Smith /*@
591dce8aebaSBarry Smith   DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFRead()`
5925edff71fSBarry Smith 
5935edff71fSBarry Smith   Not Collective
5945edff71fSBarry Smith 
595d8d19677SJose E. Roman   Input Parameters:
59672fd0fbdSBarry Smith + da    - the `DMDA`
5970af9b551SBarry Smith . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
59812b4a537SBarry Smith - array - the `array` pointer
5995edff71fSBarry Smith 
6005edff71fSBarry Smith   Level: intermediate
6015edff71fSBarry Smith 
60212b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
6034ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
6045edff71fSBarry Smith @*/
605d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOFRead(DM da, Vec vec, void *array)
606d71ae5a4SJacob Faibussowitsch {
6075edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
6085edff71fSBarry Smith 
6095edff71fSBarry Smith   PetscFunctionBegin;
6109566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
6119566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
6129566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
6135edff71fSBarry Smith 
6145edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
6159566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
6165edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
6175edff71fSBarry Smith     gxm = xm;
6185edff71fSBarry Smith     gym = ym;
6195edff71fSBarry Smith     gzm = zm;
6205edff71fSBarry Smith     gxs = xs;
6215edff71fSBarry Smith     gys = ys;
6225edff71fSBarry Smith     gzs = zs;
6235edff71fSBarry Smith   }
6245edff71fSBarry Smith 
6255edff71fSBarry Smith   if (dim == 1) {
6269566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
6275edff71fSBarry Smith   } else if (dim == 2) {
6289566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
6295edff71fSBarry Smith   } else if (dim == 3) {
6309566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
63163a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
6323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6335edff71fSBarry Smith }
6345edff71fSBarry Smith 
6351e5d2365SBarry Smith /*@C
6361e5d2365SBarry Smith   DMDAVecGetArrayDOFWrite - Returns a multiple dimension array that shares data with
63712b4a537SBarry Smith   the underlying vector and is indexed using the global or local dimensions of a `DMDA`
6381e5d2365SBarry Smith 
6391e5d2365SBarry Smith   Not Collective
6401e5d2365SBarry Smith 
641d8d19677SJose E. Roman   Input Parameters:
64272fd0fbdSBarry Smith + da  - the `DMDA`
6430af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
6441e5d2365SBarry Smith 
6451e5d2365SBarry Smith   Output Parameter:
6461e5d2365SBarry Smith . array - the array
6471e5d2365SBarry Smith 
6480f99b6f4SBarry Smith   Level: intermediate
6490f99b6f4SBarry Smith 
6501e5d2365SBarry Smith   Notes:
651dce8aebaSBarry Smith   Call `DMDAVecRestoreArrayDOFWrite()` once you have finished accessing the vector entries.
6521e5d2365SBarry Smith 
6531e5d2365SBarry Smith   In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
6541e5d2365SBarry Smith 
6550af9b551SBarry 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
6560af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
6570af9b551SBarry Smith 
65860225df5SJacob Faibussowitsch   Fortran Notes:
659*ce78bad3SBarry Smith   Use  `DMDAVecGetArrayWrite()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
6600f99b6f4SBarry 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
6610f99b6f4SBarry Smith   dimension of the `DMDA`.
6621e5d2365SBarry Smith 
6630af9b551SBarry 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
6640af9b551SBarry Smith   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
6650af9b551SBarry Smith   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
6661e5d2365SBarry Smith 
667*ce78bad3SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
66860225df5SJacob Faibussowitsch           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
6691e5d2365SBarry Smith @*/
670d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOFWrite(DM da, Vec vec, void *array)
671d71ae5a4SJacob Faibussowitsch {
6721e5d2365SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
6731e5d2365SBarry Smith 
6741e5d2365SBarry Smith   PetscFunctionBegin;
6759566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
6769566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
6779566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
6781e5d2365SBarry Smith 
6791e5d2365SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
6809566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
6811e5d2365SBarry Smith   if (N == xm * ym * zm * dof) {
6821e5d2365SBarry Smith     gxm = xm;
6831e5d2365SBarry Smith     gym = ym;
6841e5d2365SBarry Smith     gzm = zm;
6851e5d2365SBarry Smith     gxs = xs;
6861e5d2365SBarry Smith     gys = ys;
6871e5d2365SBarry Smith     gzs = zs;
68863a3b9bcSJacob 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);
6891e5d2365SBarry Smith 
6901e5d2365SBarry Smith   if (dim == 1) {
6919566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
6921e5d2365SBarry Smith   } else if (dim == 2) {
6939566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
6941e5d2365SBarry Smith   } else if (dim == 3) {
6959566063dSJacob Faibussowitsch     PetscCall(VecGetArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
69663a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
6973ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6981e5d2365SBarry Smith }
6991e5d2365SBarry Smith 
7001e5d2365SBarry Smith /*@
701dce8aebaSBarry Smith   DMDAVecRestoreArrayDOFWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFWrite()`
7021e5d2365SBarry Smith 
7031e5d2365SBarry Smith   Not Collective
7041e5d2365SBarry Smith 
705d8d19677SJose E. Roman   Input Parameters:
70672fd0fbdSBarry Smith + da    - the `DMDA`
7070af9b551SBarry Smith . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
70812b4a537SBarry Smith - array - the `array` pointer
7091e5d2365SBarry Smith 
7101e5d2365SBarry Smith   Level: intermediate
7111e5d2365SBarry Smith 
71212b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
71360225df5SJacob Faibussowitsch           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
7141e5d2365SBarry Smith @*/
715d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM da, Vec vec, void *array)
716d71ae5a4SJacob Faibussowitsch {
7171e5d2365SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
7181e5d2365SBarry Smith 
7191e5d2365SBarry Smith   PetscFunctionBegin;
7209566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
7219566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
7229566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
7231e5d2365SBarry Smith 
7241e5d2365SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
7259566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
7261e5d2365SBarry Smith   if (N == xm * ym * zm * dof) {
7271e5d2365SBarry Smith     gxm = xm;
7281e5d2365SBarry Smith     gym = ym;
7291e5d2365SBarry Smith     gzm = zm;
7301e5d2365SBarry Smith     gxs = xs;
7311e5d2365SBarry Smith     gys = ys;
7321e5d2365SBarry Smith     gzs = zs;
7331e5d2365SBarry Smith   }
7341e5d2365SBarry Smith 
7351e5d2365SBarry Smith   if (dim == 1) {
7369566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
7371e5d2365SBarry Smith   } else if (dim == 2) {
7389566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
7391e5d2365SBarry Smith   } else if (dim == 3) {
7409566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
74163a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
7423ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7431e5d2365SBarry Smith }
744