xref: /petsc/src/dm/impls/da/dagetarray.c (revision 20f4b53cbb5e9bd9ef12b76a8697d60d197cda17)
147c6ae99SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I   "petscdmda.h"   I*/
347c6ae99SBarry Smith 
40f99b6f4SBarry Smith /*MC
50f99b6f4SBarry Smith   DMDAVecGetArrayF90 - check Fortran Notes at `DMDAVecGetArray()`
671b416bcSMatthew G. Knepley 
771b416bcSMatthew G. Knepley   Level: intermediate
80f99b6f4SBarry Smith M*/
90f99b6f4SBarry Smith 
1047c6ae99SBarry Smith /*@C
11aa219208SBarry Smith    DMDAVecGetArray - Returns a multiple dimension array that shares data with
1247c6ae99SBarry Smith       the underlying vector and is indexed using the global dimensions.
1347c6ae99SBarry Smith 
14*20f4b53cSBarry Smith    Logically Collective
1547c6ae99SBarry Smith 
16d8d19677SJose E. Roman    Input Parameters:
1747c6ae99SBarry Smith +  da - the distributed array
18dce8aebaSBarry Smith -  vec - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
1947c6ae99SBarry Smith 
2047c6ae99SBarry Smith    Output Parameter:
2147c6ae99SBarry Smith .  array - the array
2247c6ae99SBarry Smith 
23dce8aebaSBarry Smith   Level: intermediate
24dce8aebaSBarry Smith 
2547c6ae99SBarry Smith    Notes:
26dce8aebaSBarry Smith     Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
2747c6ae99SBarry Smith 
2847c6ae99SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
2947c6ae99SBarry Smith 
307e57d48aSBarry Smith     If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
317e57d48aSBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
3247c6ae99SBarry Smith 
33dce8aebaSBarry Smith     appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
3447c6ae99SBarry Smith 
3595452b02SPatrick Sanan   Fortran Notes:
360f99b6f4SBarry Smith   Use `DMDAVecGetArrayF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
37dce8aebaSBarry 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
380f99b6f4SBarry Smith   dimension of the `DMDA`.
3947c6ae99SBarry Smith 
400f99b6f4SBarry 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
410f99b6f4SBarry Smith   array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
420f99b6f4SBarry Smith   `DMDAGetCorners()` for a global array or `DMDAGetGhostCorners()` for a local array.
4347c6ae99SBarry Smith 
44dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
45db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
46db781477SPatrick Sanan           `DMStagVecGetArray()`
4747c6ae99SBarry Smith @*/
48d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArray(DM da, Vec vec, void *array)
49d71ae5a4SJacob Faibussowitsch {
5047c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
5147c6ae99SBarry Smith 
5247c6ae99SBarry Smith   PetscFunctionBegin;
53a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
546db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
556db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
569566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
579566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
589566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5947c6ae99SBarry Smith 
6047c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
619566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
6247c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
6347c6ae99SBarry Smith     gxm = xm;
6447c6ae99SBarry Smith     gym = ym;
6547c6ae99SBarry Smith     gzm = zm;
6647c6ae99SBarry Smith     gxs = xs;
6747c6ae99SBarry Smith     gys = ys;
6847c6ae99SBarry Smith     gzs = zs;
6963a3b9bcSJacob 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);
7047c6ae99SBarry Smith 
7147c6ae99SBarry Smith   if (dim == 1) {
729566063dSJacob Faibussowitsch     PetscCall(VecGetArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
7347c6ae99SBarry Smith   } else if (dim == 2) {
749566063dSJacob Faibussowitsch     PetscCall(VecGetArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
7547c6ae99SBarry Smith   } else if (dim == 3) {
769566063dSJacob Faibussowitsch     PetscCall(VecGetArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
7763a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7947c6ae99SBarry Smith }
8047c6ae99SBarry Smith 
810f99b6f4SBarry Smith /*MC
820f99b6f4SBarry Smith   DMDAVecRestoreArrayF90 - check Fortran Notes at `DMDAVecRestoreArray()`
8353c0d4aeSBarry Smith 
8453c0d4aeSBarry Smith   Level: intermediate
850f99b6f4SBarry Smith M*/
860f99b6f4SBarry Smith 
8747c6ae99SBarry Smith /*@
88dce8aebaSBarry Smith    DMDAVecRestoreArray - Restores a multiple dimension array obtained with `DMDAVecGetArray()`
8947c6ae99SBarry Smith 
90*20f4b53cSBarry Smith    Logically Collective
9147c6ae99SBarry Smith 
92d8d19677SJose E. Roman    Input Parameters:
9347c6ae99SBarry Smith +  da - the distributed array
9447c6ae99SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
950f99b6f4SBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
96e5c84f05SJed Brown -  array - the array, non-NULL pointer is zeroed
9747c6ae99SBarry Smith 
9847c6ae99SBarry Smith   Level: intermediate
9947c6ae99SBarry Smith 
100dce8aebaSBarry Smith   Fortran Note:
101dce8aebaSBarry Smith   From Fortran use `DMDAVecRestoreArayF90()`
10247c6ae99SBarry Smith 
1030f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `DMDAVecRestoreArayF90()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`,
104db781477SPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
105db781477SPatrick Sanan           `DMDStagVecRestoreArray()`
10647c6ae99SBarry Smith @*/
107d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArray(DM da, Vec vec, void *array)
108d71ae5a4SJacob Faibussowitsch {
10947c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
11047c6ae99SBarry Smith 
11147c6ae99SBarry Smith   PetscFunctionBegin;
112a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
1136db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
1146db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
1159566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
1169566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
1179566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
11847c6ae99SBarry Smith 
11947c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
1209566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
12147c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
12247c6ae99SBarry Smith     gxm = xm;
12347c6ae99SBarry Smith     gym = ym;
12447c6ae99SBarry Smith     gzm = zm;
12547c6ae99SBarry Smith     gxs = xs;
12647c6ae99SBarry Smith     gys = ys;
12747c6ae99SBarry Smith     gzs = zs;
12863a3b9bcSJacob 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);
12947c6ae99SBarry Smith 
13047c6ae99SBarry Smith   if (dim == 1) {
1319566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
13247c6ae99SBarry Smith   } else if (dim == 2) {
1339566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
13447c6ae99SBarry Smith   } else if (dim == 3) {
1359566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
13663a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
1373ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13847c6ae99SBarry Smith }
13947c6ae99SBarry Smith 
1400f99b6f4SBarry Smith /*MC
1410f99b6f4SBarry Smith   DMDAVecGetArrayWriteF90 - check Fortran Notes at `DMDAVecGetArrayWrite()`
14253c0d4aeSBarry Smith 
14353c0d4aeSBarry Smith   Level: intermediate
1440f99b6f4SBarry Smith M*/
1450f99b6f4SBarry Smith 
14647c6ae99SBarry Smith /*@C
147fdc842d1SBarry Smith    DMDAVecGetArrayWrite - Returns a multiple dimension array that shares data with
148fdc842d1SBarry Smith       the underlying vector and is indexed using the global dimensions.
149fdc842d1SBarry Smith 
150*20f4b53cSBarry Smith    Logically Collective
151fdc842d1SBarry Smith 
152d8d19677SJose E. Roman    Input Parameters:
153fdc842d1SBarry Smith +  da - the distributed array
154dce8aebaSBarry Smith -  vec - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
155fdc842d1SBarry Smith 
156fdc842d1SBarry Smith    Output Parameter:
157fdc842d1SBarry Smith .  array - the array
158fdc842d1SBarry Smith 
159dce8aebaSBarry Smith   Level: intermediate
160dce8aebaSBarry Smith 
161fdc842d1SBarry Smith    Notes:
162dce8aebaSBarry Smith     Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
163fdc842d1SBarry Smith 
164fdc842d1SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
165fdc842d1SBarry Smith 
166dce8aebaSBarry Smith     If vec is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
167fdc842d1SBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
168dce8aebaSBarry Smith     appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
169fdc842d1SBarry Smith 
170fdc842d1SBarry Smith    Fortran Notes:
1710f99b6f4SBarry Smith    From Fortran use `DMDAVecGetArrayWriteF90()` 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 
1750f99b6f4SBarry 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
1760f99b6f4SBarry Smith    array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
1770f99b6f4SBarry Smith    `DMDAGetCorners()` for a global array or `DMDAGetGhostCorners()` for a local array.
178fdc842d1SBarry Smith 
179dce8aebaSBarry Smith   Developer Note:
180dce8aebaSBarry Smith   This has code duplication with `DMDAVecGetArray()` and `DMDAVecGetArrayRead()`
181fdc842d1SBarry Smith 
182dce8aebaSBarry Smith .seealso: `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);
192fdc842d1SBarry Smith   PetscValidPointer(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
2230f99b6f4SBarry Smith   DMDAVecRestoreArrayWriteF90 - check 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 
231*20f4b53cSBarry Smith    Logically Collective
232fdc842d1SBarry Smith 
233d8d19677SJose E. Roman    Input Parameters:
234fdc842d1SBarry Smith +  da - the distributed array
235fdc842d1SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
236dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
237fdc842d1SBarry Smith -  array - the array, non-NULL pointer is zeroed
238fdc842d1SBarry Smith 
239fdc842d1SBarry Smith   Level: intermediate
240fdc842d1SBarry Smith 
241dce8aebaSBarry Smith   Fortran Note:
2420f99b6f4SBarry Smith   Use `DMDAVecRestoreArayWriteF90()`
243fdc842d1SBarry Smith 
244dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayWrite()`,
245db781477SPatrick Sanan           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
246fdc842d1SBarry Smith @*/
247d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayWrite(DM da, Vec vec, void *array)
248d71ae5a4SJacob Faibussowitsch {
249fdc842d1SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
250fdc842d1SBarry Smith 
251fdc842d1SBarry Smith   PetscFunctionBegin;
252fdc842d1SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
253fdc842d1SBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
254fdc842d1SBarry Smith   PetscValidPointer(array, 3);
2551bb6d2a8SBarry Smith   if (da->localSection) {
2569566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(vec, (PetscScalar **)array));
2573ba16761SJacob Faibussowitsch     PetscFunctionReturn(PETSC_SUCCESS);
258fdc842d1SBarry Smith   }
2599566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
2609566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
2619566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
262fdc842d1SBarry Smith 
263fdc842d1SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
2649566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
265fdc842d1SBarry Smith   if (N == xm * ym * zm * dof) {
266fdc842d1SBarry Smith     gxm = xm;
267fdc842d1SBarry Smith     gym = ym;
268fdc842d1SBarry Smith     gzm = zm;
269fdc842d1SBarry Smith     gxs = xs;
270fdc842d1SBarry Smith     gys = ys;
271fdc842d1SBarry Smith     gzs = zs;
27263a3b9bcSJacob 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);
273fdc842d1SBarry Smith 
274fdc842d1SBarry Smith   if (dim == 1) {
2759566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
276fdc842d1SBarry Smith   } else if (dim == 2) {
2779566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
278fdc842d1SBarry Smith   } else if (dim == 3) {
2799566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
28063a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
2813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
282fdc842d1SBarry Smith }
283fdc842d1SBarry Smith 
284fdc842d1SBarry Smith /*@C
285aa219208SBarry Smith    DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
28647c6ae99SBarry Smith       the underlying vector and is indexed using the global dimensions.
28747c6ae99SBarry Smith 
288*20f4b53cSBarry Smith    Logically Collective
28947c6ae99SBarry Smith 
290d8d19677SJose E. Roman    Input Parameters:
29147c6ae99SBarry Smith +  da - the distributed array
29247c6ae99SBarry Smith -  vec - the vector, either a vector the same size as one obtained with
293dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
29447c6ae99SBarry Smith 
29547c6ae99SBarry Smith    Output Parameter:
29647c6ae99SBarry Smith .  array - the array
29747c6ae99SBarry Smith 
298dce8aebaSBarry Smith   Level: intermediate
299dce8aebaSBarry Smith 
30047c6ae99SBarry Smith    Notes:
301dce8aebaSBarry Smith     Call `DMDAVecRestoreArrayDOF()` once you have finished accessing the vector entries.
30247c6ae99SBarry Smith 
30347c6ae99SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
30447c6ae99SBarry Smith 
3050f99b6f4SBarry Smith    Fortran Notes:
3060f99b6f4SBarry Smith     Use `DMDAVecGetArrayF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
3070f99b6f4SBarry 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
3080f99b6f4SBarry Smith    dimension of the `DMDA`.
3090f99b6f4SBarry Smith 
3100f99b6f4SBarry 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
3110f99b6f4SBarry Smith    array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
3120f99b6f4SBarry Smith    `DMDAGetCorners()` for a global array or `DMDAGetGhostCorners()` for a local array.
3131b82215eSBarry Smith 
314dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecRestoreArrayDOF()`,
3154ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, `DMDAVecGetArrayDOFRead()`
31647c6ae99SBarry Smith @*/
317d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOF(DM da, Vec vec, void *array)
318d71ae5a4SJacob Faibussowitsch {
31947c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
32047c6ae99SBarry Smith 
32147c6ae99SBarry Smith   PetscFunctionBegin;
3229566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
3239566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
3249566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
32547c6ae99SBarry Smith 
32647c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
3279566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
32847c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
32947c6ae99SBarry Smith     gxm = xm;
33047c6ae99SBarry Smith     gym = ym;
33147c6ae99SBarry Smith     gzm = zm;
33247c6ae99SBarry Smith     gxs = xs;
33347c6ae99SBarry Smith     gys = ys;
33447c6ae99SBarry Smith     gzs = zs;
33563a3b9bcSJacob 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);
33647c6ae99SBarry Smith 
33747c6ae99SBarry Smith   if (dim == 1) {
3389566063dSJacob Faibussowitsch     PetscCall(VecGetArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
33947c6ae99SBarry Smith   } else if (dim == 2) {
3409566063dSJacob Faibussowitsch     PetscCall(VecGetArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
34147c6ae99SBarry Smith   } else if (dim == 3) {
3429566063dSJacob Faibussowitsch     PetscCall(VecGetArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
34363a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
3443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
34547c6ae99SBarry Smith }
34647c6ae99SBarry Smith 
34747c6ae99SBarry Smith /*@
348dce8aebaSBarry Smith    DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOF()`
34947c6ae99SBarry Smith 
350*20f4b53cSBarry Smith    Logically Collective
35147c6ae99SBarry Smith 
352d8d19677SJose E. Roman    Input Parameters:
35347c6ae99SBarry Smith +  da - the distributed array
35447c6ae99SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
355dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
35647c6ae99SBarry Smith -  array - the array
35747c6ae99SBarry Smith 
35847c6ae99SBarry Smith   Level: intermediate
35947c6ae99SBarry Smith 
3600f99b6f4SBarry Smith   Fortran Note:
3610f99b6f4SBarry Smith   Use `DMDAVecRestoreArayF90()`
3620f99b6f4SBarry Smith 
363dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
3644ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
36547c6ae99SBarry Smith @*/
366d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOF(DM da, Vec vec, void *array)
367d71ae5a4SJacob Faibussowitsch {
36847c6ae99SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
36947c6ae99SBarry Smith 
37047c6ae99SBarry Smith   PetscFunctionBegin;
3719566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
3729566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
3739566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
37447c6ae99SBarry Smith 
37547c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
3769566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
37747c6ae99SBarry Smith   if (N == xm * ym * zm * dof) {
37847c6ae99SBarry Smith     gxm = xm;
37947c6ae99SBarry Smith     gym = ym;
38047c6ae99SBarry Smith     gzm = zm;
38147c6ae99SBarry Smith     gxs = xs;
38247c6ae99SBarry Smith     gys = ys;
38347c6ae99SBarry Smith     gzs = zs;
38447c6ae99SBarry Smith   }
38547c6ae99SBarry Smith 
38647c6ae99SBarry Smith   if (dim == 1) {
3879566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
38847c6ae99SBarry Smith   } else if (dim == 2) {
3899566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
39047c6ae99SBarry Smith   } else if (dim == 3) {
3919566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
39263a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
3933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
39447c6ae99SBarry Smith }
39547c6ae99SBarry Smith 
3960f99b6f4SBarry Smith /*MC
3970f99b6f4SBarry Smith   DMDAVecGetArrayReadF90 - check Fortran Notes at `DMDAVecGetArrayRead()`
39853c0d4aeSBarry Smith 
39953c0d4aeSBarry Smith   Level: intermediate
4000f99b6f4SBarry Smith M*/
4010f99b6f4SBarry Smith 
4025edff71fSBarry Smith /*@C
4035edff71fSBarry Smith    DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
4045edff71fSBarry Smith       the underlying vector and is indexed using the global dimensions.
4055edff71fSBarry Smith 
406*20f4b53cSBarry Smith    Not Collective
4075edff71fSBarry Smith 
408d8d19677SJose E. Roman    Input Parameters:
4095edff71fSBarry Smith +  da - the distributed array
410dce8aebaSBarry Smith -  vec - the vector, either a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
4115edff71fSBarry Smith 
4125edff71fSBarry Smith    Output Parameter:
4135edff71fSBarry Smith .  array - the array
4145edff71fSBarry Smith 
415dce8aebaSBarry Smith   Level: intermediate
416dce8aebaSBarry Smith 
4175edff71fSBarry Smith    Notes:
418dce8aebaSBarry Smith     Call `DMDAVecRestoreArrayRead()` once you have finished accessing the vector entries.
4195edff71fSBarry Smith 
4205edff71fSBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
4215edff71fSBarry Smith 
422dce8aebaSBarry Smith     If vec is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
4235edff71fSBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
424dce8aebaSBarry Smith     appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
4255edff71fSBarry Smith 
42695452b02SPatrick Sanan   Fortran Notes:
4270f99b6f4SBarry Smith   Use `DMDAVecGetArrayReadF90()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
428dce8aebaSBarry 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
4290f99b6f4SBarry Smith   dimension of the `DMDA`.
4300f99b6f4SBarry Smith 
4310f99b6f4SBarry 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
4325edff71fSBarry Smith   array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
4330f99b6f4SBarry Smith   `DMDAGetCorners()` for a global array or `DMDAGetGhostCorners()` for a local array.
4345edff71fSBarry Smith 
4350f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAVecGetArrayReadF90()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayRead()`, `DMDAVecRestoreArrayDOF()`
436db781477SPatrick Sanan           `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
437db781477SPatrick Sanan           `DMStagVecGetArrayRead()`
4385edff71fSBarry Smith @*/
439d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayRead(DM da, Vec vec, void *array)
440d71ae5a4SJacob Faibussowitsch {
4415edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
4425edff71fSBarry Smith 
4435edff71fSBarry Smith   PetscFunctionBegin;
444a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
4455edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
4465edff71fSBarry Smith   PetscValidPointer(array, 3);
4479566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
4489566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
4499566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
4505edff71fSBarry Smith 
4515edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
4529566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
4535edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
4545edff71fSBarry Smith     gxm = xm;
4555edff71fSBarry Smith     gym = ym;
4565edff71fSBarry Smith     gzm = zm;
4575edff71fSBarry Smith     gxs = xs;
4585edff71fSBarry Smith     gys = ys;
4595edff71fSBarry Smith     gzs = zs;
46063a3b9bcSJacob 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);
4615edff71fSBarry Smith 
4625edff71fSBarry Smith   if (dim == 1) {
4639566063dSJacob Faibussowitsch     PetscCall(VecGetArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
4645edff71fSBarry Smith   } else if (dim == 2) {
4659566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
4665edff71fSBarry Smith   } else if (dim == 3) {
4679566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
46863a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
4693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4705edff71fSBarry Smith }
4715edff71fSBarry Smith 
4720f99b6f4SBarry Smith /*MC
4730f99b6f4SBarry Smith   DMDAVecRestoreArrayReadF90 - check Fortran Notes at `DMDAVecRestoreArrayRead()`
47453c0d4aeSBarry Smith 
47553c0d4aeSBarry Smith   Level: intermediate
4760f99b6f4SBarry Smith M*/
4770f99b6f4SBarry Smith 
4785edff71fSBarry Smith /*@
479dce8aebaSBarry Smith    DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayRead()`
4805edff71fSBarry Smith 
481*20f4b53cSBarry Smith    Not Collective
4825edff71fSBarry Smith 
483d8d19677SJose E. Roman    Input Parameters:
4845edff71fSBarry Smith +  da - the distributed array
4855edff71fSBarry Smith .  vec - the vector, either a vector the same size as one obtained with
486dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
4875edff71fSBarry Smith -  array - the array, non-NULL pointer is zeroed
4885edff71fSBarry Smith 
4895edff71fSBarry Smith   Level: intermediate
4905edff71fSBarry Smith 
491dce8aebaSBarry Smith   Fortran Note:
4920f99b6f4SBarry Smith   Use `DMDAVecRestoreArrayReadF90()`
4935edff71fSBarry Smith 
4940f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAVecRestoreArrayReadF90()`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayRead()`,
495db781477SPatrick Sanan           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`,
496db781477SPatrick Sanan           `DMStagVecRestoreArrayRead()`
4975edff71fSBarry Smith @*/
498d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayRead(DM da, Vec vec, void *array)
499d71ae5a4SJacob Faibussowitsch {
5005edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
5015edff71fSBarry Smith 
5025edff71fSBarry Smith   PetscFunctionBegin;
503a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
5045edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
5055edff71fSBarry Smith   PetscValidPointer(array, 3);
5069566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
5079566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
5089566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5095edff71fSBarry Smith 
5105edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
5119566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
5125edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
5135edff71fSBarry Smith     gxm = xm;
5145edff71fSBarry Smith     gym = ym;
5155edff71fSBarry Smith     gzm = zm;
5165edff71fSBarry Smith     gxs = xs;
5175edff71fSBarry Smith     gys = ys;
5185edff71fSBarry Smith     gzs = zs;
51963a3b9bcSJacob 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);
5205edff71fSBarry Smith 
5215edff71fSBarry Smith   if (dim == 1) {
5229566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
5235edff71fSBarry Smith   } else if (dim == 2) {
5249566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
5255edff71fSBarry Smith   } else if (dim == 3) {
5269566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
52763a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
5283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5295edff71fSBarry Smith }
5305edff71fSBarry Smith 
5315edff71fSBarry Smith /*@C
5325edff71fSBarry Smith    DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
5335edff71fSBarry Smith       the underlying vector and is indexed using the global dimensions.
5345edff71fSBarry Smith 
5355edff71fSBarry Smith    Not Collective
5365edff71fSBarry Smith 
537d8d19677SJose E. Roman    Input Parameters:
5385edff71fSBarry Smith +  da - the distributed array
5395edff71fSBarry Smith -  vec - the vector, either a vector the same size as one obtained with
540dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
5415edff71fSBarry Smith 
5425edff71fSBarry Smith    Output Parameter:
5435edff71fSBarry Smith .  array - the array
5445edff71fSBarry Smith 
545dce8aebaSBarry Smith   Level: intermediate
546dce8aebaSBarry Smith 
5475edff71fSBarry Smith    Notes:
548dce8aebaSBarry Smith    Call `DMDAVecRestoreArrayDOFRead()` once you have finished accessing the vector entries.
5495edff71fSBarry Smith 
5505edff71fSBarry Smith    In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
5515edff71fSBarry Smith 
552dce8aebaSBarry Smith    Fortran Note:
5530f99b6f4SBarry Smith    Use  `DMDAVecGetArrayReadF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
5540f99b6f4SBarry 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
5550f99b6f4SBarry Smith    dimension of the `DMDA`.
5560f99b6f4SBarry Smith 
5570f99b6f4SBarry 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
5580f99b6f4SBarry Smith    array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
5590f99b6f4SBarry Smith    `DMDAGetCorners()` for a global array or `DMDAGetGhostCorners()` for a local array.
5605edff71fSBarry Smith 
561dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
5624ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
5635edff71fSBarry Smith @*/
564d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOFRead(DM da, Vec vec, void *array)
565d71ae5a4SJacob Faibussowitsch {
5665edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
5675edff71fSBarry Smith 
5685edff71fSBarry Smith   PetscFunctionBegin;
5699566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
5709566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
5719566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5725edff71fSBarry Smith 
5735edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
5749566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
5755edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
5765edff71fSBarry Smith     gxm = xm;
5775edff71fSBarry Smith     gym = ym;
5785edff71fSBarry Smith     gzm = zm;
5795edff71fSBarry Smith     gxs = xs;
5805edff71fSBarry Smith     gys = ys;
5815edff71fSBarry Smith     gzs = zs;
58263a3b9bcSJacob 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);
5835edff71fSBarry Smith 
5845edff71fSBarry Smith   if (dim == 1) {
5859566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
5865edff71fSBarry Smith   } else if (dim == 2) {
5879566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
5885edff71fSBarry Smith   } else if (dim == 3) {
5899566063dSJacob Faibussowitsch     PetscCall(VecGetArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
59063a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
5913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
5925edff71fSBarry Smith }
5935edff71fSBarry Smith 
5945edff71fSBarry Smith /*@
595dce8aebaSBarry Smith    DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFRead()`
5965edff71fSBarry Smith 
5975edff71fSBarry Smith    Not Collective
5985edff71fSBarry Smith 
599d8d19677SJose E. Roman    Input Parameters:
6005edff71fSBarry Smith +  da - the distributed array
6015edff71fSBarry Smith .  vec - the vector, either a vector the same size as one obtained with
602dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
6035edff71fSBarry Smith -  array - the array
6045edff71fSBarry Smith 
6055edff71fSBarry Smith   Level: intermediate
6065edff71fSBarry Smith 
6070f99b6f4SBarry Smith   Fortran Note:
6080f99b6f4SBarry Smith   Use `DMDAVecRestoreArrayReadF90()`
6090f99b6f4SBarry Smith 
610dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
6114ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
6125edff71fSBarry Smith @*/
613d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOFRead(DM da, Vec vec, void *array)
614d71ae5a4SJacob Faibussowitsch {
6155edff71fSBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
6165edff71fSBarry Smith 
6175edff71fSBarry Smith   PetscFunctionBegin;
6189566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
6199566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
6209566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
6215edff71fSBarry Smith 
6225edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
6239566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
6245edff71fSBarry Smith   if (N == xm * ym * zm * dof) {
6255edff71fSBarry Smith     gxm = xm;
6265edff71fSBarry Smith     gym = ym;
6275edff71fSBarry Smith     gzm = zm;
6285edff71fSBarry Smith     gxs = xs;
6295edff71fSBarry Smith     gys = ys;
6305edff71fSBarry Smith     gzs = zs;
6315edff71fSBarry Smith   }
6325edff71fSBarry Smith 
6335edff71fSBarry Smith   if (dim == 1) {
6349566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
6355edff71fSBarry Smith   } else if (dim == 2) {
6369566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
6375edff71fSBarry Smith   } else if (dim == 3) {
6389566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
63963a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
6403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6415edff71fSBarry Smith }
6425edff71fSBarry Smith 
6431e5d2365SBarry Smith /*@C
6441e5d2365SBarry Smith    DMDAVecGetArrayDOFWrite - Returns a multiple dimension array that shares data with
6451e5d2365SBarry Smith       the underlying vector and is indexed using the global dimensions.
6461e5d2365SBarry Smith 
6471e5d2365SBarry Smith    Not Collective
6481e5d2365SBarry Smith 
649d8d19677SJose E. Roman    Input Parameters:
6501e5d2365SBarry Smith +  da - the distributed array
6511e5d2365SBarry Smith -  vec - the vector, either a vector the same size as one obtained with
652dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
6531e5d2365SBarry Smith 
6541e5d2365SBarry Smith    Output Parameter:
6551e5d2365SBarry Smith .  array - the array
6561e5d2365SBarry Smith 
6570f99b6f4SBarry Smith    Level: intermediate
6580f99b6f4SBarry Smith 
6591e5d2365SBarry Smith    Notes:
660dce8aebaSBarry Smith     Call `DMDAVecRestoreArrayDOFWrite()` once you have finished accessing the vector entries.
6611e5d2365SBarry Smith 
6621e5d2365SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
6631e5d2365SBarry Smith 
6640f99b6f4SBarry Smith    Fortran Note:
6650f99b6f4SBarry Smith    Use  `DMDAVecGetArrayWriteF90()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
6660f99b6f4SBarry 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
6670f99b6f4SBarry Smith    dimension of the `DMDA`.
6681e5d2365SBarry Smith 
6690f99b6f4SBarry 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
6700f99b6f4SBarry Smith    array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
6710f99b6f4SBarry Smith    `DMDAGetCorners()` for a global array or `DMDAGetGhostCorners()` for a local array.
6721e5d2365SBarry Smith 
6730f99b6f4SBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `DMDAVecRestoreArrayWriteF90()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
6744ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
6751e5d2365SBarry Smith @*/
676d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOFWrite(DM da, Vec vec, void *array)
677d71ae5a4SJacob Faibussowitsch {
6781e5d2365SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
6791e5d2365SBarry Smith 
6801e5d2365SBarry Smith   PetscFunctionBegin;
6819566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
6829566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
6839566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
6841e5d2365SBarry Smith 
6851e5d2365SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
6869566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
6871e5d2365SBarry Smith   if (N == xm * ym * zm * dof) {
6881e5d2365SBarry Smith     gxm = xm;
6891e5d2365SBarry Smith     gym = ym;
6901e5d2365SBarry Smith     gzm = zm;
6911e5d2365SBarry Smith     gxs = xs;
6921e5d2365SBarry Smith     gys = ys;
6931e5d2365SBarry Smith     gzs = zs;
69463a3b9bcSJacob 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);
6951e5d2365SBarry Smith 
6961e5d2365SBarry Smith   if (dim == 1) {
6979566063dSJacob Faibussowitsch     PetscCall(VecGetArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
6981e5d2365SBarry Smith   } else if (dim == 2) {
6999566063dSJacob Faibussowitsch     PetscCall(VecGetArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
7001e5d2365SBarry Smith   } else if (dim == 3) {
7019566063dSJacob Faibussowitsch     PetscCall(VecGetArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
70263a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
7033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7041e5d2365SBarry Smith }
7051e5d2365SBarry Smith 
7061e5d2365SBarry Smith /*@
707dce8aebaSBarry Smith    DMDAVecRestoreArrayDOFWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFWrite()`
7081e5d2365SBarry Smith 
7091e5d2365SBarry Smith    Not Collective
7101e5d2365SBarry Smith 
711d8d19677SJose E. Roman    Input Parameters:
7121e5d2365SBarry Smith +  da - the distributed array
7131e5d2365SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
714dce8aebaSBarry Smith          `DMCreateGlobalVector()` or `DMCreateLocalVector()`
7151e5d2365SBarry Smith -  array - the array
7161e5d2365SBarry Smith 
7171e5d2365SBarry Smith   Level: intermediate
7181e5d2365SBarry Smith 
7190f99b6f4SBarry Smith   Fortran Note:
7200f99b6f4SBarry Smith   Use `DMDAVecRestoreArrayWriteF90()`
7210f99b6f4SBarry Smith 
722dce8aebaSBarry Smith .seealso: `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
7234ab966ffSPatrick Sanan           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
7241e5d2365SBarry Smith @*/
725d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM da, Vec vec, void *array)
726d71ae5a4SJacob Faibussowitsch {
7271e5d2365SBarry Smith   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
7281e5d2365SBarry Smith 
7291e5d2365SBarry Smith   PetscFunctionBegin;
7309566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
7319566063dSJacob Faibussowitsch   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
7329566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
7331e5d2365SBarry Smith 
7341e5d2365SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
7359566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(vec, &N));
7361e5d2365SBarry Smith   if (N == xm * ym * zm * dof) {
7371e5d2365SBarry Smith     gxm = xm;
7381e5d2365SBarry Smith     gym = ym;
7391e5d2365SBarry Smith     gzm = zm;
7401e5d2365SBarry Smith     gxs = xs;
7411e5d2365SBarry Smith     gys = ys;
7421e5d2365SBarry Smith     gzs = zs;
7431e5d2365SBarry Smith   }
7441e5d2365SBarry Smith 
7451e5d2365SBarry Smith   if (dim == 1) {
7469566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
7471e5d2365SBarry Smith   } else if (dim == 2) {
7489566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
7491e5d2365SBarry Smith   } else if (dim == 3) {
7509566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
75163a3b9bcSJacob Faibussowitsch   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
7523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7531e5d2365SBarry Smith }
754