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