xref: /petsc/src/dm/impls/da/dagetarray.c (revision 9895aa37ac365bac650f6bd8bf977519f7222510)
1 
2 #include <petsc-private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
3 
4 #undef __FUNCT__
5 #define __FUNCT__ "DMDAVecGetArray"
6 /*@C
7    DMDAVecGetArray - Returns a multiple dimension array that shares data with
8       the underlying vector and is indexed using the global dimensions.
9 
10    Collective on Vec
11 
12    Input Parameter:
13 +  da - the distributed array
14 -  vec - the vector, either a vector the same size as one obtained with
15          DMCreateGlobalVector() or DMCreateLocalVector()
16 
17    Output Parameter:
18 .  array - the array
19 
20    Notes:
21     Call DMDAVecRestoreArray() once you have finished accessing the vector entries.
22 
23     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
24 
25     If vec is a local vector (obtained with DMCreateLocalVector() etc) then they ghost point locations are accessable. If it is
26     a global vector then the ghost points are not accessable. Of course with the local vector you will have had to do the
27 
28     appropriate DMLocalToGlobalBegin() and DMLocalToGlobalEnd() to have correct values in the ghost locations.
29 
30   Fortran Notes: From Fortran use DMDAVecGetArrayF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
31        dimension. For a DMDA created with a dof of 1 use the dimension of the DMDA, for a DMDA created with a dof greater than 1 use one more than the
32        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
33        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
34        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include finclude/petscdmda.h90 to access this routine.
35 
36   Due to bugs in the compiler DMDAVecGetArrayF90() does not work with gfortran versions before 4.5
37 
38   Level: intermediate
39 
40 .keywords: distributed array, get, corners, nodes, local indices, coordinates
41 
42 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF()
43           DMDAVecGetArrayDOF()
44 @*/
45 PetscErrorCode  DMDAVecGetArray(DM da,Vec vec,void *array)
46 {
47   PetscErrorCode ierr;
48   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
49 
50   PetscFunctionBegin;
51   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
52   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
53   PetscValidPointer(array, 3);
54   if (da->defaultSection) {
55     ierr = VecGetArray(vec,(PetscScalar**)array);CHKERRQ(ierr);
56     PetscFunctionReturn(0);
57   }
58   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
59   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
60   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
61 
62   /* Handle case where user passes in global vector as opposed to local */
63   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
64   if (N == xm*ym*zm*dof) {
65     gxm = xm;
66     gym = ym;
67     gzm = zm;
68     gxs = xs;
69     gys = ys;
70     gzs = zs;
71   } 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);
72 
73   if (dim == 1) {
74     ierr = VecGetArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
75   } else if (dim == 2) {
76     ierr = VecGetArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
77   } else if (dim == 3) {
78     ierr = VecGetArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
79   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
80   PetscFunctionReturn(0);
81 }
82 
83 #undef __FUNCT__
84 #define __FUNCT__ "DMDAVecRestoreArray"
85 /*@
86    DMDAVecRestoreArray - Restores a multiple dimension array obtained with DMDAVecGetArray()
87 
88    Collective on Vec
89 
90    Input Parameter:
91 +  da - the distributed array
92 .  vec - the vector, either a vector the same size as one obtained with
93          DMCreateGlobalVector() or DMCreateLocalVector()
94 -  array - the array
95 
96   Level: intermediate
97 
98   Fortran Notes: From Fortran use DMDAVecRestoreArrayF90()
99 
100 .keywords: distributed array, get, corners, nodes, local indices, coordinates
101 
102 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray()
103 @*/
104 PetscErrorCode  DMDAVecRestoreArray(DM da,Vec vec,void *array)
105 {
106   PetscErrorCode ierr;
107   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
108 
109   PetscFunctionBegin;
110   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
111   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
112   PetscValidPointer(array, 3);
113   if (da->defaultSection) {
114     ierr = VecRestoreArray(vec,(PetscScalar**)array);CHKERRQ(ierr);
115     PetscFunctionReturn(0);
116   }
117   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
118   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
119   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
120 
121   /* Handle case where user passes in global vector as opposed to local */
122   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
123   if (N == xm*ym*zm*dof) {
124     gxm = xm;
125     gym = ym;
126     gzm = zm;
127     gxs = xs;
128     gys = ys;
129     gzs = zs;
130   } 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);
131 
132   if (dim == 1) {
133     ierr = VecRestoreArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
134   } else if (dim == 2) {
135     ierr = VecRestoreArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
136   } else if (dim == 3) {
137     ierr = VecRestoreArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
138   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
139   PetscFunctionReturn(0);
140 }
141 
142 #undef __FUNCT__
143 #define __FUNCT__ "DMDAVecGetArrayDOF"
144 /*@C
145    DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
146       the underlying vector and is indexed using the global dimensions.
147 
148    Not Collective
149 
150    Input Parameter:
151 +  da - the distributed array
152 -  vec - the vector, either a vector the same size as one obtained with
153          DMCreateGlobalVector() or DMCreateLocalVector()
154 
155    Output Parameter:
156 .  array - the array
157 
158    Notes:
159     Call DMDAVecRestoreArrayDOF() once you have finished accessing the vector entries.
160 
161     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
162 
163     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayF90() and declare your array with one higher dimension,
164     see src/dm/examples/tutorials/ex11f90.F
165 
166   Level: intermediate
167 
168 .keywords: distributed array, get, corners, nodes, local indices, coordinates
169 
170 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF()
171 @*/
172 PetscErrorCode  DMDAVecGetArrayDOF(DM da,Vec vec,void *array)
173 {
174   PetscErrorCode ierr;
175   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
176 
177   PetscFunctionBegin;
178   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
179   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
180   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
181 
182   /* Handle case where user passes in global vector as opposed to local */
183   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
184   if (N == xm*ym*zm*dof) {
185     gxm = xm;
186     gym = ym;
187     gzm = zm;
188     gxs = xs;
189     gys = ys;
190     gzs = zs;
191   } 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);
192 
193   if (dim == 1) {
194     ierr = VecGetArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
195   } else if (dim == 2) {
196     ierr = VecGetArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
197   } else if (dim == 3) {
198     ierr = VecGetArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
199   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
200   PetscFunctionReturn(0);
201 }
202 
203 #undef __FUNCT__
204 #define __FUNCT__ "DMDAVecRestoreArrayDOF"
205 /*@
206    DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with DMDAVecGetArrayDOF()
207 
208    Not Collective
209 
210    Input Parameter:
211 +  da - the distributed array
212 .  vec - the vector, either a vector the same size as one obtained with
213          DMCreateGlobalVector() or DMCreateLocalVector()
214 -  array - the array
215 
216   Level: intermediate
217 
218 .keywords: distributed array, get, corners, nodes, local indices, coordinates
219 
220 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF()
221 @*/
222 PetscErrorCode  DMDAVecRestoreArrayDOF(DM da,Vec vec,void *array)
223 {
224   PetscErrorCode ierr;
225   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
226 
227   PetscFunctionBegin;
228   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
229   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
230   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
231 
232   /* Handle case where user passes in global vector as opposed to local */
233   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
234   if (N == xm*ym*zm*dof) {
235     gxm = xm;
236     gym = ym;
237     gzm = zm;
238     gxs = xs;
239     gys = ys;
240     gzs = zs;
241   }
242 
243   if (dim == 1) {
244     ierr = VecRestoreArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
245   } else if (dim == 2) {
246     ierr = VecRestoreArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
247   } else if (dim == 3) {
248     ierr = VecRestoreArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
249   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
250   PetscFunctionReturn(0);
251 }
252 
253 
254 
255 
256 
257 
258 
259 
260 
261 
262 
263 
264 
265