xref: /petsc/src/dm/impls/da/dagetarray.c (revision 1b82215ed1dc7986f1727aad36942f6eeb755622)
147c6ae99SBarry Smith 
2c6db04a5SJed Brown #include <petscdmda.h>    /*I   "petscdmda.h"   I*/
347c6ae99SBarry Smith 
447c6ae99SBarry Smith #undef __FUNCT__
5aa219208SBarry Smith #define __FUNCT__ "DMDAVecGetArray"
647c6ae99SBarry Smith /*@C
7aa219208SBarry Smith    DMDAVecGetArray - Returns a multiple dimension array that shares data with
847c6ae99SBarry Smith       the underlying vector and is indexed using the global dimensions.
947c6ae99SBarry Smith 
1047c6ae99SBarry Smith    Not Collective
1147c6ae99SBarry Smith 
1247c6ae99SBarry Smith    Input Parameter:
1347c6ae99SBarry Smith +  da - the distributed array
1447c6ae99SBarry Smith -  vec - the vector, either a vector the same size as one obtained with
15564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
1647c6ae99SBarry Smith 
1747c6ae99SBarry Smith    Output Parameter:
1847c6ae99SBarry Smith .  array - the array
1947c6ae99SBarry Smith 
2047c6ae99SBarry Smith    Notes:
21aa219208SBarry Smith     Call DMDAVecRestoreArray() once you have finished accessing the vector entries.
2247c6ae99SBarry Smith 
2347c6ae99SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
2447c6ae99SBarry Smith 
25564755cdSBarry Smith     If vec is a local vector (obtained with DMCreateLocalVector() etc) then they ghost point locations are accessable. If it is
2647c6ae99SBarry Smith     a global vector then the ghost points are not accessable. Of course with the local vector you will have had to do the
2747c6ae99SBarry Smith 
289a42bb27SBarry Smith     appropriate DMLocalToGlobalBegin() and DMLocalToGlobalEnd() to have correct values in the ghost locations.
2947c6ae99SBarry Smith 
30aa219208SBarry Smith   Fortran Notes: From Fortran use DMDAVecGetArrayF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
31aa219208SBarry Smith        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
32aa219208SBarry Smith        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
3347c6ae99SBarry Smith        array(1:dof,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
343c48a1e8SJed Brown        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include finclude/petscdmda.h90 to access this routine.
3547c6ae99SBarry Smith 
36596f5af7SJed Brown   Due to bugs in the compiler DMDAVecGetArrayF90() does not work with gfortran versions before 4.5
3747c6ae99SBarry Smith 
3847c6ae99SBarry Smith   Level: intermediate
3947c6ae99SBarry Smith 
4047c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
4147c6ae99SBarry Smith 
42aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF()
436db82c94SMatthew G Knepley           DMDAVecGetArrayDOF()
4447c6ae99SBarry Smith @*/
457087cfbeSBarry Smith PetscErrorCode  DMDAVecGetArray(DM da,Vec vec,void *array)
4647c6ae99SBarry Smith {
4747c6ae99SBarry Smith   PetscErrorCode ierr;
4847c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
4947c6ae99SBarry Smith 
5047c6ae99SBarry Smith   PetscFunctionBegin;
516db82c94SMatthew G Knepley   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
526db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
536db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
54aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
55aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
561321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
5747c6ae99SBarry Smith 
5847c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
5947c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
6047c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
6147c6ae99SBarry Smith     gxm = xm;
6247c6ae99SBarry Smith     gym = ym;
6347c6ae99SBarry Smith     gzm = zm;
6447c6ae99SBarry Smith     gxs = xs;
6547c6ae99SBarry Smith     gys = ys;
6647c6ae99SBarry Smith     gzs = zs;
6730729d88SBarry Smith   } 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);
6847c6ae99SBarry Smith 
6947c6ae99SBarry Smith   if (dim == 1) {
7047c6ae99SBarry Smith     ierr = VecGetArray1d(vec,gxm*dof,gxs*dof,(PetscScalar **)array);CHKERRQ(ierr);
7147c6ae99SBarry Smith   } else if (dim == 2) {
7247c6ae99SBarry Smith     ierr = VecGetArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
7347c6ae99SBarry Smith   } else if (dim == 3) {
7447c6ae99SBarry Smith     ierr = VecGetArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
7530729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
7647c6ae99SBarry Smith 
7747c6ae99SBarry Smith   PetscFunctionReturn(0);
7847c6ae99SBarry Smith }
7947c6ae99SBarry Smith 
8047c6ae99SBarry Smith #undef __FUNCT__
81aa219208SBarry Smith #define __FUNCT__ "DMDAVecRestoreArray"
8247c6ae99SBarry Smith /*@
83aa219208SBarry Smith    DMDAVecRestoreArray - Restores a multiple dimension array obtained with DMDAVecGetArray()
8447c6ae99SBarry Smith 
8547c6ae99SBarry Smith    Not Collective
8647c6ae99SBarry Smith 
8747c6ae99SBarry Smith    Input Parameter:
8847c6ae99SBarry Smith +  da - the distributed array
8947c6ae99SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
90564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
9147c6ae99SBarry Smith -  array - the array
9247c6ae99SBarry Smith 
9347c6ae99SBarry Smith   Level: intermediate
9447c6ae99SBarry Smith 
95aa219208SBarry Smith   Fortran Notes: From Fortran use DMDAVecRestoreArrayF90()
9647c6ae99SBarry Smith 
9747c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
9847c6ae99SBarry Smith 
99aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray()
10047c6ae99SBarry Smith @*/
1017087cfbeSBarry Smith PetscErrorCode  DMDAVecRestoreArray(DM da,Vec vec,void *array)
10247c6ae99SBarry Smith {
10347c6ae99SBarry Smith   PetscErrorCode ierr;
10447c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
10547c6ae99SBarry Smith 
10647c6ae99SBarry Smith   PetscFunctionBegin;
1076db82c94SMatthew G Knepley   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
1086db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
1096db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
110aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
111aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
1121321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
11347c6ae99SBarry Smith 
11447c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
11547c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
11647c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
11747c6ae99SBarry Smith     gxm = xm;
11847c6ae99SBarry Smith     gym = ym;
11947c6ae99SBarry Smith     gzm = zm;
12047c6ae99SBarry Smith     gxs = xs;
12147c6ae99SBarry Smith     gys = ys;
12247c6ae99SBarry Smith     gzs = zs;
12330729d88SBarry Smith   } 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);
12447c6ae99SBarry Smith 
12547c6ae99SBarry Smith   if (dim == 1) {
12647c6ae99SBarry Smith     ierr = VecRestoreArray1d(vec,gxm*dof,gxs*dof,(PetscScalar **)array);CHKERRQ(ierr);
12747c6ae99SBarry Smith   } else if (dim == 2) {
12847c6ae99SBarry Smith     ierr = VecRestoreArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
12947c6ae99SBarry Smith   } else if (dim == 3) {
13047c6ae99SBarry Smith     ierr = VecRestoreArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
13130729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
13247c6ae99SBarry Smith   PetscFunctionReturn(0);
13347c6ae99SBarry Smith }
13447c6ae99SBarry Smith 
13547c6ae99SBarry Smith #undef __FUNCT__
136aa219208SBarry Smith #define __FUNCT__ "DMDAVecGetArrayDOF"
13747c6ae99SBarry Smith /*@C
138aa219208SBarry Smith    DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
13947c6ae99SBarry Smith       the underlying vector and is indexed using the global dimensions.
14047c6ae99SBarry Smith 
14147c6ae99SBarry Smith    Not Collective
14247c6ae99SBarry Smith 
14347c6ae99SBarry Smith    Input Parameter:
14447c6ae99SBarry Smith +  da - the distributed array
14547c6ae99SBarry Smith -  vec - the vector, either a vector the same size as one obtained with
146564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
14747c6ae99SBarry Smith 
14847c6ae99SBarry Smith    Output Parameter:
14947c6ae99SBarry Smith .  array - the array
15047c6ae99SBarry Smith 
15147c6ae99SBarry Smith    Notes:
152aa219208SBarry Smith     Call DMDAVecRestoreArrayDOF() once you have finished accessing the vector entries.
15347c6ae99SBarry Smith 
15447c6ae99SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
15547c6ae99SBarry Smith 
156*1b82215eSBarry Smith     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayF90() and declare your array with one higher dimension,
157*1b82215eSBarry Smith     see src/dm/examples/tutorials/ex11f90.F
158*1b82215eSBarry Smith 
15947c6ae99SBarry Smith   Level: intermediate
16047c6ae99SBarry Smith 
16147c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
16247c6ae99SBarry Smith 
163aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF()
16447c6ae99SBarry Smith @*/
1657087cfbeSBarry Smith PetscErrorCode  DMDAVecGetArrayDOF(DM da,Vec vec,void *array)
16647c6ae99SBarry Smith {
16747c6ae99SBarry Smith   PetscErrorCode ierr;
16847c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
16947c6ae99SBarry Smith 
17047c6ae99SBarry Smith   PetscFunctionBegin;
171aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
172aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
1731321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
17447c6ae99SBarry Smith 
17547c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
17647c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
17747c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
17847c6ae99SBarry Smith     gxm = xm;
17947c6ae99SBarry Smith     gym = ym;
18047c6ae99SBarry Smith     gzm = zm;
18147c6ae99SBarry Smith     gxs = xs;
18247c6ae99SBarry Smith     gys = ys;
18347c6ae99SBarry Smith     gzs = zs;
18430729d88SBarry Smith   } 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);
18547c6ae99SBarry Smith 
18647c6ae99SBarry Smith   if (dim == 1) {
18747c6ae99SBarry Smith     ierr = VecGetArray2d(vec,gxm,dof,gxs,0,(PetscScalar ***)array);CHKERRQ(ierr);
18847c6ae99SBarry Smith   } else if (dim == 2) {
18947c6ae99SBarry Smith     ierr = VecGetArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
19047c6ae99SBarry Smith   } else if (dim == 3) {
19147c6ae99SBarry Smith     ierr = VecGetArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
19230729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
19347c6ae99SBarry Smith   PetscFunctionReturn(0);
19447c6ae99SBarry Smith }
19547c6ae99SBarry Smith 
19647c6ae99SBarry Smith #undef __FUNCT__
197aa219208SBarry Smith #define __FUNCT__ "DMDAVecRestoreArrayDOF"
19847c6ae99SBarry Smith /*@
199aa219208SBarry Smith    DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with DMDAVecGetArrayDOF()
20047c6ae99SBarry Smith 
20147c6ae99SBarry Smith    Not Collective
20247c6ae99SBarry Smith 
20347c6ae99SBarry Smith    Input Parameter:
20447c6ae99SBarry Smith +  da - the distributed array
20547c6ae99SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
206564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
20747c6ae99SBarry Smith -  array - the array
20847c6ae99SBarry Smith 
20947c6ae99SBarry Smith   Level: intermediate
21047c6ae99SBarry Smith 
21147c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
21247c6ae99SBarry Smith 
213aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF()
21447c6ae99SBarry Smith @*/
2157087cfbeSBarry Smith PetscErrorCode  DMDAVecRestoreArrayDOF(DM da,Vec vec,void *array)
21647c6ae99SBarry Smith {
21747c6ae99SBarry Smith   PetscErrorCode ierr;
21847c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
21947c6ae99SBarry Smith 
22047c6ae99SBarry Smith   PetscFunctionBegin;
221aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
222aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
2231321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
22447c6ae99SBarry Smith 
22547c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
22647c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
22747c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
22847c6ae99SBarry Smith     gxm = xm;
22947c6ae99SBarry Smith     gym = ym;
23047c6ae99SBarry Smith     gzm = zm;
23147c6ae99SBarry Smith     gxs = xs;
23247c6ae99SBarry Smith     gys = ys;
23347c6ae99SBarry Smith     gzs = zs;
23447c6ae99SBarry Smith   }
23547c6ae99SBarry Smith 
23647c6ae99SBarry Smith   if (dim == 1) {
23747c6ae99SBarry Smith     ierr = VecRestoreArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
23847c6ae99SBarry Smith   } else if (dim == 2) {
23947c6ae99SBarry Smith     ierr = VecRestoreArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
24047c6ae99SBarry Smith   } else if (dim == 3) {
24147c6ae99SBarry Smith     ierr = VecRestoreArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
24230729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
24347c6ae99SBarry Smith   PetscFunctionReturn(0);
24447c6ae99SBarry Smith }
24547c6ae99SBarry Smith 
24647c6ae99SBarry Smith 
24747c6ae99SBarry Smith 
24847c6ae99SBarry Smith 
24947c6ae99SBarry Smith 
25047c6ae99SBarry Smith 
25147c6ae99SBarry Smith 
25247c6ae99SBarry Smith 
25347c6ae99SBarry Smith 
25447c6ae99SBarry Smith 
25547c6ae99SBarry Smith 
25647c6ae99SBarry Smith 
25747c6ae99SBarry Smith 
258