xref: /petsc/src/dm/impls/da/dagetarray.c (revision 1e5d23653b4749164fe87aa5a60b30836473cd8d)
147c6ae99SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
347c6ae99SBarry Smith 
447c6ae99SBarry Smith /*@C
5aa219208SBarry Smith    DMDAVecGetArray - Returns a multiple dimension array that shares data with
647c6ae99SBarry Smith       the underlying vector and is indexed using the global dimensions.
747c6ae99SBarry Smith 
8d083f849SBarry Smith    Logically collective on da
947c6ae99SBarry Smith 
1047c6ae99SBarry Smith    Input Parameter:
1147c6ae99SBarry Smith +  da - the distributed array
127e57d48aSBarry Smith -  vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
1347c6ae99SBarry Smith 
1447c6ae99SBarry Smith    Output Parameter:
1547c6ae99SBarry Smith .  array - the array
1647c6ae99SBarry Smith 
1747c6ae99SBarry Smith    Notes:
18aa219208SBarry Smith     Call DMDAVecRestoreArray() once you have finished accessing the vector entries.
1947c6ae99SBarry Smith 
2047c6ae99SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
2147c6ae99SBarry Smith 
227e57d48aSBarry Smith     If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
237e57d48aSBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
2447c6ae99SBarry Smith 
257e57d48aSBarry Smith     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
2647c6ae99SBarry Smith 
2795452b02SPatrick Sanan   Fortran Notes:
2895452b02SPatrick Sanan     From Fortran use DMDAVecGetArrayF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
29aa219208SBarry 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
30aa219208SBarry Smith        dimension of the DMDA. The order of the indices is array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) (when dof is 1) otherwise
31bef27881SBarry Smith        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
32af0996ceSBarry Smith        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
3347c6ae99SBarry Smith 
34596f5af7SJed Brown   Due to bugs in the compiler DMDAVecGetArrayF90() does not work with gfortran versions before 4.5
3547c6ae99SBarry Smith 
3647c6ae99SBarry Smith   Level: intermediate
3747c6ae99SBarry Smith 
38aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF()
393bd220d7SPatrick Sanan           DMDAVecGetArrayDOF(), DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(),
403bd220d7SPatrick Sanan           DMStagVecGetArray()
4147c6ae99SBarry Smith @*/
427087cfbeSBarry Smith PetscErrorCode  DMDAVecGetArray(DM da,Vec vec,void *array)
4347c6ae99SBarry Smith {
4447c6ae99SBarry Smith   PetscErrorCode ierr;
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);
506db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
51aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
52aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
53ea78f98cSLisandro Dalcin   ierr = DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
5447c6ae99SBarry Smith 
5547c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
5647c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
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;
6430729d88SBarry Smith   } else if (N != gxm*gym*gzm*dof) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMDA local sizes %D %D\n",N,xm*ym*zm*dof,gxm*gym*gzm*dof);
6547c6ae99SBarry Smith 
6647c6ae99SBarry Smith   if (dim == 1) {
6747c6ae99SBarry Smith     ierr = VecGetArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
6847c6ae99SBarry Smith   } else if (dim == 2) {
6947c6ae99SBarry Smith     ierr = VecGetArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
7047c6ae99SBarry Smith   } else if (dim == 3) {
7147c6ae99SBarry Smith     ierr = VecGetArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
7230729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
7347c6ae99SBarry Smith   PetscFunctionReturn(0);
7447c6ae99SBarry Smith }
7547c6ae99SBarry Smith 
7647c6ae99SBarry Smith /*@
77aa219208SBarry Smith    DMDAVecRestoreArray - Restores a multiple dimension array obtained with DMDAVecGetArray()
7847c6ae99SBarry Smith 
79d083f849SBarry Smith    Logically collective on da
8047c6ae99SBarry Smith 
8147c6ae99SBarry Smith    Input Parameter:
8247c6ae99SBarry Smith +  da - the distributed array
8347c6ae99SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
84564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
85e5c84f05SJed Brown -  array - the array, non-NULL pointer is zeroed
8647c6ae99SBarry Smith 
8747c6ae99SBarry Smith   Level: intermediate
8847c6ae99SBarry Smith 
8995452b02SPatrick Sanan   Fortran Notes:
90fdc842d1SBarry Smith     From Fortran use DMDAVecRestoreArayF90()
9147c6ae99SBarry Smith 
92fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(),
933bd220d7SPatrick Sanan           DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(),
943bd220d7SPatrick Sanan           DMDStagVecRestoreArray()
9547c6ae99SBarry Smith @*/
967087cfbeSBarry Smith PetscErrorCode  DMDAVecRestoreArray(DM da,Vec vec,void *array)
9747c6ae99SBarry Smith {
9847c6ae99SBarry Smith   PetscErrorCode ierr;
9947c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
10047c6ae99SBarry Smith 
10147c6ae99SBarry Smith   PetscFunctionBegin;
102a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
1036db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
1046db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
105aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
106aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
107ea78f98cSLisandro Dalcin   ierr = DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
10847c6ae99SBarry Smith 
10947c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
11047c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
11147c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
11247c6ae99SBarry Smith     gxm = xm;
11347c6ae99SBarry Smith     gym = ym;
11447c6ae99SBarry Smith     gzm = zm;
11547c6ae99SBarry Smith     gxs = xs;
11647c6ae99SBarry Smith     gys = ys;
11747c6ae99SBarry Smith     gzs = zs;
11830729d88SBarry Smith   } else if (N != gxm*gym*gzm*dof) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMDA local sizes %D %D\n",N,xm*ym*zm*dof,gxm*gym*gzm*dof);
11947c6ae99SBarry Smith 
12047c6ae99SBarry Smith   if (dim == 1) {
12147c6ae99SBarry Smith     ierr = VecRestoreArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
12247c6ae99SBarry Smith   } else if (dim == 2) {
12347c6ae99SBarry Smith     ierr = VecRestoreArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
12447c6ae99SBarry Smith   } else if (dim == 3) {
12547c6ae99SBarry Smith     ierr = VecRestoreArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
12630729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
12747c6ae99SBarry Smith   PetscFunctionReturn(0);
12847c6ae99SBarry Smith }
12947c6ae99SBarry Smith 
13047c6ae99SBarry Smith /*@C
131fdc842d1SBarry Smith    DMDAVecGetArrayWrite - Returns a multiple dimension array that shares data with
132fdc842d1SBarry Smith       the underlying vector and is indexed using the global dimensions.
133fdc842d1SBarry Smith 
134fdc842d1SBarry Smith    Logically collective on Vec
135fdc842d1SBarry Smith 
136fdc842d1SBarry Smith    Input Parameter:
137fdc842d1SBarry Smith +  da - the distributed array
138fdc842d1SBarry Smith -  vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
139fdc842d1SBarry Smith 
140fdc842d1SBarry Smith    Output Parameter:
141fdc842d1SBarry Smith .  array - the array
142fdc842d1SBarry Smith 
143fdc842d1SBarry Smith    Notes:
144fdc842d1SBarry Smith     Call DMDAVecRestoreArray() once you have finished accessing the vector entries.
145fdc842d1SBarry Smith 
146fdc842d1SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
147fdc842d1SBarry Smith 
148fdc842d1SBarry Smith     If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
149fdc842d1SBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
150fdc842d1SBarry Smith 
151fdc842d1SBarry Smith     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
152fdc842d1SBarry Smith 
153fdc842d1SBarry Smith   Fortran Notes:
154fdc842d1SBarry Smith     From Fortran use DMDAVecGetArrayF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
155fdc842d1SBarry 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
156fdc842d1SBarry Smith        dimension of the DMDA. The order of the indices is array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) (when dof is 1) otherwise
157fdc842d1SBarry Smith        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
158fdc842d1SBarry Smith        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
159fdc842d1SBarry Smith 
160fdc842d1SBarry Smith   Due to bugs in the compiler DMDAVecGetArrayF90() does not work with gfortran versions before 4.5
161fdc842d1SBarry Smith 
162fdc842d1SBarry Smith   Level: intermediate
163fdc842d1SBarry Smith 
164fdc842d1SBarry Smith   Developer Notes: This has code duplication with DMDAVecGetArray() and DMDAVecGetArrayRead()
165fdc842d1SBarry Smith 
166fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArrayWrite(), DMDAVecRestoreArrayDOF()
167fdc842d1SBarry Smith           DMDAVecGetArrayDOF(), DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead()
168fdc842d1SBarry Smith @*/
169fdc842d1SBarry Smith PetscErrorCode  DMDAVecGetArrayWrite(DM da,Vec vec,void *array)
170fdc842d1SBarry Smith {
171fdc842d1SBarry Smith   PetscErrorCode ierr;
172fdc842d1SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
173fdc842d1SBarry Smith 
174fdc842d1SBarry Smith   PetscFunctionBegin;
175fdc842d1SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
176fdc842d1SBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
177fdc842d1SBarry Smith   PetscValidPointer(array, 3);
1781bb6d2a8SBarry Smith   if (da->localSection) {
179fdc842d1SBarry Smith     ierr = VecGetArrayWrite(vec,(PetscScalar**)array);CHKERRQ(ierr);
180fdc842d1SBarry Smith     PetscFunctionReturn(0);
181fdc842d1SBarry Smith   }
182fdc842d1SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
183fdc842d1SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
184ea78f98cSLisandro Dalcin   ierr = DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
185fdc842d1SBarry Smith 
186fdc842d1SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
187fdc842d1SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
188fdc842d1SBarry Smith   if (N == xm*ym*zm*dof) {
189fdc842d1SBarry Smith     gxm = xm;
190fdc842d1SBarry Smith     gym = ym;
191fdc842d1SBarry Smith     gzm = zm;
192fdc842d1SBarry Smith     gxs = xs;
193fdc842d1SBarry Smith     gys = ys;
194fdc842d1SBarry Smith     gzs = zs;
195fdc842d1SBarry Smith   } else if (N != gxm*gym*gzm*dof) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMDA local sizes %D %D\n",N,xm*ym*zm*dof,gxm*gym*gzm*dof);
196fdc842d1SBarry Smith 
197fdc842d1SBarry Smith   if (dim == 1) {
198fdc842d1SBarry Smith     ierr = VecGetArray1dWrite(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
199fdc842d1SBarry Smith   } else if (dim == 2) {
200fdc842d1SBarry Smith     ierr = VecGetArray2dWrite(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
201fdc842d1SBarry Smith   } else if (dim == 3) {
202fdc842d1SBarry Smith     ierr = VecGetArray3dWrite(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
203fdc842d1SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
204fdc842d1SBarry Smith   PetscFunctionReturn(0);
205fdc842d1SBarry Smith }
206fdc842d1SBarry Smith 
207fdc842d1SBarry Smith /*@
208fdc842d1SBarry Smith    DMDAVecRestoreArrayWrite - Restores a multiple dimension array obtained with DMDAVecGetArrayWrite()
209fdc842d1SBarry Smith 
210fdc842d1SBarry Smith    Logically collective on Vec
211fdc842d1SBarry Smith 
212fdc842d1SBarry Smith    Input Parameter:
213fdc842d1SBarry Smith +  da - the distributed array
214fdc842d1SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
215fdc842d1SBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
216fdc842d1SBarry Smith -  array - the array, non-NULL pointer is zeroed
217fdc842d1SBarry Smith 
218fdc842d1SBarry Smith   Level: intermediate
219fdc842d1SBarry Smith 
220fdc842d1SBarry Smith   Fortran Notes:
221fdc842d1SBarry Smith     From Fortran use DMDAVecRestoreArayF90()
222fdc842d1SBarry Smith 
223fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArrayWrite(),
224fdc842d1SBarry Smith           DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead()
225fdc842d1SBarry Smith @*/
226fdc842d1SBarry Smith PetscErrorCode  DMDAVecRestoreArrayWrite(DM da,Vec vec,void *array)
227fdc842d1SBarry Smith {
228fdc842d1SBarry Smith   PetscErrorCode ierr;
229fdc842d1SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
230fdc842d1SBarry Smith 
231fdc842d1SBarry Smith   PetscFunctionBegin;
232fdc842d1SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
233fdc842d1SBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
234fdc842d1SBarry Smith   PetscValidPointer(array, 3);
2351bb6d2a8SBarry Smith   if (da->localSection) {
236fdc842d1SBarry Smith     ierr = VecRestoreArray(vec,(PetscScalar**)array);CHKERRQ(ierr);
237fdc842d1SBarry Smith     PetscFunctionReturn(0);
238fdc842d1SBarry Smith   }
239fdc842d1SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
240fdc842d1SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
241ea78f98cSLisandro Dalcin   ierr = DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
242fdc842d1SBarry Smith 
243fdc842d1SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
244fdc842d1SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
245fdc842d1SBarry Smith   if (N == xm*ym*zm*dof) {
246fdc842d1SBarry Smith     gxm = xm;
247fdc842d1SBarry Smith     gym = ym;
248fdc842d1SBarry Smith     gzm = zm;
249fdc842d1SBarry Smith     gxs = xs;
250fdc842d1SBarry Smith     gys = ys;
251fdc842d1SBarry Smith     gzs = zs;
252fdc842d1SBarry Smith   } else if (N != gxm*gym*gzm*dof) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMDA local sizes %D %D\n",N,xm*ym*zm*dof,gxm*gym*gzm*dof);
253fdc842d1SBarry Smith 
254fdc842d1SBarry Smith   if (dim == 1) {
255fdc842d1SBarry Smith     ierr = VecRestoreArray1dWrite(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
256fdc842d1SBarry Smith   } else if (dim == 2) {
257fdc842d1SBarry Smith     ierr = VecRestoreArray2dWrite(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
258fdc842d1SBarry Smith   } else if (dim == 3) {
259fdc842d1SBarry Smith     ierr = VecRestoreArray3dWrite(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
260fdc842d1SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
261fdc842d1SBarry Smith   PetscFunctionReturn(0);
262fdc842d1SBarry Smith }
263fdc842d1SBarry Smith 
264fdc842d1SBarry Smith /*@C
265aa219208SBarry Smith    DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
26647c6ae99SBarry Smith       the underlying vector and is indexed using the global dimensions.
26747c6ae99SBarry Smith 
2682ddbf755SMatthew G. Knepley    Logically collective
26947c6ae99SBarry Smith 
27047c6ae99SBarry Smith    Input Parameter:
27147c6ae99SBarry Smith +  da - the distributed array
27247c6ae99SBarry Smith -  vec - the vector, either a vector the same size as one obtained with
273564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
27447c6ae99SBarry Smith 
27547c6ae99SBarry Smith    Output Parameter:
27647c6ae99SBarry Smith .  array - the array
27747c6ae99SBarry Smith 
27847c6ae99SBarry Smith    Notes:
279aa219208SBarry Smith     Call DMDAVecRestoreArrayDOF() once you have finished accessing the vector entries.
28047c6ae99SBarry Smith 
28147c6ae99SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
28247c6ae99SBarry Smith 
2831b82215eSBarry Smith     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayF90() and declare your array with one higher dimension,
284c4762a1bSJed Brown     see src/dm/tutorials/ex11f90.F
2851b82215eSBarry Smith 
28647c6ae99SBarry Smith   Level: intermediate
28747c6ae99SBarry Smith 
288fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF(),
289856ac3e9SJed Brown           DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), DMStagVecGetArrayDOF(), DMDAVecGetArrayDOFRead()
29047c6ae99SBarry Smith @*/
2917087cfbeSBarry Smith PetscErrorCode  DMDAVecGetArrayDOF(DM da,Vec vec,void *array)
29247c6ae99SBarry Smith {
29347c6ae99SBarry Smith   PetscErrorCode ierr;
29447c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
29547c6ae99SBarry Smith 
29647c6ae99SBarry Smith   PetscFunctionBegin;
297aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
298aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
299ea78f98cSLisandro Dalcin   ierr = DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
30047c6ae99SBarry Smith 
30147c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
30247c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
30347c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
30447c6ae99SBarry Smith     gxm = xm;
30547c6ae99SBarry Smith     gym = ym;
30647c6ae99SBarry Smith     gzm = zm;
30747c6ae99SBarry Smith     gxs = xs;
30847c6ae99SBarry Smith     gys = ys;
30947c6ae99SBarry Smith     gzs = zs;
31030729d88SBarry Smith   } else if (N != gxm*gym*gzm*dof) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMDA local sizes %D %D\n",N,xm*ym*zm*dof,gxm*gym*gzm*dof);
31147c6ae99SBarry Smith 
31247c6ae99SBarry Smith   if (dim == 1) {
31347c6ae99SBarry Smith     ierr = VecGetArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
31447c6ae99SBarry Smith   } else if (dim == 2) {
31547c6ae99SBarry Smith     ierr = VecGetArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
31647c6ae99SBarry Smith   } else if (dim == 3) {
31747c6ae99SBarry Smith     ierr = VecGetArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
31830729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
31947c6ae99SBarry Smith   PetscFunctionReturn(0);
32047c6ae99SBarry Smith }
32147c6ae99SBarry Smith 
32247c6ae99SBarry Smith /*@
323aa219208SBarry Smith    DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with DMDAVecGetArrayDOF()
32447c6ae99SBarry Smith 
3252ddbf755SMatthew G. Knepley    Logically collective
32647c6ae99SBarry Smith 
32747c6ae99SBarry Smith    Input Parameter:
32847c6ae99SBarry Smith +  da - the distributed array
32947c6ae99SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
330564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
33147c6ae99SBarry Smith -  array - the array
33247c6ae99SBarry Smith 
33347c6ae99SBarry Smith   Level: intermediate
33447c6ae99SBarry Smith 
335fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF(),
336856ac3e9SJed Brown           DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), DMStagVecRestoreArrayDOF(), DMDAVecRestoreArrayDOFRead()
33747c6ae99SBarry Smith @*/
3387087cfbeSBarry Smith PetscErrorCode  DMDAVecRestoreArrayDOF(DM da,Vec vec,void *array)
33947c6ae99SBarry Smith {
34047c6ae99SBarry Smith   PetscErrorCode ierr;
34147c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
34247c6ae99SBarry Smith 
34347c6ae99SBarry Smith   PetscFunctionBegin;
344aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
345aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
346ea78f98cSLisandro Dalcin   ierr = DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
34747c6ae99SBarry Smith 
34847c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
34947c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
35047c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
35147c6ae99SBarry Smith     gxm = xm;
35247c6ae99SBarry Smith     gym = ym;
35347c6ae99SBarry Smith     gzm = zm;
35447c6ae99SBarry Smith     gxs = xs;
35547c6ae99SBarry Smith     gys = ys;
35647c6ae99SBarry Smith     gzs = zs;
35747c6ae99SBarry Smith   }
35847c6ae99SBarry Smith 
35947c6ae99SBarry Smith   if (dim == 1) {
36047c6ae99SBarry Smith     ierr = VecRestoreArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
36147c6ae99SBarry Smith   } else if (dim == 2) {
36247c6ae99SBarry Smith     ierr = VecRestoreArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
36347c6ae99SBarry Smith   } else if (dim == 3) {
36447c6ae99SBarry Smith     ierr = VecRestoreArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
36530729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
36647c6ae99SBarry Smith   PetscFunctionReturn(0);
36747c6ae99SBarry Smith }
36847c6ae99SBarry Smith 
3695edff71fSBarry Smith /*@C
3705edff71fSBarry Smith    DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
3715edff71fSBarry Smith       the underlying vector and is indexed using the global dimensions.
3725edff71fSBarry Smith 
3732ddbf755SMatthew G. Knepley    Not collective
3745edff71fSBarry Smith 
3755edff71fSBarry Smith    Input Parameter:
3765edff71fSBarry Smith +  da - the distributed array
3775edff71fSBarry Smith -  vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
3785edff71fSBarry Smith 
3795edff71fSBarry Smith    Output Parameter:
3805edff71fSBarry Smith .  array - the array
3815edff71fSBarry Smith 
3825edff71fSBarry Smith    Notes:
3835edff71fSBarry Smith     Call DMDAVecRestoreArrayRead() once you have finished accessing the vector entries.
3845edff71fSBarry Smith 
3855edff71fSBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
3865edff71fSBarry Smith 
3875edff71fSBarry Smith     If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
3885edff71fSBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
3895edff71fSBarry Smith 
3905edff71fSBarry Smith     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
3915edff71fSBarry Smith 
39295452b02SPatrick Sanan   Fortran Notes:
39395452b02SPatrick Sanan     From Fortran use DMDAVecGetArrayReadF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
3945edff71fSBarry 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
3955edff71fSBarry Smith        dimension of the DMDA. The order of the indices is array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) (when dof is 1) otherwise
3965edff71fSBarry Smith        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
397af0996ceSBarry Smith        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
3985edff71fSBarry Smith 
3994ebdaf59SBarry Smith   Due to bugs in the compiler DMDAVecGetArrayReadF90() does not work with gfortran versions before 4.5
4005edff71fSBarry Smith 
4015edff71fSBarry Smith   Level: intermediate
4025edff71fSBarry Smith 
403fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArrayRead(), DMDAVecRestoreArrayDOF()
4043bd220d7SPatrick Sanan           DMDAVecGetArrayDOF(), DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(),
4053bd220d7SPatrick Sanan           DMStagVecGetArrayRead()
4065edff71fSBarry Smith @*/
4075edff71fSBarry Smith PetscErrorCode  DMDAVecGetArrayRead(DM da,Vec vec,void *array)
4085edff71fSBarry Smith {
4095edff71fSBarry Smith   PetscErrorCode ierr;
4105edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
4115edff71fSBarry Smith 
4125edff71fSBarry Smith   PetscFunctionBegin;
413a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
4145edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
4155edff71fSBarry Smith   PetscValidPointer(array, 3);
4165edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
4175edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
418ea78f98cSLisandro Dalcin   ierr = DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
4195edff71fSBarry Smith 
4205edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
4215edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
4225edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
4235edff71fSBarry Smith     gxm = xm;
4245edff71fSBarry Smith     gym = ym;
4255edff71fSBarry Smith     gzm = zm;
4265edff71fSBarry Smith     gxs = xs;
4275edff71fSBarry Smith     gys = ys;
4285edff71fSBarry Smith     gzs = zs;
4295edff71fSBarry Smith   } else if (N != gxm*gym*gzm*dof) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMDA local sizes %D %D\n",N,xm*ym*zm*dof,gxm*gym*gzm*dof);
4305edff71fSBarry Smith 
4315edff71fSBarry Smith   if (dim == 1) {
4325edff71fSBarry Smith     ierr = VecGetArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
4335edff71fSBarry Smith   } else if (dim == 2) {
4345edff71fSBarry Smith     ierr = VecGetArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
4355edff71fSBarry Smith   } else if (dim == 3) {
4365edff71fSBarry Smith     ierr = VecGetArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
4375edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
4385edff71fSBarry Smith   PetscFunctionReturn(0);
4395edff71fSBarry Smith }
4405edff71fSBarry Smith 
4415edff71fSBarry Smith /*@
4425edff71fSBarry Smith    DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with DMDAVecGetArrayRead()
4435edff71fSBarry Smith 
4442ddbf755SMatthew G. Knepley    Not collective
4455edff71fSBarry Smith 
4465edff71fSBarry Smith    Input Parameter:
4475edff71fSBarry Smith +  da - the distributed array
4485edff71fSBarry Smith .  vec - the vector, either a vector the same size as one obtained with
4495edff71fSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
4505edff71fSBarry Smith -  array - the array, non-NULL pointer is zeroed
4515edff71fSBarry Smith 
4525edff71fSBarry Smith   Level: intermediate
4535edff71fSBarry Smith 
45495452b02SPatrick Sanan   Fortran Notes:
45595452b02SPatrick Sanan     From Fortran use DMDAVecRestoreArrayReadF90()
4565edff71fSBarry Smith 
457fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArrayRead(),
4583bd220d7SPatrick Sanan           DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(),
4593bd220d7SPatrick Sanan           DMStagVecRestoreArrayRead()
4605edff71fSBarry Smith @*/
4615edff71fSBarry Smith PetscErrorCode  DMDAVecRestoreArrayRead(DM da,Vec vec,void *array)
4625edff71fSBarry Smith {
4635edff71fSBarry Smith   PetscErrorCode ierr;
4645edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
4655edff71fSBarry Smith 
4665edff71fSBarry Smith   PetscFunctionBegin;
467a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
4685edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
4695edff71fSBarry Smith   PetscValidPointer(array, 3);
4705edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
4715edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
472ea78f98cSLisandro Dalcin   ierr = DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
4735edff71fSBarry Smith 
4745edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
4755edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
4765edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
4775edff71fSBarry Smith     gxm = xm;
4785edff71fSBarry Smith     gym = ym;
4795edff71fSBarry Smith     gzm = zm;
4805edff71fSBarry Smith     gxs = xs;
4815edff71fSBarry Smith     gys = ys;
4825edff71fSBarry Smith     gzs = zs;
4835edff71fSBarry Smith   } else if (N != gxm*gym*gzm*dof) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMDA local sizes %D %D\n",N,xm*ym*zm*dof,gxm*gym*gzm*dof);
4845edff71fSBarry Smith 
4855edff71fSBarry Smith   if (dim == 1) {
4865edff71fSBarry Smith     ierr = VecRestoreArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
4875edff71fSBarry Smith   } else if (dim == 2) {
4885edff71fSBarry Smith     ierr = VecRestoreArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
4895edff71fSBarry Smith   } else if (dim == 3) {
4905edff71fSBarry Smith     ierr = VecRestoreArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
4915edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
4925edff71fSBarry Smith   PetscFunctionReturn(0);
4935edff71fSBarry Smith }
4945edff71fSBarry Smith 
4955edff71fSBarry Smith /*@C
4965edff71fSBarry Smith    DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
4975edff71fSBarry Smith       the underlying vector and is indexed using the global dimensions.
4985edff71fSBarry Smith 
4995edff71fSBarry Smith    Not Collective
5005edff71fSBarry Smith 
5015edff71fSBarry Smith    Input Parameter:
5025edff71fSBarry Smith +  da - the distributed array
5035edff71fSBarry Smith -  vec - the vector, either a vector the same size as one obtained with
5045edff71fSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
5055edff71fSBarry Smith 
5065edff71fSBarry Smith    Output Parameter:
5075edff71fSBarry Smith .  array - the array
5085edff71fSBarry Smith 
5095edff71fSBarry Smith    Notes:
5105edff71fSBarry Smith     Call DMDAVecRestoreArrayDOFRead() once you have finished accessing the vector entries.
5115edff71fSBarry Smith 
5125edff71fSBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
5135edff71fSBarry Smith 
5144ebdaf59SBarry Smith     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayReadF90() and declare your array with one higher dimension,
515c4762a1bSJed Brown     see src/dm/tutorials/ex11f90.F
5165edff71fSBarry Smith 
5175edff71fSBarry Smith   Level: intermediate
5185edff71fSBarry Smith 
519856ac3e9SJed Brown .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(),
5203bd220d7SPatrick Sanan           DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), DMStagVecGetArrayDOFRead()
5215edff71fSBarry Smith @*/
5225edff71fSBarry Smith PetscErrorCode  DMDAVecGetArrayDOFRead(DM da,Vec vec,void *array)
5235edff71fSBarry Smith {
5245edff71fSBarry Smith   PetscErrorCode ierr;
5255edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
5265edff71fSBarry Smith 
5275edff71fSBarry Smith   PetscFunctionBegin;
5285edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
5295edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
530ea78f98cSLisandro Dalcin   ierr = DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
5315edff71fSBarry Smith 
5325edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
5335edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
5345edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
5355edff71fSBarry Smith     gxm = xm;
5365edff71fSBarry Smith     gym = ym;
5375edff71fSBarry Smith     gzm = zm;
5385edff71fSBarry Smith     gxs = xs;
5395edff71fSBarry Smith     gys = ys;
5405edff71fSBarry Smith     gzs = zs;
5415edff71fSBarry Smith   } else if (N != gxm*gym*gzm*dof) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMDA local sizes %D %D\n",N,xm*ym*zm*dof,gxm*gym*gzm*dof);
5425edff71fSBarry Smith 
5435edff71fSBarry Smith   if (dim == 1) {
5445edff71fSBarry Smith     ierr = VecGetArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
5455edff71fSBarry Smith   } else if (dim == 2) {
5465edff71fSBarry Smith     ierr = VecGetArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
5475edff71fSBarry Smith   } else if (dim == 3) {
5485edff71fSBarry Smith     ierr = VecGetArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
5495edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
5505edff71fSBarry Smith   PetscFunctionReturn(0);
5515edff71fSBarry Smith }
5525edff71fSBarry Smith 
5535edff71fSBarry Smith /*@
5545edff71fSBarry Smith    DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with DMDAVecGetArrayDOFRead()
5555edff71fSBarry Smith 
5565edff71fSBarry Smith    Not Collective
5575edff71fSBarry Smith 
5585edff71fSBarry Smith    Input Parameter:
5595edff71fSBarry Smith +  da - the distributed array
5605edff71fSBarry Smith .  vec - the vector, either a vector the same size as one obtained with
5615edff71fSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
5625edff71fSBarry Smith -  array - the array
5635edff71fSBarry Smith 
5645edff71fSBarry Smith   Level: intermediate
5655edff71fSBarry Smith 
566fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF(),
5673bd220d7SPatrick Sanan           DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), DMStagVecRestoreArrayDOFRead()
5685edff71fSBarry Smith @*/
5695edff71fSBarry Smith PetscErrorCode  DMDAVecRestoreArrayDOFRead(DM da,Vec vec,void *array)
5705edff71fSBarry Smith {
5715edff71fSBarry Smith   PetscErrorCode ierr;
5725edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
5735edff71fSBarry Smith 
5745edff71fSBarry Smith   PetscFunctionBegin;
5755edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
5765edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
577ea78f98cSLisandro Dalcin   ierr = DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
5785edff71fSBarry Smith 
5795edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
5805edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
5815edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
5825edff71fSBarry Smith     gxm = xm;
5835edff71fSBarry Smith     gym = ym;
5845edff71fSBarry Smith     gzm = zm;
5855edff71fSBarry Smith     gxs = xs;
5865edff71fSBarry Smith     gys = ys;
5875edff71fSBarry Smith     gzs = zs;
5885edff71fSBarry Smith   }
5895edff71fSBarry Smith 
5905edff71fSBarry Smith   if (dim == 1) {
5915edff71fSBarry Smith     ierr = VecRestoreArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
5925edff71fSBarry Smith   } else if (dim == 2) {
5935edff71fSBarry Smith     ierr = VecRestoreArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
5945edff71fSBarry Smith   } else if (dim == 3) {
5955edff71fSBarry Smith     ierr = VecRestoreArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
5965edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
5975edff71fSBarry Smith   PetscFunctionReturn(0);
5985edff71fSBarry Smith }
5995edff71fSBarry Smith 
600*1e5d2365SBarry Smith /*@C
601*1e5d2365SBarry Smith    DMDAVecGetArrayDOFWrite - Returns a multiple dimension array that shares data with
602*1e5d2365SBarry Smith       the underlying vector and is indexed using the global dimensions.
603*1e5d2365SBarry Smith 
604*1e5d2365SBarry Smith    Not Collective
605*1e5d2365SBarry Smith 
606*1e5d2365SBarry Smith    Input Parameter:
607*1e5d2365SBarry Smith +  da - the distributed array
608*1e5d2365SBarry Smith -  vec - the vector, either a vector the same size as one obtained with
609*1e5d2365SBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
610*1e5d2365SBarry Smith 
611*1e5d2365SBarry Smith    Output Parameter:
612*1e5d2365SBarry Smith .  array - the array
613*1e5d2365SBarry Smith 
614*1e5d2365SBarry Smith    Notes:
615*1e5d2365SBarry Smith     Call DMDAVecRestoreArrayDOFWrite() once you have finished accessing the vector entries.
616*1e5d2365SBarry Smith 
617*1e5d2365SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
618*1e5d2365SBarry Smith 
619*1e5d2365SBarry Smith     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayWriteF90() and declare your array with one higher dimension,
620*1e5d2365SBarry Smith     see src/dm/tutorials/ex11f90.F
621*1e5d2365SBarry Smith 
622*1e5d2365SBarry Smith   Level: intermediate
623*1e5d2365SBarry Smith 
624*1e5d2365SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(),
625*1e5d2365SBarry Smith           DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMStagVecGetArrayDOFWrite(),
626*1e5d2365SBarry Smith           DMStagVecRestoreArrayDOFRead(), DMStagVecRestoreArrayDOFRead()
627*1e5d2365SBarry Smith @*/
628*1e5d2365SBarry Smith PetscErrorCode  DMDAVecGetArrayDOFWrite(DM da,Vec vec,void *array)
629*1e5d2365SBarry Smith {
630*1e5d2365SBarry Smith   PetscErrorCode ierr;
631*1e5d2365SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
632*1e5d2365SBarry Smith 
633*1e5d2365SBarry Smith   PetscFunctionBegin;
634*1e5d2365SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
635*1e5d2365SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
636*1e5d2365SBarry Smith   ierr = DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
637*1e5d2365SBarry Smith 
638*1e5d2365SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
639*1e5d2365SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
640*1e5d2365SBarry Smith   if (N == xm*ym*zm*dof) {
641*1e5d2365SBarry Smith     gxm = xm;
642*1e5d2365SBarry Smith     gym = ym;
643*1e5d2365SBarry Smith     gzm = zm;
644*1e5d2365SBarry Smith     gxs = xs;
645*1e5d2365SBarry Smith     gys = ys;
646*1e5d2365SBarry Smith     gzs = zs;
647*1e5d2365SBarry Smith   } else if (N != gxm*gym*gzm*dof) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %D is not compatible with DMDA local sizes %D %D\n",N,xm*ym*zm*dof,gxm*gym*gzm*dof);
648*1e5d2365SBarry Smith 
649*1e5d2365SBarry Smith   if (dim == 1) {
650*1e5d2365SBarry Smith     ierr = VecGetArray2dWrite(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
651*1e5d2365SBarry Smith   } else if (dim == 2) {
652*1e5d2365SBarry Smith     ierr = VecGetArray3dWrite(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
653*1e5d2365SBarry Smith   } else if (dim == 3) {
654*1e5d2365SBarry Smith     ierr = VecGetArray4dWrite(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
655*1e5d2365SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
656*1e5d2365SBarry Smith   PetscFunctionReturn(0);
657*1e5d2365SBarry Smith }
658*1e5d2365SBarry Smith 
659*1e5d2365SBarry Smith /*@
660*1e5d2365SBarry Smith    DMDAVecRestoreArrayDOFWrite - Restores a multiple dimension array obtained with DMDAVecGetArrayDOFWrite()
661*1e5d2365SBarry Smith 
662*1e5d2365SBarry Smith    Not Collective
663*1e5d2365SBarry Smith 
664*1e5d2365SBarry Smith    Input Parameter:
665*1e5d2365SBarry Smith +  da - the distributed array
666*1e5d2365SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
667*1e5d2365SBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
668*1e5d2365SBarry Smith -  array - the array
669*1e5d2365SBarry Smith 
670*1e5d2365SBarry Smith   Level: intermediate
671*1e5d2365SBarry Smith 
672*1e5d2365SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF(),
673*1e5d2365SBarry Smith           DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMStagVecRestoreArrayDOFWrite(),
674*1e5d2365SBarry Smith           DMStagVecRestoreArrayDOFRead(), DMStagVecRestoreArrayDOFRead()
675*1e5d2365SBarry Smith @*/
676*1e5d2365SBarry Smith PetscErrorCode  DMDAVecRestoreArrayDOFWrite(DM da,Vec vec,void *array)
677*1e5d2365SBarry Smith {
678*1e5d2365SBarry Smith   PetscErrorCode ierr;
679*1e5d2365SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
680*1e5d2365SBarry Smith 
681*1e5d2365SBarry Smith   PetscFunctionBegin;
682*1e5d2365SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
683*1e5d2365SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
684*1e5d2365SBarry Smith   ierr = DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
685*1e5d2365SBarry Smith 
686*1e5d2365SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
687*1e5d2365SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
688*1e5d2365SBarry Smith   if (N == xm*ym*zm*dof) {
689*1e5d2365SBarry Smith     gxm = xm;
690*1e5d2365SBarry Smith     gym = ym;
691*1e5d2365SBarry Smith     gzm = zm;
692*1e5d2365SBarry Smith     gxs = xs;
693*1e5d2365SBarry Smith     gys = ys;
694*1e5d2365SBarry Smith     gzs = zs;
695*1e5d2365SBarry Smith   }
696*1e5d2365SBarry Smith 
697*1e5d2365SBarry Smith   if (dim == 1) {
698*1e5d2365SBarry Smith     ierr = VecRestoreArray2dWrite(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
699*1e5d2365SBarry Smith   } else if (dim == 2) {
700*1e5d2365SBarry Smith     ierr = VecRestoreArray3dWrite(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
701*1e5d2365SBarry Smith   } else if (dim == 3) {
702*1e5d2365SBarry Smith     ierr = VecRestoreArray4dWrite(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
703*1e5d2365SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
704*1e5d2365SBarry Smith   PetscFunctionReturn(0);
705*1e5d2365SBarry Smith }
706*1e5d2365SBarry Smith 
707