xref: /petsc/src/dm/impls/da/dagetarray.c (revision d083f849a86f1f43e18d534ee43954e2786cb29a)
147c6ae99SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.h>    /*I   "petscdmda.h"   I*/
347c6ae99SBarry Smith 
447c6ae99SBarry Smith /*@C
5aa219208SBarry Smith    DMDAVecGetArray - Returns a multiple dimension array that shares data with
647c6ae99SBarry Smith       the underlying vector and is indexed using the global dimensions.
747c6ae99SBarry Smith 
8*d083f849SBarry Smith    Logically collective on da
947c6ae99SBarry Smith 
1047c6ae99SBarry Smith    Input Parameter:
1147c6ae99SBarry Smith +  da - the distributed array
127e57d48aSBarry Smith -  vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
1347c6ae99SBarry Smith 
1447c6ae99SBarry Smith    Output Parameter:
1547c6ae99SBarry Smith .  array - the array
1647c6ae99SBarry Smith 
1747c6ae99SBarry Smith    Notes:
18aa219208SBarry Smith     Call DMDAVecRestoreArray() once you have finished accessing the vector entries.
1947c6ae99SBarry Smith 
2047c6ae99SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
2147c6ae99SBarry Smith 
227e57d48aSBarry Smith     If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
237e57d48aSBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
2447c6ae99SBarry Smith 
257e57d48aSBarry Smith     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
2647c6ae99SBarry Smith 
2795452b02SPatrick Sanan   Fortran Notes:
2895452b02SPatrick Sanan     From Fortran use DMDAVecGetArrayF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
29aa219208SBarry 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
30aa219208SBarry 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
31bef27881SBarry Smith        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
32af0996ceSBarry Smith        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
3347c6ae99SBarry Smith 
34596f5af7SJed Brown   Due to bugs in the compiler DMDAVecGetArrayF90() does not work with gfortran versions before 4.5
3547c6ae99SBarry Smith 
3647c6ae99SBarry Smith   Level: intermediate
3747c6ae99SBarry Smith 
38aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF()
396db82c94SMatthew G Knepley           DMDAVecGetArrayDOF()
4047c6ae99SBarry Smith @*/
417087cfbeSBarry Smith PetscErrorCode  DMDAVecGetArray(DM da,Vec vec,void *array)
4247c6ae99SBarry Smith {
4347c6ae99SBarry Smith   PetscErrorCode ierr;
4447c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
4547c6ae99SBarry Smith 
4647c6ae99SBarry Smith   PetscFunctionBegin;
47a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
486db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
496db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
50aa1993deSMatthew G Knepley   if (da->defaultSection) {
51aa1993deSMatthew G Knepley     ierr = VecGetArray(vec,(PetscScalar**)array);CHKERRQ(ierr);
52aa1993deSMatthew G Knepley     PetscFunctionReturn(0);
53aa1993deSMatthew G Knepley   }
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   PetscFunctionReturn(0);
7747c6ae99SBarry Smith }
7847c6ae99SBarry Smith 
7947c6ae99SBarry Smith /*@
80aa219208SBarry Smith    DMDAVecRestoreArray - Restores a multiple dimension array obtained with DMDAVecGetArray()
8147c6ae99SBarry Smith 
82*d083f849SBarry Smith    Logically collective on da
8347c6ae99SBarry Smith 
8447c6ae99SBarry Smith    Input Parameter:
8547c6ae99SBarry Smith +  da - the distributed array
8647c6ae99SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
87564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
88e5c84f05SJed Brown -  array - the array, non-NULL pointer is zeroed
8947c6ae99SBarry Smith 
9047c6ae99SBarry Smith   Level: intermediate
9147c6ae99SBarry Smith 
9295452b02SPatrick Sanan   Fortran Notes:
9395452b02SPatrick Sanan     From Fortran use DMDAVecRestoreArrayF90()
9447c6ae99SBarry Smith 
95aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray()
9647c6ae99SBarry Smith @*/
977087cfbeSBarry Smith PetscErrorCode  DMDAVecRestoreArray(DM da,Vec vec,void *array)
9847c6ae99SBarry Smith {
9947c6ae99SBarry Smith   PetscErrorCode ierr;
10047c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
10147c6ae99SBarry Smith 
10247c6ae99SBarry Smith   PetscFunctionBegin;
103a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
1046db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
1056db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
106aa1993deSMatthew G Knepley   if (da->defaultSection) {
107aa1993deSMatthew G Knepley     ierr = VecRestoreArray(vec,(PetscScalar**)array);CHKERRQ(ierr);
108aa1993deSMatthew G Knepley     PetscFunctionReturn(0);
109aa1993deSMatthew G Knepley   }
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 /*@C
136aa219208SBarry Smith    DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
13747c6ae99SBarry Smith       the underlying vector and is indexed using the global dimensions.
13847c6ae99SBarry Smith 
1392ddbf755SMatthew G. Knepley    Logically collective
14047c6ae99SBarry Smith 
14147c6ae99SBarry Smith    Input Parameter:
14247c6ae99SBarry Smith +  da - the distributed array
14347c6ae99SBarry Smith -  vec - the vector, either a vector the same size as one obtained with
144564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
14547c6ae99SBarry Smith 
14647c6ae99SBarry Smith    Output Parameter:
14747c6ae99SBarry Smith .  array - the array
14847c6ae99SBarry Smith 
14947c6ae99SBarry Smith    Notes:
150aa219208SBarry Smith     Call DMDAVecRestoreArrayDOF() once you have finished accessing the vector entries.
15147c6ae99SBarry Smith 
15247c6ae99SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
15347c6ae99SBarry Smith 
1541b82215eSBarry Smith     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayF90() and declare your array with one higher dimension,
1551b82215eSBarry Smith     see src/dm/examples/tutorials/ex11f90.F
1561b82215eSBarry Smith 
15747c6ae99SBarry Smith   Level: intermediate
15847c6ae99SBarry Smith 
159aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF()
16047c6ae99SBarry Smith @*/
1617087cfbeSBarry Smith PetscErrorCode  DMDAVecGetArrayDOF(DM da,Vec vec,void *array)
16247c6ae99SBarry Smith {
16347c6ae99SBarry Smith   PetscErrorCode ierr;
16447c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
16547c6ae99SBarry Smith 
16647c6ae99SBarry Smith   PetscFunctionBegin;
167aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
168aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
1691321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
17047c6ae99SBarry Smith 
17147c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
17247c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
17347c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
17447c6ae99SBarry Smith     gxm = xm;
17547c6ae99SBarry Smith     gym = ym;
17647c6ae99SBarry Smith     gzm = zm;
17747c6ae99SBarry Smith     gxs = xs;
17847c6ae99SBarry Smith     gys = ys;
17947c6ae99SBarry Smith     gzs = zs;
18030729d88SBarry 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);
18147c6ae99SBarry Smith 
18247c6ae99SBarry Smith   if (dim == 1) {
18347c6ae99SBarry Smith     ierr = VecGetArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
18447c6ae99SBarry Smith   } else if (dim == 2) {
18547c6ae99SBarry Smith     ierr = VecGetArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
18647c6ae99SBarry Smith   } else if (dim == 3) {
18747c6ae99SBarry Smith     ierr = VecGetArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
18830729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
18947c6ae99SBarry Smith   PetscFunctionReturn(0);
19047c6ae99SBarry Smith }
19147c6ae99SBarry Smith 
19247c6ae99SBarry Smith /*@
193aa219208SBarry Smith    DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with DMDAVecGetArrayDOF()
19447c6ae99SBarry Smith 
1952ddbf755SMatthew G. Knepley    Logically collective
19647c6ae99SBarry Smith 
19747c6ae99SBarry Smith    Input Parameter:
19847c6ae99SBarry Smith +  da - the distributed array
19947c6ae99SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
200564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
20147c6ae99SBarry Smith -  array - the array
20247c6ae99SBarry Smith 
20347c6ae99SBarry Smith   Level: intermediate
20447c6ae99SBarry Smith 
205aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF()
20647c6ae99SBarry Smith @*/
2077087cfbeSBarry Smith PetscErrorCode  DMDAVecRestoreArrayDOF(DM da,Vec vec,void *array)
20847c6ae99SBarry Smith {
20947c6ae99SBarry Smith   PetscErrorCode ierr;
21047c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
21147c6ae99SBarry Smith 
21247c6ae99SBarry Smith   PetscFunctionBegin;
213aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
214aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
2151321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
21647c6ae99SBarry Smith 
21747c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
21847c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
21947c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
22047c6ae99SBarry Smith     gxm = xm;
22147c6ae99SBarry Smith     gym = ym;
22247c6ae99SBarry Smith     gzm = zm;
22347c6ae99SBarry Smith     gxs = xs;
22447c6ae99SBarry Smith     gys = ys;
22547c6ae99SBarry Smith     gzs = zs;
22647c6ae99SBarry Smith   }
22747c6ae99SBarry Smith 
22847c6ae99SBarry Smith   if (dim == 1) {
22947c6ae99SBarry Smith     ierr = VecRestoreArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
23047c6ae99SBarry Smith   } else if (dim == 2) {
23147c6ae99SBarry Smith     ierr = VecRestoreArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
23247c6ae99SBarry Smith   } else if (dim == 3) {
23347c6ae99SBarry Smith     ierr = VecRestoreArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
23430729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
23547c6ae99SBarry Smith   PetscFunctionReturn(0);
23647c6ae99SBarry Smith }
23747c6ae99SBarry Smith 
2385edff71fSBarry Smith /*@C
2395edff71fSBarry Smith    DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
2405edff71fSBarry Smith       the underlying vector and is indexed using the global dimensions.
2415edff71fSBarry Smith 
2422ddbf755SMatthew G. Knepley    Not collective
2435edff71fSBarry Smith 
2445edff71fSBarry Smith    Input Parameter:
2455edff71fSBarry Smith +  da - the distributed array
2465edff71fSBarry Smith -  vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
2475edff71fSBarry Smith 
2485edff71fSBarry Smith    Output Parameter:
2495edff71fSBarry Smith .  array - the array
2505edff71fSBarry Smith 
2515edff71fSBarry Smith    Notes:
2525edff71fSBarry Smith     Call DMDAVecRestoreArrayRead() once you have finished accessing the vector entries.
2535edff71fSBarry Smith 
2545edff71fSBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
2555edff71fSBarry Smith 
2565edff71fSBarry Smith     If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
2575edff71fSBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
2585edff71fSBarry Smith 
2595edff71fSBarry Smith     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
2605edff71fSBarry Smith 
26195452b02SPatrick Sanan   Fortran Notes:
26295452b02SPatrick Sanan     From Fortran use DMDAVecGetArrayReadF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
2635edff71fSBarry 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
2645edff71fSBarry 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
2655edff71fSBarry Smith        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
266af0996ceSBarry Smith        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
2675edff71fSBarry Smith 
2684ebdaf59SBarry Smith   Due to bugs in the compiler DMDAVecGetArrayReadF90() does not work with gfortran versions before 4.5
2695edff71fSBarry Smith 
2705edff71fSBarry Smith   Level: intermediate
2715edff71fSBarry Smith 
2725edff71fSBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF()
2735edff71fSBarry Smith           DMDAVecGetArrayDOF()
2745edff71fSBarry Smith @*/
2755edff71fSBarry Smith PetscErrorCode  DMDAVecGetArrayRead(DM da,Vec vec,void *array)
2765edff71fSBarry Smith {
2775edff71fSBarry Smith   PetscErrorCode ierr;
2785edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
2795edff71fSBarry Smith 
2805edff71fSBarry Smith   PetscFunctionBegin;
281a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
2825edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
2835edff71fSBarry Smith   PetscValidPointer(array, 3);
2845edff71fSBarry Smith   if (da->defaultSection) {
2855edff71fSBarry Smith     ierr = VecGetArrayRead(vec,(const PetscScalar**)array);CHKERRQ(ierr);
2865edff71fSBarry Smith     PetscFunctionReturn(0);
2875edff71fSBarry Smith   }
2885edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
2895edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
2905edff71fSBarry Smith   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
2915edff71fSBarry Smith 
2925edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
2935edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
2945edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
2955edff71fSBarry Smith     gxm = xm;
2965edff71fSBarry Smith     gym = ym;
2975edff71fSBarry Smith     gzm = zm;
2985edff71fSBarry Smith     gxs = xs;
2995edff71fSBarry Smith     gys = ys;
3005edff71fSBarry Smith     gzs = zs;
3015edff71fSBarry 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);
3025edff71fSBarry Smith 
3035edff71fSBarry Smith   if (dim == 1) {
3045edff71fSBarry Smith     ierr = VecGetArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
3055edff71fSBarry Smith   } else if (dim == 2) {
3065edff71fSBarry Smith     ierr = VecGetArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
3075edff71fSBarry Smith   } else if (dim == 3) {
3085edff71fSBarry Smith     ierr = VecGetArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
3095edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
3105edff71fSBarry Smith   PetscFunctionReturn(0);
3115edff71fSBarry Smith }
3125edff71fSBarry Smith 
3135edff71fSBarry Smith /*@
3145edff71fSBarry Smith    DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with DMDAVecGetArrayRead()
3155edff71fSBarry Smith 
3162ddbf755SMatthew G. Knepley    Not collective
3175edff71fSBarry Smith 
3185edff71fSBarry Smith    Input Parameter:
3195edff71fSBarry Smith +  da - the distributed array
3205edff71fSBarry Smith .  vec - the vector, either a vector the same size as one obtained with
3215edff71fSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
3225edff71fSBarry Smith -  array - the array, non-NULL pointer is zeroed
3235edff71fSBarry Smith 
3245edff71fSBarry Smith   Level: intermediate
3255edff71fSBarry Smith 
32695452b02SPatrick Sanan   Fortran Notes:
32795452b02SPatrick Sanan     From Fortran use DMDAVecRestoreArrayReadF90()
3285edff71fSBarry Smith 
3295edff71fSBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray()
3305edff71fSBarry Smith @*/
3315edff71fSBarry Smith PetscErrorCode  DMDAVecRestoreArrayRead(DM da,Vec vec,void *array)
3325edff71fSBarry Smith {
3335edff71fSBarry Smith   PetscErrorCode ierr;
3345edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
3355edff71fSBarry Smith 
3365edff71fSBarry Smith   PetscFunctionBegin;
337a9a02de4SBarry Smith   PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA);
3385edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
3395edff71fSBarry Smith   PetscValidPointer(array, 3);
3405edff71fSBarry Smith   if (da->defaultSection) {
3415edff71fSBarry Smith     ierr = VecRestoreArrayRead(vec,(const PetscScalar**)array);CHKERRQ(ierr);
3425edff71fSBarry Smith     PetscFunctionReturn(0);
3435edff71fSBarry Smith   }
3445edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
3455edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
3465edff71fSBarry Smith   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
3475edff71fSBarry Smith 
3485edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
3495edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
3505edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
3515edff71fSBarry Smith     gxm = xm;
3525edff71fSBarry Smith     gym = ym;
3535edff71fSBarry Smith     gzm = zm;
3545edff71fSBarry Smith     gxs = xs;
3555edff71fSBarry Smith     gys = ys;
3565edff71fSBarry Smith     gzs = zs;
3575edff71fSBarry 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);
3585edff71fSBarry Smith 
3595edff71fSBarry Smith   if (dim == 1) {
3605edff71fSBarry Smith     ierr = VecRestoreArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
3615edff71fSBarry Smith   } else if (dim == 2) {
3625edff71fSBarry Smith     ierr = VecRestoreArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
3635edff71fSBarry Smith   } else if (dim == 3) {
3645edff71fSBarry Smith     ierr = VecRestoreArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
3655edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
3665edff71fSBarry Smith   PetscFunctionReturn(0);
3675edff71fSBarry Smith }
3685edff71fSBarry Smith 
3695edff71fSBarry Smith /*@C
3705edff71fSBarry Smith    DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
3715edff71fSBarry Smith       the underlying vector and is indexed using the global dimensions.
3725edff71fSBarry Smith 
3735edff71fSBarry Smith    Not Collective
3745edff71fSBarry Smith 
3755edff71fSBarry Smith    Input Parameter:
3765edff71fSBarry Smith +  da - the distributed array
3775edff71fSBarry Smith -  vec - the vector, either a vector the same size as one obtained with
3785edff71fSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
3795edff71fSBarry Smith 
3805edff71fSBarry Smith    Output Parameter:
3815edff71fSBarry Smith .  array - the array
3825edff71fSBarry Smith 
3835edff71fSBarry Smith    Notes:
3845edff71fSBarry Smith     Call DMDAVecRestoreArrayDOFRead() once you have finished accessing the vector entries.
3855edff71fSBarry Smith 
3865edff71fSBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
3875edff71fSBarry Smith 
3884ebdaf59SBarry Smith     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayReadF90() and declare your array with one higher dimension,
3895edff71fSBarry Smith     see src/dm/examples/tutorials/ex11f90.F
3905edff71fSBarry Smith 
3915edff71fSBarry Smith   Level: intermediate
3925edff71fSBarry Smith 
3935edff71fSBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF()
3945edff71fSBarry Smith @*/
3955edff71fSBarry Smith PetscErrorCode  DMDAVecGetArrayDOFRead(DM da,Vec vec,void *array)
3965edff71fSBarry Smith {
3975edff71fSBarry Smith   PetscErrorCode ierr;
3985edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
3995edff71fSBarry Smith 
4005edff71fSBarry Smith   PetscFunctionBegin;
4015edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
4025edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
4035edff71fSBarry Smith   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
4045edff71fSBarry Smith 
4055edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
4065edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
4075edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
4085edff71fSBarry Smith     gxm = xm;
4095edff71fSBarry Smith     gym = ym;
4105edff71fSBarry Smith     gzm = zm;
4115edff71fSBarry Smith     gxs = xs;
4125edff71fSBarry Smith     gys = ys;
4135edff71fSBarry Smith     gzs = zs;
4145edff71fSBarry 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);
4155edff71fSBarry Smith 
4165edff71fSBarry Smith   if (dim == 1) {
4175edff71fSBarry Smith     ierr = VecGetArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
4185edff71fSBarry Smith   } else if (dim == 2) {
4195edff71fSBarry Smith     ierr = VecGetArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
4205edff71fSBarry Smith   } else if (dim == 3) {
4215edff71fSBarry Smith     ierr = VecGetArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
4225edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
4235edff71fSBarry Smith   PetscFunctionReturn(0);
4245edff71fSBarry Smith }
4255edff71fSBarry Smith 
4265edff71fSBarry Smith /*@
4275edff71fSBarry Smith    DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with DMDAVecGetArrayDOFRead()
4285edff71fSBarry Smith 
4295edff71fSBarry Smith    Not Collective
4305edff71fSBarry Smith 
4315edff71fSBarry Smith    Input Parameter:
4325edff71fSBarry Smith +  da - the distributed array
4335edff71fSBarry Smith .  vec - the vector, either a vector the same size as one obtained with
4345edff71fSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
4355edff71fSBarry Smith -  array - the array
4365edff71fSBarry Smith 
4375edff71fSBarry Smith   Level: intermediate
4385edff71fSBarry Smith 
4395edff71fSBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF()
4405edff71fSBarry Smith @*/
4415edff71fSBarry Smith PetscErrorCode  DMDAVecRestoreArrayDOFRead(DM da,Vec vec,void *array)
4425edff71fSBarry Smith {
4435edff71fSBarry Smith   PetscErrorCode ierr;
4445edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
4455edff71fSBarry Smith 
4465edff71fSBarry Smith   PetscFunctionBegin;
4475edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
4485edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
4495edff71fSBarry Smith   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
4505edff71fSBarry Smith 
4515edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
4525edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
4535edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
4545edff71fSBarry Smith     gxm = xm;
4555edff71fSBarry Smith     gym = ym;
4565edff71fSBarry Smith     gzm = zm;
4575edff71fSBarry Smith     gxs = xs;
4585edff71fSBarry Smith     gys = ys;
4595edff71fSBarry Smith     gzs = zs;
4605edff71fSBarry Smith   }
4615edff71fSBarry Smith 
4625edff71fSBarry Smith   if (dim == 1) {
4635edff71fSBarry Smith     ierr = VecRestoreArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
4645edff71fSBarry Smith   } else if (dim == 2) {
4655edff71fSBarry Smith     ierr = VecRestoreArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
4665edff71fSBarry Smith   } else if (dim == 3) {
4675edff71fSBarry Smith     ierr = VecRestoreArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
4685edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
4695edff71fSBarry Smith   PetscFunctionReturn(0);
4705edff71fSBarry Smith }
4715edff71fSBarry Smith 
47247c6ae99SBarry Smith 
47347c6ae99SBarry Smith 
47447c6ae99SBarry Smith 
47547c6ae99SBarry Smith 
47647c6ae99SBarry Smith 
47747c6ae99SBarry Smith 
47847c6ae99SBarry Smith 
47947c6ae99SBarry Smith 
48047c6ae99SBarry Smith 
48147c6ae99SBarry Smith 
48247c6ae99SBarry Smith 
48347c6ae99SBarry Smith 
484