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