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