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