1af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/
247c6ae99SBarry Smith
347c6ae99SBarry Smith /*@C
4aa219208SBarry Smith DMDAVecGetArray - Returns a multiple dimension array that shares data with
512b4a537SBarry Smith the underlying vector and is indexed using the global or local dimensions of a `DMDA`.
647c6ae99SBarry Smith
720f4b53cSBarry Smith Logically Collective
847c6ae99SBarry Smith
9d8d19677SJose E. Roman Input Parameters:
1072fd0fbdSBarry Smith + da - the `DMDA`
110af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
1247c6ae99SBarry Smith
1347c6ae99SBarry Smith Output Parameter:
1447c6ae99SBarry Smith . array - the array
1547c6ae99SBarry Smith
16dce8aebaSBarry Smith Level: intermediate
17dce8aebaSBarry Smith
1847c6ae99SBarry Smith Notes:
19dce8aebaSBarry Smith Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
2047c6ae99SBarry Smith
2147c6ae99SBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
2247c6ae99SBarry Smith
230af9b551SBarry Smith If `vec` is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
240af9b551SBarry Smith a global vector then the ghost points are not accessible. Of course, with a local vector you will have had to do the
25dce8aebaSBarry Smith appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
2647c6ae99SBarry Smith
270af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
280af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
290af9b551SBarry Smith
3095452b02SPatrick Sanan Fortran Notes:
31*ce78bad3SBarry Smith Use `DMDAVecGetArray()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
32dce8aebaSBarry 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
330f99b6f4SBarry Smith dimension of the `DMDA`.
3447c6ae99SBarry Smith
350af9b551SBarry 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
360af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
370af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
3847c6ae99SBarry Smith
3912b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
40db781477SPatrick Sanan `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
41db781477SPatrick Sanan `DMStagVecGetArray()`
4247c6ae99SBarry Smith @*/
DMDAVecGetArray(DM da,Vec vec,void * array)43d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArray(DM da, Vec vec, void *array)
44d71ae5a4SJacob Faibussowitsch {
4547c6ae99SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
4647c6ae99SBarry Smith
4747c6ae99SBarry Smith PetscFunctionBegin;
48a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
496db82c94SMatthew G Knepley PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
504f572ea9SToby Isaac PetscAssertPointer(array, 3);
519566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
529566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
539566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5447c6ae99SBarry Smith
5547c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */
569566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N));
5747c6ae99SBarry Smith if (N == xm * ym * zm * dof) {
5847c6ae99SBarry Smith gxm = xm;
5947c6ae99SBarry Smith gym = ym;
6047c6ae99SBarry Smith gzm = zm;
6147c6ae99SBarry Smith gxs = xs;
6247c6ae99SBarry Smith gys = ys;
6347c6ae99SBarry Smith gzs = zs;
6463a3b9bcSJacob 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);
6547c6ae99SBarry Smith
6647c6ae99SBarry Smith if (dim == 1) {
679566063dSJacob Faibussowitsch PetscCall(VecGetArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
6847c6ae99SBarry Smith } else if (dim == 2) {
699566063dSJacob Faibussowitsch PetscCall(VecGetArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
7047c6ae99SBarry Smith } else if (dim == 3) {
719566063dSJacob Faibussowitsch PetscCall(VecGetArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
7263a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7447c6ae99SBarry Smith }
7547c6ae99SBarry Smith
7647c6ae99SBarry Smith /*@
77dce8aebaSBarry Smith DMDAVecRestoreArray - Restores a multiple dimension array obtained with `DMDAVecGetArray()`
7847c6ae99SBarry Smith
7920f4b53cSBarry Smith Logically Collective
8047c6ae99SBarry Smith
81d8d19677SJose E. Roman Input Parameters:
8272fd0fbdSBarry Smith + da - the `DMDA`
830af9b551SBarry Smith . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
8412b4a537SBarry Smith - array - the `array` pointer
8547c6ae99SBarry Smith
8647c6ae99SBarry Smith Level: intermediate
8747c6ae99SBarry Smith
88*ce78bad3SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`,
89db781477SPatrick Sanan `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
90f8d70eaaSPierre Jolivet `DMStagVecRestoreArray()`
9147c6ae99SBarry Smith @*/
DMDAVecRestoreArray(DM da,Vec vec,void * array)92d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArray(DM da, Vec vec, void *array)
93d71ae5a4SJacob Faibussowitsch {
9447c6ae99SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
9547c6ae99SBarry Smith
9647c6ae99SBarry Smith PetscFunctionBegin;
97a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
986db82c94SMatthew G Knepley PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
994f572ea9SToby Isaac PetscAssertPointer(array, 3);
1009566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
1019566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
1029566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
10347c6ae99SBarry Smith
10447c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */
1059566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N));
10647c6ae99SBarry Smith if (N == xm * ym * zm * dof) {
10747c6ae99SBarry Smith gxm = xm;
10847c6ae99SBarry Smith gym = ym;
10947c6ae99SBarry Smith gzm = zm;
11047c6ae99SBarry Smith gxs = xs;
11147c6ae99SBarry Smith gys = ys;
11247c6ae99SBarry Smith gzs = zs;
11363a3b9bcSJacob 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);
11447c6ae99SBarry Smith
11547c6ae99SBarry Smith if (dim == 1) {
1169566063dSJacob Faibussowitsch PetscCall(VecRestoreArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
11747c6ae99SBarry Smith } else if (dim == 2) {
1189566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
11947c6ae99SBarry Smith } else if (dim == 3) {
1209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
12163a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
1223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
12347c6ae99SBarry Smith }
12447c6ae99SBarry Smith
12547c6ae99SBarry Smith /*@C
126fdc842d1SBarry Smith DMDAVecGetArrayWrite - Returns a multiple dimension array that shares data with
12712b4a537SBarry Smith the underlying vector and is indexed using the global or local dimensions of a `DMDA`.
128fdc842d1SBarry Smith
12920f4b53cSBarry Smith Logically Collective
130fdc842d1SBarry Smith
131d8d19677SJose E. Roman Input Parameters:
13272fd0fbdSBarry Smith + da - the `DMDA`
1330af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
134fdc842d1SBarry Smith
135fdc842d1SBarry Smith Output Parameter:
136fdc842d1SBarry Smith . array - the array
137fdc842d1SBarry Smith
138dce8aebaSBarry Smith Level: intermediate
139dce8aebaSBarry Smith
140fdc842d1SBarry Smith Notes:
141dce8aebaSBarry Smith Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
142fdc842d1SBarry Smith
143fdc842d1SBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
144fdc842d1SBarry Smith
1450af9b551SBarry Smith if `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
146fdc842d1SBarry Smith a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
147dce8aebaSBarry Smith appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
148fdc842d1SBarry Smith
1490af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
1500af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
1510af9b551SBarry Smith
152fdc842d1SBarry Smith Fortran Notes:
153*ce78bad3SBarry Smith Use `DMDAVecGetArrayWrite()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
154dce8aebaSBarry 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
1550f99b6f4SBarry Smith dimension of the `DMDA`.
156fdc842d1SBarry Smith
1570af9b551SBarry 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
1580af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
1590af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
160fdc842d1SBarry Smith
16112b4a537SBarry Smith Developer Note:
162dce8aebaSBarry Smith This has code duplication with `DMDAVecGetArray()` and `DMDAVecGetArrayRead()`
163fdc842d1SBarry Smith
16412b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecRestoreArrayDOF()`
165db781477SPatrick Sanan `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
166fdc842d1SBarry Smith @*/
DMDAVecGetArrayWrite(DM da,Vec vec,void * array)167d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayWrite(DM da, Vec vec, void *array)
168d71ae5a4SJacob Faibussowitsch {
169fdc842d1SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
170fdc842d1SBarry Smith
171fdc842d1SBarry Smith PetscFunctionBegin;
172fdc842d1SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
173fdc842d1SBarry Smith PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
1744f572ea9SToby Isaac PetscAssertPointer(array, 3);
1751bb6d2a8SBarry Smith if (da->localSection) {
1769566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(vec, (PetscScalar **)array));
1773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
178fdc842d1SBarry Smith }
1799566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
1809566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
1819566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
182fdc842d1SBarry Smith
183fdc842d1SBarry Smith /* Handle case where user passes in global vector as opposed to local */
1849566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N));
185fdc842d1SBarry Smith if (N == xm * ym * zm * dof) {
186fdc842d1SBarry Smith gxm = xm;
187fdc842d1SBarry Smith gym = ym;
188fdc842d1SBarry Smith gzm = zm;
189fdc842d1SBarry Smith gxs = xs;
190fdc842d1SBarry Smith gys = ys;
191fdc842d1SBarry Smith gzs = zs;
19263a3b9bcSJacob 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);
193fdc842d1SBarry Smith
194fdc842d1SBarry Smith if (dim == 1) {
1959566063dSJacob Faibussowitsch PetscCall(VecGetArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
196fdc842d1SBarry Smith } else if (dim == 2) {
1979566063dSJacob Faibussowitsch PetscCall(VecGetArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
198fdc842d1SBarry Smith } else if (dim == 3) {
1999566063dSJacob Faibussowitsch PetscCall(VecGetArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
20063a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
2013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
202fdc842d1SBarry Smith }
203fdc842d1SBarry Smith
204fdc842d1SBarry Smith /*@
205dce8aebaSBarry Smith DMDAVecRestoreArrayWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayWrite()`
206fdc842d1SBarry Smith
20720f4b53cSBarry Smith Logically Collective
208fdc842d1SBarry Smith
209d8d19677SJose E. Roman Input Parameters:
21072fd0fbdSBarry Smith + da - the `DMDA`
2110af9b551SBarry Smith . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
21212b4a537SBarry Smith - array - the `array` pointer
213fdc842d1SBarry Smith
214fdc842d1SBarry Smith Level: intermediate
215fdc842d1SBarry Smith
21612b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayWrite()`,
217db781477SPatrick Sanan `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
218fdc842d1SBarry Smith @*/
DMDAVecRestoreArrayWrite(DM da,Vec vec,void * array)219d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayWrite(DM da, Vec vec, void *array)
220d71ae5a4SJacob Faibussowitsch {
221fdc842d1SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
222fdc842d1SBarry Smith
223fdc842d1SBarry Smith PetscFunctionBegin;
224fdc842d1SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
225fdc842d1SBarry Smith PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
2264f572ea9SToby Isaac PetscAssertPointer(array, 3);
2271bb6d2a8SBarry Smith if (da->localSection) {
2289566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vec, (PetscScalar **)array));
2293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
230fdc842d1SBarry Smith }
2319566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
2329566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
2339566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
234fdc842d1SBarry Smith
235fdc842d1SBarry Smith /* Handle case where user passes in global vector as opposed to local */
2369566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N));
237fdc842d1SBarry Smith if (N == xm * ym * zm * dof) {
238fdc842d1SBarry Smith gxm = xm;
239fdc842d1SBarry Smith gym = ym;
240fdc842d1SBarry Smith gzm = zm;
241fdc842d1SBarry Smith gxs = xs;
242fdc842d1SBarry Smith gys = ys;
243fdc842d1SBarry Smith gzs = zs;
24463a3b9bcSJacob 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);
245fdc842d1SBarry Smith
246fdc842d1SBarry Smith if (dim == 1) {
2479566063dSJacob Faibussowitsch PetscCall(VecRestoreArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
248fdc842d1SBarry Smith } else if (dim == 2) {
2499566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
250fdc842d1SBarry Smith } else if (dim == 3) {
2519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
25263a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
2533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
254fdc842d1SBarry Smith }
255fdc842d1SBarry Smith
256fdc842d1SBarry Smith /*@C
257aa219208SBarry Smith DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
25812b4a537SBarry Smith the underlying vector and is indexed using the global or local dimensions of a `DMDA`
25947c6ae99SBarry Smith
26020f4b53cSBarry Smith Logically Collective
26147c6ae99SBarry Smith
262d8d19677SJose E. Roman Input Parameters:
26372fd0fbdSBarry Smith + da - the `DMDA`
2640af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
26547c6ae99SBarry Smith
26647c6ae99SBarry Smith Output Parameter:
26712b4a537SBarry Smith . array - the `array` pointer
26847c6ae99SBarry Smith
269dce8aebaSBarry Smith Level: intermediate
270dce8aebaSBarry Smith
27147c6ae99SBarry Smith Notes:
272dce8aebaSBarry Smith Call `DMDAVecRestoreArrayDOF()` once you have finished accessing the vector entries.
27347c6ae99SBarry Smith
2740af9b551SBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]
2750af9b551SBarry Smith
2760af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1][0:ndof-1]` where the values are obtained from
2770af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
27847c6ae99SBarry Smith
2790f99b6f4SBarry Smith Fortran Notes:
280*ce78bad3SBarry Smith Use `DMDAVecGetArray()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
2810f99b6f4SBarry 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
2820f99b6f4SBarry Smith dimension of the `DMDA`.
2830f99b6f4SBarry Smith
2840af9b551SBarry Smith The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when ndof is 1) otherwise
2850af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
2860af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
2871b82215eSBarry Smith
28812b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecRestoreArrayDOF()`,
2894ab966ffSPatrick Sanan `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, `DMDAVecGetArrayDOFRead()`
29047c6ae99SBarry Smith @*/
DMDAVecGetArrayDOF(DM da,Vec vec,void * array)291d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOF(DM da, Vec vec, void *array)
292d71ae5a4SJacob Faibussowitsch {
29347c6ae99SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
29447c6ae99SBarry Smith
29547c6ae99SBarry Smith PetscFunctionBegin;
2969566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
2979566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
2989566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
29947c6ae99SBarry Smith
30047c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */
3019566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N));
30247c6ae99SBarry Smith if (N == xm * ym * zm * dof) {
30347c6ae99SBarry Smith gxm = xm;
30447c6ae99SBarry Smith gym = ym;
30547c6ae99SBarry Smith gzm = zm;
30647c6ae99SBarry Smith gxs = xs;
30747c6ae99SBarry Smith gys = ys;
30847c6ae99SBarry Smith gzs = zs;
30963a3b9bcSJacob 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);
31047c6ae99SBarry Smith
31147c6ae99SBarry Smith if (dim == 1) {
3129566063dSJacob Faibussowitsch PetscCall(VecGetArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
31347c6ae99SBarry Smith } else if (dim == 2) {
3149566063dSJacob Faibussowitsch PetscCall(VecGetArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
31547c6ae99SBarry Smith } else if (dim == 3) {
3169566063dSJacob Faibussowitsch PetscCall(VecGetArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
31763a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
3183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
31947c6ae99SBarry Smith }
32047c6ae99SBarry Smith
32147c6ae99SBarry Smith /*@
322dce8aebaSBarry Smith DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOF()`
32347c6ae99SBarry Smith
32420f4b53cSBarry Smith Logically Collective
32547c6ae99SBarry Smith
326d8d19677SJose E. Roman Input Parameters:
32772fd0fbdSBarry Smith + da - the `DMDA`
3280af9b551SBarry Smith . vec - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
32912b4a537SBarry Smith - array - the `array` point
33047c6ae99SBarry Smith
33147c6ae99SBarry Smith Level: intermediate
33247c6ae99SBarry Smith
33312b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
3344ab966ffSPatrick Sanan `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
33547c6ae99SBarry Smith @*/
DMDAVecRestoreArrayDOF(DM da,Vec vec,void * array)336d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOF(DM da, Vec vec, void *array)
337d71ae5a4SJacob Faibussowitsch {
33847c6ae99SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
33947c6ae99SBarry Smith
34047c6ae99SBarry Smith PetscFunctionBegin;
3419566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
3429566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
3439566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
34447c6ae99SBarry Smith
34547c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */
3469566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N));
34747c6ae99SBarry Smith if (N == xm * ym * zm * dof) {
34847c6ae99SBarry Smith gxm = xm;
34947c6ae99SBarry Smith gym = ym;
35047c6ae99SBarry Smith gzm = zm;
35147c6ae99SBarry Smith gxs = xs;
35247c6ae99SBarry Smith gys = ys;
35347c6ae99SBarry Smith gzs = zs;
35447c6ae99SBarry Smith }
35547c6ae99SBarry Smith
35647c6ae99SBarry Smith if (dim == 1) {
3579566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
35847c6ae99SBarry Smith } else if (dim == 2) {
3599566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
36047c6ae99SBarry Smith } else if (dim == 3) {
3619566063dSJacob Faibussowitsch PetscCall(VecRestoreArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
36263a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
3633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
36447c6ae99SBarry Smith }
36547c6ae99SBarry Smith
3665edff71fSBarry Smith /*@C
3675edff71fSBarry Smith DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
36812b4a537SBarry Smith the underlying vector and is indexed using the global or local dimensions of a `DMDA`.
3695edff71fSBarry Smith
37020f4b53cSBarry Smith Not Collective
3715edff71fSBarry Smith
372d8d19677SJose E. Roman Input Parameters:
37372fd0fbdSBarry Smith + da - the `DMDA`
3740af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
3755edff71fSBarry Smith
3765edff71fSBarry Smith Output Parameter:
3775edff71fSBarry Smith . array - the array
3785edff71fSBarry Smith
379dce8aebaSBarry Smith Level: intermediate
380dce8aebaSBarry Smith
3815edff71fSBarry Smith Notes:
382dce8aebaSBarry Smith Call `DMDAVecRestoreArrayRead()` once you have finished accessing the vector entries.
3835edff71fSBarry Smith
3845edff71fSBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
3855edff71fSBarry Smith
3860af9b551SBarry Smith If `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
3875edff71fSBarry Smith a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
388dce8aebaSBarry Smith appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
3895edff71fSBarry Smith
3900af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
3910af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
3920af9b551SBarry Smith
39395452b02SPatrick Sanan Fortran Notes:
394*ce78bad3SBarry Smith Use `DMDAVecGetArrayRead()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
395dce8aebaSBarry 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
3960f99b6f4SBarry Smith dimension of the `DMDA`.
3970f99b6f4SBarry Smith
3980af9b551SBarry 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
3990af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
4000af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
4015edff71fSBarry Smith
402*ce78bad3SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`,
40342747ad1SJacob Faibussowitsch `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayRead()`,
40442747ad1SJacob Faibussowitsch `DMDAVecRestoreArrayDOF()`, `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`,
40542747ad1SJacob Faibussowitsch `DMDAVecRestoreArray()`, `DMStagVecGetArrayRead()`
4065edff71fSBarry Smith @*/
DMDAVecGetArrayRead(DM da,Vec vec,void * array)407d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayRead(DM da, Vec vec, void *array)
408d71ae5a4SJacob Faibussowitsch {
4095edff71fSBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
4105edff71fSBarry Smith
4115edff71fSBarry Smith PetscFunctionBegin;
412a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
4135edff71fSBarry Smith PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
4144f572ea9SToby Isaac PetscAssertPointer(array, 3);
4159566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
4169566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
4179566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
4185edff71fSBarry Smith
4195edff71fSBarry Smith /* Handle case where user passes in global vector as opposed to local */
4209566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N));
4215edff71fSBarry Smith if (N == xm * ym * zm * dof) {
4225edff71fSBarry Smith gxm = xm;
4235edff71fSBarry Smith gym = ym;
4245edff71fSBarry Smith gzm = zm;
4255edff71fSBarry Smith gxs = xs;
4265edff71fSBarry Smith gys = ys;
4275edff71fSBarry Smith gzs = zs;
42863a3b9bcSJacob 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);
4295edff71fSBarry Smith
4305edff71fSBarry Smith if (dim == 1) {
4319566063dSJacob Faibussowitsch PetscCall(VecGetArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
4325edff71fSBarry Smith } else if (dim == 2) {
4339566063dSJacob Faibussowitsch PetscCall(VecGetArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
4345edff71fSBarry Smith } else if (dim == 3) {
4359566063dSJacob Faibussowitsch PetscCall(VecGetArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
43663a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
4373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4385edff71fSBarry Smith }
4395edff71fSBarry Smith
4405edff71fSBarry Smith /*@
441dce8aebaSBarry Smith DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayRead()`
4425edff71fSBarry Smith
44320f4b53cSBarry Smith Not Collective
4445edff71fSBarry Smith
445d8d19677SJose E. Roman Input Parameters:
44672fd0fbdSBarry Smith + da - the `DMDA`
4470af9b551SBarry Smith . vec - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
44812b4a537SBarry Smith - array - the `array` pointer
4495edff71fSBarry Smith
4505edff71fSBarry Smith Level: intermediate
4515edff71fSBarry Smith
452*ce78bad3SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayRead()`,
453db781477SPatrick Sanan `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`,
454db781477SPatrick Sanan `DMStagVecRestoreArrayRead()`
4555edff71fSBarry Smith @*/
DMDAVecRestoreArrayRead(DM da,Vec vec,void * array)456d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayRead(DM da, Vec vec, void *array)
457d71ae5a4SJacob Faibussowitsch {
4585edff71fSBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
4595edff71fSBarry Smith
4605edff71fSBarry Smith PetscFunctionBegin;
461a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA);
4625edff71fSBarry Smith PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
4634f572ea9SToby Isaac PetscAssertPointer(array, 3);
4649566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
4659566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
4669566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
4675edff71fSBarry Smith
4685edff71fSBarry Smith /* Handle case where user passes in global vector as opposed to local */
4699566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N));
4705edff71fSBarry Smith if (N == xm * ym * zm * dof) {
4715edff71fSBarry Smith gxm = xm;
4725edff71fSBarry Smith gym = ym;
4735edff71fSBarry Smith gzm = zm;
4745edff71fSBarry Smith gxs = xs;
4755edff71fSBarry Smith gys = ys;
4765edff71fSBarry Smith gzs = zs;
47763a3b9bcSJacob 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);
4785edff71fSBarry Smith
4795edff71fSBarry Smith if (dim == 1) {
4809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
4815edff71fSBarry Smith } else if (dim == 2) {
4829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
4835edff71fSBarry Smith } else if (dim == 3) {
4849566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
48563a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
4863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
4875edff71fSBarry Smith }
4885edff71fSBarry Smith
4895edff71fSBarry Smith /*@C
4905edff71fSBarry Smith DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
49112b4a537SBarry Smith the underlying vector and is indexed using the global or local dimensions of a `DMDA`
4925edff71fSBarry Smith
4935edff71fSBarry Smith Not Collective
4945edff71fSBarry Smith
495d8d19677SJose E. Roman Input Parameters:
49672fd0fbdSBarry Smith + da - the `DMDA`
4970af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
4985edff71fSBarry Smith
4995edff71fSBarry Smith Output Parameter:
5005edff71fSBarry Smith . array - the array
5015edff71fSBarry Smith
502dce8aebaSBarry Smith Level: intermediate
503dce8aebaSBarry Smith
5045edff71fSBarry Smith Notes:
505dce8aebaSBarry Smith Call `DMDAVecRestoreArrayDOFRead()` once you have finished accessing the vector entries.
5065edff71fSBarry Smith
5075edff71fSBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
5085edff71fSBarry Smith
5090af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
5100af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
5110af9b551SBarry Smith
51260225df5SJacob Faibussowitsch Fortran Notes:
513*ce78bad3SBarry Smith Use `DMDAVecGetArrayRead()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
5140f99b6f4SBarry 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
5150f99b6f4SBarry Smith dimension of the `DMDA`.
5160f99b6f4SBarry Smith
5170af9b551SBarry 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
5180af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
5190af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
5205edff71fSBarry Smith
52112b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
5224ab966ffSPatrick Sanan `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
5235edff71fSBarry Smith @*/
DMDAVecGetArrayDOFRead(DM da,Vec vec,void * array)524d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOFRead(DM da, Vec vec, void *array)
525d71ae5a4SJacob Faibussowitsch {
5265edff71fSBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
5275edff71fSBarry Smith
5285edff71fSBarry Smith PetscFunctionBegin;
5299566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
5309566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
5319566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5325edff71fSBarry Smith
5335edff71fSBarry Smith /* Handle case where user passes in global vector as opposed to local */
5349566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N));
5355edff71fSBarry Smith if (N == xm * ym * zm * dof) {
5365edff71fSBarry Smith gxm = xm;
5375edff71fSBarry Smith gym = ym;
5385edff71fSBarry Smith gzm = zm;
5395edff71fSBarry Smith gxs = xs;
5405edff71fSBarry Smith gys = ys;
5415edff71fSBarry Smith gzs = zs;
54263a3b9bcSJacob 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);
5435edff71fSBarry Smith
5445edff71fSBarry Smith if (dim == 1) {
5459566063dSJacob Faibussowitsch PetscCall(VecGetArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
5465edff71fSBarry Smith } else if (dim == 2) {
5479566063dSJacob Faibussowitsch PetscCall(VecGetArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
5485edff71fSBarry Smith } else if (dim == 3) {
5499566063dSJacob Faibussowitsch PetscCall(VecGetArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
55063a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
5513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5525edff71fSBarry Smith }
5535edff71fSBarry Smith
5545edff71fSBarry Smith /*@
555dce8aebaSBarry Smith DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFRead()`
5565edff71fSBarry Smith
5575edff71fSBarry Smith Not Collective
5585edff71fSBarry Smith
559d8d19677SJose E. Roman Input Parameters:
56072fd0fbdSBarry Smith + da - the `DMDA`
5610af9b551SBarry Smith . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
56212b4a537SBarry Smith - array - the `array` pointer
5635edff71fSBarry Smith
5645edff71fSBarry Smith Level: intermediate
5655edff71fSBarry Smith
56612b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
5674ab966ffSPatrick Sanan `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
5685edff71fSBarry Smith @*/
DMDAVecRestoreArrayDOFRead(DM da,Vec vec,void * array)569d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOFRead(DM da, Vec vec, void *array)
570d71ae5a4SJacob Faibussowitsch {
5715edff71fSBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
5725edff71fSBarry Smith
5735edff71fSBarry Smith PetscFunctionBegin;
5749566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
5759566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
5769566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
5775edff71fSBarry Smith
5785edff71fSBarry Smith /* Handle case where user passes in global vector as opposed to local */
5799566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N));
5805edff71fSBarry Smith if (N == xm * ym * zm * dof) {
5815edff71fSBarry Smith gxm = xm;
5825edff71fSBarry Smith gym = ym;
5835edff71fSBarry Smith gzm = zm;
5845edff71fSBarry Smith gxs = xs;
5855edff71fSBarry Smith gys = ys;
5865edff71fSBarry Smith gzs = zs;
5875edff71fSBarry Smith }
5885edff71fSBarry Smith
5895edff71fSBarry Smith if (dim == 1) {
5909566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
5915edff71fSBarry Smith } else if (dim == 2) {
5929566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
5935edff71fSBarry Smith } else if (dim == 3) {
5949566063dSJacob Faibussowitsch PetscCall(VecRestoreArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
59563a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
5963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
5975edff71fSBarry Smith }
5985edff71fSBarry Smith
5991e5d2365SBarry Smith /*@C
6001e5d2365SBarry Smith DMDAVecGetArrayDOFWrite - Returns a multiple dimension array that shares data with
60112b4a537SBarry Smith the underlying vector and is indexed using the global or local dimensions of a `DMDA`
6021e5d2365SBarry Smith
6031e5d2365SBarry Smith Not Collective
6041e5d2365SBarry Smith
605d8d19677SJose E. Roman Input Parameters:
60672fd0fbdSBarry Smith + da - the `DMDA`
6070af9b551SBarry Smith - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
6081e5d2365SBarry Smith
6091e5d2365SBarry Smith Output Parameter:
6101e5d2365SBarry Smith . array - the array
6111e5d2365SBarry Smith
6120f99b6f4SBarry Smith Level: intermediate
6130f99b6f4SBarry Smith
6141e5d2365SBarry Smith Notes:
615dce8aebaSBarry Smith Call `DMDAVecRestoreArrayDOFWrite()` once you have finished accessing the vector entries.
6161e5d2365SBarry Smith
6171e5d2365SBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
6181e5d2365SBarry Smith
6190af9b551SBarry Smith The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1][0:dof-1]` where the values are obtained from
6200af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
6210af9b551SBarry Smith
62260225df5SJacob Faibussowitsch Fortran Notes:
623*ce78bad3SBarry Smith Use `DMDAVecGetArrayWrite()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
6240f99b6f4SBarry 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
6250f99b6f4SBarry Smith dimension of the `DMDA`.
6261e5d2365SBarry Smith
6270af9b551SBarry 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
6280af9b551SBarry Smith `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
6290af9b551SBarry Smith `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
6301e5d2365SBarry Smith
631*ce78bad3SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
63260225df5SJacob Faibussowitsch `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
6331e5d2365SBarry Smith @*/
DMDAVecGetArrayDOFWrite(DM da,Vec vec,void * array)634d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGetArrayDOFWrite(DM da, Vec vec, void *array)
635d71ae5a4SJacob Faibussowitsch {
6361e5d2365SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
6371e5d2365SBarry Smith
6381e5d2365SBarry Smith PetscFunctionBegin;
6399566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
6409566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
6419566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
6421e5d2365SBarry Smith
6431e5d2365SBarry Smith /* Handle case where user passes in global vector as opposed to local */
6449566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N));
6451e5d2365SBarry Smith if (N == xm * ym * zm * dof) {
6461e5d2365SBarry Smith gxm = xm;
6471e5d2365SBarry Smith gym = ym;
6481e5d2365SBarry Smith gzm = zm;
6491e5d2365SBarry Smith gxs = xs;
6501e5d2365SBarry Smith gys = ys;
6511e5d2365SBarry Smith gzs = zs;
65263a3b9bcSJacob 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);
6531e5d2365SBarry Smith
6541e5d2365SBarry Smith if (dim == 1) {
6559566063dSJacob Faibussowitsch PetscCall(VecGetArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
6561e5d2365SBarry Smith } else if (dim == 2) {
6579566063dSJacob Faibussowitsch PetscCall(VecGetArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
6581e5d2365SBarry Smith } else if (dim == 3) {
6599566063dSJacob Faibussowitsch PetscCall(VecGetArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
66063a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
6613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
6621e5d2365SBarry Smith }
6631e5d2365SBarry Smith
6641e5d2365SBarry Smith /*@
665dce8aebaSBarry Smith DMDAVecRestoreArrayDOFWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFWrite()`
6661e5d2365SBarry Smith
6671e5d2365SBarry Smith Not Collective
6681e5d2365SBarry Smith
669d8d19677SJose E. Roman Input Parameters:
67072fd0fbdSBarry Smith + da - the `DMDA`
6710af9b551SBarry Smith . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
67212b4a537SBarry Smith - array - the `array` pointer
6731e5d2365SBarry Smith
6741e5d2365SBarry Smith Level: intermediate
6751e5d2365SBarry Smith
67612b4a537SBarry Smith .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
67760225df5SJacob Faibussowitsch `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
6781e5d2365SBarry Smith @*/
DMDAVecRestoreArrayDOFWrite(DM da,Vec vec,void * array)679d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM da, Vec vec, void *array)
680d71ae5a4SJacob Faibussowitsch {
6811e5d2365SBarry Smith PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
6821e5d2365SBarry Smith
6831e5d2365SBarry Smith PetscFunctionBegin;
6849566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
6859566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
6869566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
6871e5d2365SBarry Smith
6881e5d2365SBarry Smith /* Handle case where user passes in global vector as opposed to local */
6899566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec, &N));
6901e5d2365SBarry Smith if (N == xm * ym * zm * dof) {
6911e5d2365SBarry Smith gxm = xm;
6921e5d2365SBarry Smith gym = ym;
6931e5d2365SBarry Smith gzm = zm;
6941e5d2365SBarry Smith gxs = xs;
6951e5d2365SBarry Smith gys = ys;
6961e5d2365SBarry Smith gzs = zs;
6971e5d2365SBarry Smith }
6981e5d2365SBarry Smith
6991e5d2365SBarry Smith if (dim == 1) {
7009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
7011e5d2365SBarry Smith } else if (dim == 2) {
7029566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
7031e5d2365SBarry Smith } else if (dim == 3) {
7049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
70563a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
7063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
7071e5d2365SBarry Smith }
708