xref: /petsc/src/dm/impls/da/dagetarray.c (revision 0f55b88d1ac75fa99e70a53aff055236a607e7f1)
1 
2 #include <petscdmda.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    Not Collective
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(1:dof,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   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
55   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
56   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
57 
58   /* Handle case where user passes in global vector as opposed to local */
59   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
60   if (N == xm*ym*zm*dof) {
61     gxm = xm;
62     gym = ym;
63     gzm = zm;
64     gxs = xs;
65     gys = ys;
66     gzs = zs;
67   } 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);
68 
69   if (dim == 1) {
70     ierr = VecGetArray1d(vec,gxm*dof,gxs*dof,(PetscScalar **)array);CHKERRQ(ierr);
71   } else if (dim == 2) {
72     ierr = VecGetArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
73   } else if (dim == 3) {
74     ierr = VecGetArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
75   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
76 
77   PetscFunctionReturn(0);
78 }
79 
80 #undef __FUNCT__
81 #define __FUNCT__ "DMDAVecRestoreArray"
82 /*@
83    DMDAVecRestoreArray - Restores a multiple dimension array obtained with DMDAVecGetArray()
84 
85    Not Collective
86 
87    Input Parameter:
88 +  da - the distributed array
89 .  vec - the vector, either a vector the same size as one obtained with
90          DMCreateGlobalVector() or DMCreateLocalVector()
91 -  array - the array
92 
93   Level: intermediate
94 
95   Fortran Notes: From Fortran use DMDAVecRestoreArrayF90()
96 
97 .keywords: distributed array, get, corners, nodes, local indices, coordinates
98 
99 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray()
100 @*/
101 PetscErrorCode  DMDAVecRestoreArray(DM da,Vec vec,void *array)
102 {
103   PetscErrorCode ierr;
104   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
105 
106   PetscFunctionBegin;
107   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
108   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
109   PetscValidPointer(array, 3);
110   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
111   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
112   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
113 
114   /* Handle case where user passes in global vector as opposed to local */
115   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
116   if (N == xm*ym*zm*dof) {
117     gxm = xm;
118     gym = ym;
119     gzm = zm;
120     gxs = xs;
121     gys = ys;
122     gzs = zs;
123   } 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);
124 
125   if (dim == 1) {
126     ierr = VecRestoreArray1d(vec,gxm*dof,gxs*dof,(PetscScalar **)array);CHKERRQ(ierr);
127   } else if (dim == 2) {
128     ierr = VecRestoreArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
129   } else if (dim == 3) {
130     ierr = VecRestoreArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
131   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
132   PetscFunctionReturn(0);
133 }
134 
135 #undef __FUNCT__
136 #define __FUNCT__ "DMDAVecGetArrayDOF"
137 /*@C
138    DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
139       the underlying vector and is indexed using the global dimensions.
140 
141    Not Collective
142 
143    Input Parameter:
144 +  da - the distributed array
145 -  vec - the vector, either a vector the same size as one obtained with
146          DMCreateGlobalVector() or DMCreateLocalVector()
147 
148    Output Parameter:
149 .  array - the array
150 
151    Notes:
152     Call DMDAVecRestoreArrayDOF() once you have finished accessing the vector entries.
153 
154     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
155 
156   Level: intermediate
157 
158 .keywords: distributed array, get, corners, nodes, local indices, coordinates
159 
160 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF()
161 @*/
162 PetscErrorCode  DMDAVecGetArrayDOF(DM da,Vec vec,void *array)
163 {
164   PetscErrorCode ierr;
165   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
166 
167   PetscFunctionBegin;
168   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
169   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
170   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
171 
172   /* Handle case where user passes in global vector as opposed to local */
173   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
174   if (N == xm*ym*zm*dof) {
175     gxm = xm;
176     gym = ym;
177     gzm = zm;
178     gxs = xs;
179     gys = ys;
180     gzs = zs;
181   } 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);
182 
183   if (dim == 1) {
184     ierr = VecGetArray2d(vec,gxm,dof,gxs,0,(PetscScalar ***)array);CHKERRQ(ierr);
185   } else if (dim == 2) {
186     ierr = VecGetArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
187   } else if (dim == 3) {
188     ierr = VecGetArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
189   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
190   PetscFunctionReturn(0);
191 }
192 
193 #undef __FUNCT__
194 #define __FUNCT__ "DMDAVecRestoreArrayDOF"
195 /*@
196    DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with DMDAVecGetArrayDOF()
197 
198    Not Collective
199 
200    Input Parameter:
201 +  da - the distributed array
202 .  vec - the vector, either a vector the same size as one obtained with
203          DMCreateGlobalVector() or DMCreateLocalVector()
204 -  array - the array
205 
206   Level: intermediate
207 
208 .keywords: distributed array, get, corners, nodes, local indices, coordinates
209 
210 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF()
211 @*/
212 PetscErrorCode  DMDAVecRestoreArrayDOF(DM da,Vec vec,void *array)
213 {
214   PetscErrorCode ierr;
215   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
216 
217   PetscFunctionBegin;
218   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
219   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
220   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
221 
222   /* Handle case where user passes in global vector as opposed to local */
223   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
224   if (N == xm*ym*zm*dof) {
225     gxm = xm;
226     gym = ym;
227     gzm = zm;
228     gxs = xs;
229     gys = ys;
230     gzs = zs;
231   }
232 
233   if (dim == 1) {
234     ierr = VecRestoreArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
235   } else if (dim == 2) {
236     ierr = VecRestoreArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
237   } else if (dim == 3) {
238     ierr = VecRestoreArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
239   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
240   PetscFunctionReturn(0);
241 }
242 
243 
244 
245 
246 
247 
248 
249 
250 
251 
252 
253 
254 
255