xref: /petsc/src/dm/impls/da/dagetarray.c (revision dce8aeba1c9b69b19f651c53d8a6b674bd7e9cbd)
147c6ae99SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I   "petscdmda.h"   I*/
347c6ae99SBarry Smith 
447c6ae99SBarry Smith /*@C
5aa219208SBarry Smith    DMDAVecGetArray - Returns a multiple dimension array that shares data with
647c6ae99SBarry Smith       the underlying vector and is indexed using the global dimensions.
747c6ae99SBarry Smith 
8d083f849SBarry Smith    Logically collective on da
947c6ae99SBarry Smith 
10d8d19677SJose E. Roman    Input Parameters:
1147c6ae99SBarry Smith +  da - the distributed array
12*dce8aebaSBarry Smith -  vec - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
1347c6ae99SBarry Smith 
1447c6ae99SBarry Smith    Output Parameter:
1547c6ae99SBarry Smith .  array - the array
1647c6ae99SBarry Smith 
17*dce8aebaSBarry Smith   Level: intermediate
18*dce8aebaSBarry Smith 
1947c6ae99SBarry Smith    Notes:
20*dce8aebaSBarry Smith     Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
2147c6ae99SBarry Smith 
2247c6ae99SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
2347c6ae99SBarry Smith 
247e57d48aSBarry Smith     If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
257e57d48aSBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
2647c6ae99SBarry Smith 
27*dce8aebaSBarry Smith     appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
2847c6ae99SBarry Smith 
2995452b02SPatrick Sanan   Fortran Notes:
30*dce8aebaSBarry Smith     From Fortran use `DMDAVecGetArrayF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
31*dce8aebaSBarry 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
32*dce8aebaSBarry Smith        dimension of the `DMDA`. The order of the indices is array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) (when dof is 1) otherwise
33bef27881SBarry Smith        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
34af0996ceSBarry Smith        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
3547c6ae99SBarry Smith 
36*dce8aebaSBarry Smith   Due to bugs in the compiler `DMDAVecGetArrayF90()` does not work with gfortran versions before 4.5
3747c6ae99SBarry Smith 
38*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
39db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
40db781477SPatrick Sanan           `DMStagVecGetArray()`
4147c6ae99SBarry Smith @*/
42d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArray(DM da, Vec vec, void *array)
43d71ae5a4SJacob Faibussowitsch {
4447c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
4547c6ae99SBarry Smith 
4647c6ae99SBarry Smith   PetscFunctionBegin;
47a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
486db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
496db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
509566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
519566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
529566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5347c6ae99SBarry Smith 
5447c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
559566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
5647c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
5747c6ae99SBarry Smith     gxm = xm;
5847c6ae99SBarry Smith     gym = ym;
5947c6ae99SBarry Smith     gzm = zm;
6047c6ae99SBarry Smith     gxs = xs;
6147c6ae99SBarry Smith     gys = ys;
6247c6ae99SBarry Smith     gzs = zs;
6363a3b9bcSJacob 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);
6447c6ae99SBarry Smith 
6547c6ae99SBarry Smith   if (dim == 1) {
669566063dSJacob Faibussowitsch     PetscCall(VecGetArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
6747c6ae99SBarry Smith   } else if (dim == 2) {
689566063dSJacob Faibussowitsch     PetscCall(VecGetArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
6947c6ae99SBarry Smith   } else if (dim == 3) {
709566063dSJacob Faibussowitsch     PetscCall(VecGetArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
7163a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
7247c6ae99SBarry Smith   PetscFunctionReturn(0);
7347c6ae99SBarry Smith }
7447c6ae99SBarry Smith 
7547c6ae99SBarry Smith /*@
76*dce8aebaSBarry Smith    DMDAVecRestoreArray - Restores a multiple dimension array obtained with `DMDAVecGetArray()`
7747c6ae99SBarry Smith 
78d083f849SBarry Smith    Logically collective on da
7947c6ae99SBarry Smith 
80d8d19677SJose E. Roman    Input Parameters:
8147c6ae99SBarry Smith +  da - the distributed array
8247c6ae99SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
83564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
84e5c84f05SJed Brown -  array - the array, non-NULL pointer is zeroed
8547c6ae99SBarry Smith 
8647c6ae99SBarry Smith   Level: intermediate
8747c6ae99SBarry Smith 
88*dce8aebaSBarry Smith   Fortran Note:
89*dce8aebaSBarry Smith     From Fortran use `DMDAVecRestoreArayF90()`
9047c6ae99SBarry Smith 
91*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`,
92db781477SPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
93db781477SPatrick Sanan           `DMDStagVecRestoreArray()`
9447c6ae99SBarry Smith @*/
95d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArray(DM da, Vec vec, void *array)
96d71ae5a4SJacob Faibussowitsch {
9747c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
9847c6ae99SBarry Smith 
9947c6ae99SBarry Smith   PetscFunctionBegin;
100a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
1016db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
1026db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
1039566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
1049566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
1059566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
10647c6ae99SBarry Smith 
10747c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
1089566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
10947c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
11047c6ae99SBarry Smith     gxm = xm;
11147c6ae99SBarry Smith     gym = ym;
11247c6ae99SBarry Smith     gzm = zm;
11347c6ae99SBarry Smith     gxs = xs;
11447c6ae99SBarry Smith     gys = ys;
11547c6ae99SBarry Smith     gzs = zs;
11663a3b9bcSJacob 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);
11747c6ae99SBarry Smith 
11847c6ae99SBarry Smith   if (dim == 1) {
1199566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
12047c6ae99SBarry Smith   } else if (dim == 2) {
1219566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
12247c6ae99SBarry Smith   } else if (dim == 3) {
1239566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
12463a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
12547c6ae99SBarry Smith   PetscFunctionReturn(0);
12647c6ae99SBarry Smith }
12747c6ae99SBarry Smith 
12847c6ae99SBarry Smith /*@C
129fdc842d1SBarry Smith    DMDAVecGetArrayWrite - Returns a multiple dimension array that shares data with
130fdc842d1SBarry Smith       the underlying vector and is indexed using the global dimensions.
131fdc842d1SBarry Smith 
132fdc842d1SBarry Smith    Logically collective on Vec
133fdc842d1SBarry Smith 
134d8d19677SJose E. Roman    Input Parameters:
135fdc842d1SBarry Smith +  da - the distributed array
136*dce8aebaSBarry Smith -  vec - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
137fdc842d1SBarry Smith 
138fdc842d1SBarry Smith    Output Parameter:
139fdc842d1SBarry Smith .  array - the array
140fdc842d1SBarry Smith 
141*dce8aebaSBarry Smith   Level: intermediate
142*dce8aebaSBarry Smith 
143fdc842d1SBarry Smith    Notes:
144*dce8aebaSBarry Smith     Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
145fdc842d1SBarry Smith 
146fdc842d1SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
147fdc842d1SBarry Smith 
148*dce8aebaSBarry Smith     If vec is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
149fdc842d1SBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
150fdc842d1SBarry Smith 
151*dce8aebaSBarry Smith     appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
152fdc842d1SBarry Smith 
153fdc842d1SBarry Smith   Fortran Notes:
154*dce8aebaSBarry Smith     From Fortran use `DMDAVecGetArrayF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
155*dce8aebaSBarry 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
156*dce8aebaSBarry Smith        dimension of the `DMDA`. The order of the indices is array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) (when dof is 1) otherwise
157fdc842d1SBarry Smith        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
158*dce8aebaSBarry Smith        DMDAGetCorners() for a global array or `DMDAGetGhostCorners()` for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
159fdc842d1SBarry Smith 
160*dce8aebaSBarry Smith   Due to bugs in the compiler `DMDAVecGetArrayF90()` does not work with gfortran versions before 4.5
161fdc842d1SBarry Smith 
162*dce8aebaSBarry Smith   Developer Note:
163*dce8aebaSBarry Smith   This has code duplication with `DMDAVecGetArray()` and `DMDAVecGetArrayRead()`
164fdc842d1SBarry Smith 
165*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecRestoreArrayDOF()`
166db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
167fdc842d1SBarry Smith @*/
168d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayWrite(DM da, Vec vec, void *array)
169d71ae5a4SJacob Faibussowitsch {
170fdc842d1SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
171fdc842d1SBarry Smith 
172fdc842d1SBarry Smith   PetscFunctionBegin;
173fdc842d1SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
174fdc842d1SBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
175fdc842d1SBarry Smith   PetscValidPointer(array, 3);
1761bb6d2a8SBarry Smith   if (da->localSection) {
1779566063dSJacob Faibussowitsch     PetscCall(VecGetArrayWrite(vec, (PetscScalar **)array));
178fdc842d1SBarry Smith     PetscFunctionReturn(0);
179fdc842d1SBarry Smith   }
1809566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
1819566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
1829566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
183fdc842d1SBarry Smith 
184fdc842d1SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
1859566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
186fdc842d1SBarry Smith   if (N == xm * ym * zm * dof) {
187fdc842d1SBarry Smith     gxm = xm;
188fdc842d1SBarry Smith     gym = ym;
189fdc842d1SBarry Smith     gzm = zm;
190fdc842d1SBarry Smith     gxs = xs;
191fdc842d1SBarry Smith     gys = ys;
192fdc842d1SBarry Smith     gzs = zs;
19363a3b9bcSJacob 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);
194fdc842d1SBarry Smith 
195fdc842d1SBarry Smith   if (dim == 1) {
1969566063dSJacob Faibussowitsch     PetscCall(VecGetArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
197fdc842d1SBarry Smith   } else if (dim == 2) {
1989566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
199fdc842d1SBarry Smith   } else if (dim == 3) {
2009566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
20163a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
202fdc842d1SBarry Smith   PetscFunctionReturn(0);
203fdc842d1SBarry Smith }
204fdc842d1SBarry Smith 
205fdc842d1SBarry Smith /*@
206*dce8aebaSBarry Smith    DMDAVecRestoreArrayWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayWrite()`
207fdc842d1SBarry Smith 
208*dce8aebaSBarry Smith    Logically collective on vec
209fdc842d1SBarry Smith 
210d8d19677SJose E. Roman    Input Parameters:
211fdc842d1SBarry Smith +  da - the distributed array
212fdc842d1SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
213*dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
214fdc842d1SBarry Smith -  array - the array, non-NULL pointer is zeroed
215fdc842d1SBarry Smith 
216fdc842d1SBarry Smith   Level: intermediate
217fdc842d1SBarry Smith 
218*dce8aebaSBarry Smith   Fortran Note:
219*dce8aebaSBarry Smith     From Fortran use `DMDAVecRestoreArayF90()`
220fdc842d1SBarry Smith 
221*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayWrite()`,
222db781477SPatrick Sanan           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
223fdc842d1SBarry Smith @*/
224d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayWrite(DM da, Vec vec, void *array)
225d71ae5a4SJacob Faibussowitsch {
226fdc842d1SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
227fdc842d1SBarry Smith 
228fdc842d1SBarry Smith   PetscFunctionBegin;
229fdc842d1SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
230fdc842d1SBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
231fdc842d1SBarry Smith   PetscValidPointer(array, 3);
2321bb6d2a8SBarry Smith   if (da->localSection) {
2339566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(vec, (PetscScalar **)array));
234fdc842d1SBarry Smith     PetscFunctionReturn(0);
235fdc842d1SBarry Smith   }
2369566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
2379566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
2389566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
239fdc842d1SBarry Smith 
240fdc842d1SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
2419566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
242fdc842d1SBarry Smith   if (N == xm * ym * zm * dof) {
243fdc842d1SBarry Smith     gxm = xm;
244fdc842d1SBarry Smith     gym = ym;
245fdc842d1SBarry Smith     gzm = zm;
246fdc842d1SBarry Smith     gxs = xs;
247fdc842d1SBarry Smith     gys = ys;
248fdc842d1SBarry Smith     gzs = zs;
24963a3b9bcSJacob 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);
250fdc842d1SBarry Smith 
251fdc842d1SBarry Smith   if (dim == 1) {
2529566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
253fdc842d1SBarry Smith   } else if (dim == 2) {
2549566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
255fdc842d1SBarry Smith   } else if (dim == 3) {
2569566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
25763a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
258fdc842d1SBarry Smith   PetscFunctionReturn(0);
259fdc842d1SBarry Smith }
260fdc842d1SBarry Smith 
261fdc842d1SBarry Smith /*@C
262aa219208SBarry Smith    DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
26347c6ae99SBarry Smith       the underlying vector and is indexed using the global dimensions.
26447c6ae99SBarry Smith 
2652ddbf755SMatthew G. Knepley    Logically collective
26647c6ae99SBarry Smith 
267d8d19677SJose E. Roman    Input Parameters:
26847c6ae99SBarry Smith +  da - the distributed array
26947c6ae99SBarry Smith -  vec - the vector, either a vector the same size as one obtained with
270*dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
27147c6ae99SBarry Smith 
27247c6ae99SBarry Smith    Output Parameter:
27347c6ae99SBarry Smith .  array - the array
27447c6ae99SBarry Smith 
275*dce8aebaSBarry Smith   Level: intermediate
276*dce8aebaSBarry Smith 
27747c6ae99SBarry Smith    Notes:
278*dce8aebaSBarry Smith     Call `DMDAVecRestoreArrayDOF()` once you have finished accessing the vector entries.
27947c6ae99SBarry Smith 
28047c6ae99SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
28147c6ae99SBarry Smith 
282*dce8aebaSBarry Smith     In Fortran 90 you do not need a version of `DMDAVecRestoreArrayDOF()` just use `DMDAVecRestoreArrayF90()` and declare your array with one higher dimension,
283c4762a1bSJed Brown     see src/dm/tutorials/ex11f90.F
2841b82215eSBarry Smith 
285*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecRestoreArrayDOF()`,
2864ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, `DMDAVecGetArrayDOFRead()`
28747c6ae99SBarry Smith @*/
288d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOF(DM da, Vec vec, void *array)
289d71ae5a4SJacob Faibussowitsch {
29047c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
29147c6ae99SBarry Smith 
29247c6ae99SBarry Smith   PetscFunctionBegin;
2939566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
2949566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
2959566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
29647c6ae99SBarry Smith 
29747c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
2989566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
29947c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
30047c6ae99SBarry Smith     gxm = xm;
30147c6ae99SBarry Smith     gym = ym;
30247c6ae99SBarry Smith     gzm = zm;
30347c6ae99SBarry Smith     gxs = xs;
30447c6ae99SBarry Smith     gys = ys;
30547c6ae99SBarry Smith     gzs = zs;
30663a3b9bcSJacob 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);
30747c6ae99SBarry Smith 
30847c6ae99SBarry Smith   if (dim == 1) {
3099566063dSJacob Faibussowitsch     PetscCall(VecGetArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
31047c6ae99SBarry Smith   } else if (dim == 2) {
3119566063dSJacob Faibussowitsch     PetscCall(VecGetArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
31247c6ae99SBarry Smith   } else if (dim == 3) {
3139566063dSJacob Faibussowitsch     PetscCall(VecGetArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
31463a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
31547c6ae99SBarry Smith   PetscFunctionReturn(0);
31647c6ae99SBarry Smith }
31747c6ae99SBarry Smith 
31847c6ae99SBarry Smith /*@
319*dce8aebaSBarry Smith    DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOF()`
32047c6ae99SBarry Smith 
3212ddbf755SMatthew G. Knepley    Logically collective
32247c6ae99SBarry Smith 
323d8d19677SJose E. Roman    Input Parameters:
32447c6ae99SBarry Smith +  da - the distributed array
32547c6ae99SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
326*dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
32747c6ae99SBarry Smith -  array - the array
32847c6ae99SBarry Smith 
32947c6ae99SBarry Smith   Level: intermediate
33047c6ae99SBarry Smith 
331*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
3324ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
33347c6ae99SBarry Smith @*/
334d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOF(DM da, Vec vec, void *array)
335d71ae5a4SJacob Faibussowitsch {
33647c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
33747c6ae99SBarry Smith 
33847c6ae99SBarry Smith   PetscFunctionBegin;
3399566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
3409566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
3419566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
34247c6ae99SBarry Smith 
34347c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
3449566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
34547c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
34647c6ae99SBarry Smith     gxm = xm;
34747c6ae99SBarry Smith     gym = ym;
34847c6ae99SBarry Smith     gzm = zm;
34947c6ae99SBarry Smith     gxs = xs;
35047c6ae99SBarry Smith     gys = ys;
35147c6ae99SBarry Smith     gzs = zs;
35247c6ae99SBarry Smith   }
35347c6ae99SBarry Smith 
35447c6ae99SBarry Smith   if (dim == 1) {
3559566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
35647c6ae99SBarry Smith   } else if (dim == 2) {
3579566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
35847c6ae99SBarry Smith   } else if (dim == 3) {
3599566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
36063a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
36147c6ae99SBarry Smith   PetscFunctionReturn(0);
36247c6ae99SBarry Smith }
36347c6ae99SBarry Smith 
3645edff71fSBarry Smith /*@C
3655edff71fSBarry Smith    DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
3665edff71fSBarry Smith       the underlying vector and is indexed using the global dimensions.
3675edff71fSBarry Smith 
3682ddbf755SMatthew G. Knepley    Not collective
3695edff71fSBarry Smith 
370d8d19677SJose E. Roman    Input Parameters:
3715edff71fSBarry Smith +  da - the distributed array
372*dce8aebaSBarry Smith -  vec - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
3735edff71fSBarry Smith 
3745edff71fSBarry Smith    Output Parameter:
3755edff71fSBarry Smith .  array - the array
3765edff71fSBarry Smith 
377*dce8aebaSBarry Smith   Level: intermediate
378*dce8aebaSBarry Smith 
3795edff71fSBarry Smith    Notes:
380*dce8aebaSBarry Smith     Call `DMDAVecRestoreArrayRead()` once you have finished accessing the vector entries.
3815edff71fSBarry Smith 
3825edff71fSBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
3835edff71fSBarry Smith 
384*dce8aebaSBarry Smith     If vec is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
3855edff71fSBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
386*dce8aebaSBarry Smith     appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
3875edff71fSBarry Smith 
38895452b02SPatrick Sanan   Fortran Notes:
389*dce8aebaSBarry Smith     From Fortran use `DMDAVecGetArrayReadF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
390*dce8aebaSBarry 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
391*dce8aebaSBarry Smith        dimension of the `DMDA`. The order of the indices is array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) (when dof is 1) otherwise
3925edff71fSBarry Smith        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
393*dce8aebaSBarry Smith        DMDAGetCorners() for a global array or `DMDAGetGhostCorners()` for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
3945edff71fSBarry Smith 
395*dce8aebaSBarry Smith   Due to bugs in the compiler `DMDAVecGetArrayReadF90()` does not work with gfortran versions before 4.5
3965edff71fSBarry Smith 
397*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayRead()`, `DMDAVecRestoreArrayDOF()`
398db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
399db781477SPatrick Sanan           `DMStagVecGetArrayRead()`
4005edff71fSBarry Smith @*/
401d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayRead(DM da, Vec vec, void *array)
402d71ae5a4SJacob Faibussowitsch {
4035edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
4045edff71fSBarry Smith 
4055edff71fSBarry Smith   PetscFunctionBegin;
406a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
4075edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
4085edff71fSBarry Smith   PetscValidPointer(array, 3);
4099566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
4109566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
4119566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
4125edff71fSBarry Smith 
4135edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
4149566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
4155edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
4165edff71fSBarry Smith     gxm = xm;
4175edff71fSBarry Smith     gym = ym;
4185edff71fSBarry Smith     gzm = zm;
4195edff71fSBarry Smith     gxs = xs;
4205edff71fSBarry Smith     gys = ys;
4215edff71fSBarry Smith     gzs = zs;
42263a3b9bcSJacob 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);
4235edff71fSBarry Smith 
4245edff71fSBarry Smith   if (dim == 1) {
4259566063dSJacob Faibussowitsch     PetscCall(VecGetArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
4265edff71fSBarry Smith   } else if (dim == 2) {
4279566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
4285edff71fSBarry Smith   } else if (dim == 3) {
4299566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
43063a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
4315edff71fSBarry Smith   PetscFunctionReturn(0);
4325edff71fSBarry Smith }
4335edff71fSBarry Smith 
4345edff71fSBarry Smith /*@
435*dce8aebaSBarry Smith    DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayRead()`
4365edff71fSBarry Smith 
4372ddbf755SMatthew G. Knepley    Not collective
4385edff71fSBarry Smith 
439d8d19677SJose E. Roman    Input Parameters:
4405edff71fSBarry Smith +  da - the distributed array
4415edff71fSBarry Smith .  vec - the vector, either a vector the same size as one obtained with
442*dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
4435edff71fSBarry Smith -  array - the array, non-NULL pointer is zeroed
4445edff71fSBarry Smith 
4455edff71fSBarry Smith   Level: intermediate
4465edff71fSBarry Smith 
447*dce8aebaSBarry Smith   Fortran Note:
448*dce8aebaSBarry Smith     From Fortran use `DMDAVecRestoreArrayReadF90()`
4495edff71fSBarry Smith 
450*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayRead()`,
451db781477SPatrick Sanan           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`,
452db781477SPatrick Sanan           `DMStagVecRestoreArrayRead()`
4535edff71fSBarry Smith @*/
454d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayRead(DM da, Vec vec, void *array)
455d71ae5a4SJacob Faibussowitsch {
4565edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
4575edff71fSBarry Smith 
4585edff71fSBarry Smith   PetscFunctionBegin;
459a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
4605edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
4615edff71fSBarry Smith   PetscValidPointer(array, 3);
4629566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
4639566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
4649566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
4655edff71fSBarry Smith 
4665edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
4679566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
4685edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
4695edff71fSBarry Smith     gxm = xm;
4705edff71fSBarry Smith     gym = ym;
4715edff71fSBarry Smith     gzm = zm;
4725edff71fSBarry Smith     gxs = xs;
4735edff71fSBarry Smith     gys = ys;
4745edff71fSBarry Smith     gzs = zs;
47563a3b9bcSJacob 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);
4765edff71fSBarry Smith 
4775edff71fSBarry Smith   if (dim == 1) {
4789566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
4795edff71fSBarry Smith   } else if (dim == 2) {
4809566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
4815edff71fSBarry Smith   } else if (dim == 3) {
4829566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
48363a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
4845edff71fSBarry Smith   PetscFunctionReturn(0);
4855edff71fSBarry Smith }
4865edff71fSBarry Smith 
4875edff71fSBarry Smith /*@C
4885edff71fSBarry Smith    DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
4895edff71fSBarry Smith       the underlying vector and is indexed using the global dimensions.
4905edff71fSBarry Smith 
4915edff71fSBarry Smith    Not Collective
4925edff71fSBarry Smith 
493d8d19677SJose E. Roman    Input Parameters:
4945edff71fSBarry Smith +  da - the distributed array
4955edff71fSBarry Smith -  vec - the vector, either a vector the same size as one obtained with
496*dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
4975edff71fSBarry Smith 
4985edff71fSBarry Smith    Output Parameter:
4995edff71fSBarry Smith .  array - the array
5005edff71fSBarry Smith 
501*dce8aebaSBarry Smith   Level: intermediate
502*dce8aebaSBarry Smith 
5035edff71fSBarry Smith    Notes:
504*dce8aebaSBarry Smith     Call `DMDAVecRestoreArrayDOFRead()` once you have finished accessing the vector entries.
5055edff71fSBarry Smith 
5065edff71fSBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
5075edff71fSBarry Smith 
508*dce8aebaSBarry Smith    Fortran Note:
509*dce8aebaSBarry Smith     In Fortran you do not need a version of `DMDAVecRestoreArrayDOF()` just use  `DMDAVecRestoreArrayReadF90()` and declare your array with one higher dimension,
510c4762a1bSJed Brown     see src/dm/tutorials/ex11f90.F
5115edff71fSBarry Smith 
512*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
5134ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
5145edff71fSBarry Smith @*/
515d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOFRead(DM da, Vec vec, void *array)
516d71ae5a4SJacob Faibussowitsch {
5175edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
5185edff71fSBarry Smith 
5195edff71fSBarry Smith   PetscFunctionBegin;
5209566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
5219566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
5229566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5235edff71fSBarry Smith 
5245edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
5259566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
5265edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
5275edff71fSBarry Smith     gxm = xm;
5285edff71fSBarry Smith     gym = ym;
5295edff71fSBarry Smith     gzm = zm;
5305edff71fSBarry Smith     gxs = xs;
5315edff71fSBarry Smith     gys = ys;
5325edff71fSBarry Smith     gzs = zs;
53363a3b9bcSJacob 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);
5345edff71fSBarry Smith 
5355edff71fSBarry Smith   if (dim == 1) {
5369566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
5375edff71fSBarry Smith   } else if (dim == 2) {
5389566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
5395edff71fSBarry Smith   } else if (dim == 3) {
5409566063dSJacob Faibussowitsch     PetscCall(VecGetArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
54163a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
5425edff71fSBarry Smith   PetscFunctionReturn(0);
5435edff71fSBarry Smith }
5445edff71fSBarry Smith 
5455edff71fSBarry Smith /*@
546*dce8aebaSBarry Smith    DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFRead()`
5475edff71fSBarry Smith 
5485edff71fSBarry Smith    Not Collective
5495edff71fSBarry Smith 
550d8d19677SJose E. Roman    Input Parameters:
5515edff71fSBarry Smith +  da - the distributed array
5525edff71fSBarry Smith .  vec - the vector, either a vector the same size as one obtained with
553*dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
5545edff71fSBarry Smith -  array - the array
5555edff71fSBarry Smith 
5565edff71fSBarry Smith   Level: intermediate
5575edff71fSBarry Smith 
558*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
5594ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
5605edff71fSBarry Smith @*/
561d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOFRead(DM da, Vec vec, void *array)
562d71ae5a4SJacob Faibussowitsch {
5635edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
5645edff71fSBarry Smith 
5655edff71fSBarry Smith   PetscFunctionBegin;
5669566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
5679566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
5689566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5695edff71fSBarry Smith 
5705edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
5719566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
5725edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
5735edff71fSBarry Smith     gxm = xm;
5745edff71fSBarry Smith     gym = ym;
5755edff71fSBarry Smith     gzm = zm;
5765edff71fSBarry Smith     gxs = xs;
5775edff71fSBarry Smith     gys = ys;
5785edff71fSBarry Smith     gzs = zs;
5795edff71fSBarry Smith   }
5805edff71fSBarry Smith 
5815edff71fSBarry Smith   if (dim == 1) {
5829566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
5835edff71fSBarry Smith   } else if (dim == 2) {
5849566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
5855edff71fSBarry Smith   } else if (dim == 3) {
5869566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
58763a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
5885edff71fSBarry Smith   PetscFunctionReturn(0);
5895edff71fSBarry Smith }
5905edff71fSBarry Smith 
5911e5d2365SBarry Smith /*@C
5921e5d2365SBarry Smith    DMDAVecGetArrayDOFWrite - Returns a multiple dimension array that shares data with
5931e5d2365SBarry Smith       the underlying vector and is indexed using the global dimensions.
5941e5d2365SBarry Smith 
5951e5d2365SBarry Smith    Not Collective
5961e5d2365SBarry Smith 
597d8d19677SJose E. Roman    Input Parameters:
5981e5d2365SBarry Smith +  da - the distributed array
5991e5d2365SBarry Smith -  vec - the vector, either a vector the same size as one obtained with
600*dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
6011e5d2365SBarry Smith 
6021e5d2365SBarry Smith    Output Parameter:
6031e5d2365SBarry Smith .  array - the array
6041e5d2365SBarry Smith 
6051e5d2365SBarry Smith    Notes:
606*dce8aebaSBarry Smith     Call `DMDAVecRestoreArrayDOFWrite()` once you have finished accessing the vector entries.
6071e5d2365SBarry Smith 
6081e5d2365SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
6091e5d2365SBarry Smith 
610*dce8aebaSBarry Smith    Fortran Noe:
611*dce8aebaSBarry Smith     In Fortran you do not need a version of `DMDAVecRestoreArrayDOF()` just use  `DMDAVecRestoreArrayWriteF90()` and declare your array with one higher dimension,
6121e5d2365SBarry Smith     see src/dm/tutorials/ex11f90.F
6131e5d2365SBarry Smith 
6141e5d2365SBarry Smith   Level: intermediate
6151e5d2365SBarry Smith 
616*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
6174ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
6181e5d2365SBarry Smith @*/
619d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOFWrite(DM da, Vec vec, void *array)
620d71ae5a4SJacob Faibussowitsch {
6211e5d2365SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
6221e5d2365SBarry Smith 
6231e5d2365SBarry Smith   PetscFunctionBegin;
6249566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
6259566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
6269566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
6271e5d2365SBarry Smith 
6281e5d2365SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
6299566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
6301e5d2365SBarry Smith   if (N == xm * ym * zm * dof) {
6311e5d2365SBarry Smith     gxm = xm;
6321e5d2365SBarry Smith     gym = ym;
6331e5d2365SBarry Smith     gzm = zm;
6341e5d2365SBarry Smith     gxs = xs;
6351e5d2365SBarry Smith     gys = ys;
6361e5d2365SBarry Smith     gzs = zs;
63763a3b9bcSJacob 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);
6381e5d2365SBarry Smith 
6391e5d2365SBarry Smith   if (dim == 1) {
6409566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
6411e5d2365SBarry Smith   } else if (dim == 2) {
6429566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
6431e5d2365SBarry Smith   } else if (dim == 3) {
6449566063dSJacob Faibussowitsch     PetscCall(VecGetArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
64563a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
6461e5d2365SBarry Smith   PetscFunctionReturn(0);
6471e5d2365SBarry Smith }
6481e5d2365SBarry Smith 
6491e5d2365SBarry Smith /*@
650*dce8aebaSBarry Smith    DMDAVecRestoreArrayDOFWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFWrite()`
6511e5d2365SBarry Smith 
6521e5d2365SBarry Smith    Not Collective
6531e5d2365SBarry Smith 
654d8d19677SJose E. Roman    Input Parameters:
6551e5d2365SBarry Smith +  da - the distributed array
6561e5d2365SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
657*dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
6581e5d2365SBarry Smith -  array - the array
6591e5d2365SBarry Smith 
6601e5d2365SBarry Smith   Level: intermediate
6611e5d2365SBarry Smith 
662*dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
6634ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
6641e5d2365SBarry Smith @*/
665d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM da, Vec vec, void *array)
666d71ae5a4SJacob Faibussowitsch {
6671e5d2365SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
6681e5d2365SBarry Smith 
6691e5d2365SBarry Smith   PetscFunctionBegin;
6709566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
6719566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
6729566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
6731e5d2365SBarry Smith 
6741e5d2365SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
6759566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
6761e5d2365SBarry Smith   if (N == xm * ym * zm * dof) {
6771e5d2365SBarry Smith     gxm = xm;
6781e5d2365SBarry Smith     gym = ym;
6791e5d2365SBarry Smith     gzm = zm;
6801e5d2365SBarry Smith     gxs = xs;
6811e5d2365SBarry Smith     gys = ys;
6821e5d2365SBarry Smith     gzs = zs;
6831e5d2365SBarry Smith   }
6841e5d2365SBarry Smith 
6851e5d2365SBarry Smith   if (dim == 1) {
6869566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
6871e5d2365SBarry Smith   } else if (dim == 2) {
6889566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
6891e5d2365SBarry Smith   } else if (dim == 3) {
6909566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
69163a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
6921e5d2365SBarry Smith   PetscFunctionReturn(0);
6931e5d2365SBarry Smith }
694