xref: /petsc/src/dm/impls/da/dagetarray.c (revision 2ddbf755899deb97b28dcf235e8c5785f592fce4)
147c6ae99SBarry Smith 
2af0996ceSBarry Smith #include <petsc/private/dmdaimpl.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 
10*2ddbf755SMatthew G. Knepley    Logically collective on Vec
1147c6ae99SBarry Smith 
1247c6ae99SBarry Smith    Input Parameter:
1347c6ae99SBarry Smith +  da - the distributed array
147e57d48aSBarry Smith -  vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
1547c6ae99SBarry Smith 
1647c6ae99SBarry Smith    Output Parameter:
1747c6ae99SBarry Smith .  array - the array
1847c6ae99SBarry Smith 
1947c6ae99SBarry Smith    Notes:
20aa219208SBarry Smith     Call DMDAVecRestoreArray() once you have finished accessing the vector entries.
2147c6ae99SBarry Smith 
2247c6ae99SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
2347c6ae99SBarry Smith 
247e57d48aSBarry Smith     If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
257e57d48aSBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
2647c6ae99SBarry Smith 
277e57d48aSBarry Smith     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
2847c6ae99SBarry Smith 
29aa219208SBarry Smith   Fortran Notes: From Fortran use DMDAVecGetArrayF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
30aa219208SBarry 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
31aa219208SBarry 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
32bef27881SBarry Smith        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
33af0996ceSBarry Smith        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
3447c6ae99SBarry Smith 
35596f5af7SJed Brown   Due to bugs in the compiler DMDAVecGetArrayF90() does not work with gfortran versions before 4.5
3647c6ae99SBarry Smith 
3747c6ae99SBarry Smith   Level: intermediate
3847c6ae99SBarry Smith 
3947c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
4047c6ae99SBarry Smith 
41aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF()
426db82c94SMatthew G Knepley           DMDAVecGetArrayDOF()
4347c6ae99SBarry Smith @*/
447087cfbeSBarry Smith PetscErrorCode  DMDAVecGetArray(DM da,Vec vec,void *array)
4547c6ae99SBarry Smith {
4647c6ae99SBarry Smith   PetscErrorCode ierr;
4747c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
4847c6ae99SBarry Smith 
4947c6ae99SBarry Smith   PetscFunctionBegin;
506db82c94SMatthew G Knepley   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
516db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
526db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
53aa1993deSMatthew G Knepley   if (da->defaultSection) {
54aa1993deSMatthew G Knepley     ierr = VecGetArray(vec,(PetscScalar**)array);CHKERRQ(ierr);
55aa1993deSMatthew G Knepley     PetscFunctionReturn(0);
56aa1993deSMatthew G Knepley   }
57aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
58aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
591321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
6047c6ae99SBarry Smith 
6147c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
6247c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
6347c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
6447c6ae99SBarry Smith     gxm = xm;
6547c6ae99SBarry Smith     gym = ym;
6647c6ae99SBarry Smith     gzm = zm;
6747c6ae99SBarry Smith     gxs = xs;
6847c6ae99SBarry Smith     gys = ys;
6947c6ae99SBarry Smith     gzs = zs;
7030729d88SBarry 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);
7147c6ae99SBarry Smith 
7247c6ae99SBarry Smith   if (dim == 1) {
7347c6ae99SBarry Smith     ierr = VecGetArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
7447c6ae99SBarry Smith   } else if (dim == 2) {
7547c6ae99SBarry Smith     ierr = VecGetArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
7647c6ae99SBarry Smith   } else if (dim == 3) {
7747c6ae99SBarry Smith     ierr = VecGetArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
7830729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
7947c6ae99SBarry Smith   PetscFunctionReturn(0);
8047c6ae99SBarry Smith }
8147c6ae99SBarry Smith 
8247c6ae99SBarry Smith #undef __FUNCT__
83aa219208SBarry Smith #define __FUNCT__ "DMDAVecRestoreArray"
8447c6ae99SBarry Smith /*@
85aa219208SBarry Smith    DMDAVecRestoreArray - Restores a multiple dimension array obtained with DMDAVecGetArray()
8647c6ae99SBarry Smith 
87*2ddbf755SMatthew G. Knepley    Logically collective on Vec
8847c6ae99SBarry Smith 
8947c6ae99SBarry Smith    Input Parameter:
9047c6ae99SBarry Smith +  da - the distributed array
9147c6ae99SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
92564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
93e5c84f05SJed Brown -  array - the array, non-NULL pointer is zeroed
9447c6ae99SBarry Smith 
9547c6ae99SBarry Smith   Level: intermediate
9647c6ae99SBarry Smith 
97aa219208SBarry Smith   Fortran Notes: From Fortran use DMDAVecRestoreArrayF90()
9847c6ae99SBarry Smith 
9947c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
10047c6ae99SBarry Smith 
101aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray()
10247c6ae99SBarry Smith @*/
1037087cfbeSBarry Smith PetscErrorCode  DMDAVecRestoreArray(DM da,Vec vec,void *array)
10447c6ae99SBarry Smith {
10547c6ae99SBarry Smith   PetscErrorCode ierr;
10647c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
10747c6ae99SBarry Smith 
10847c6ae99SBarry Smith   PetscFunctionBegin;
1096db82c94SMatthew G Knepley   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
1106db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
1116db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
112aa1993deSMatthew G Knepley   if (da->defaultSection) {
113aa1993deSMatthew G Knepley     ierr = VecRestoreArray(vec,(PetscScalar**)array);CHKERRQ(ierr);
114aa1993deSMatthew G Knepley     PetscFunctionReturn(0);
115aa1993deSMatthew G Knepley   }
116aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
117aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
1181321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
11947c6ae99SBarry Smith 
12047c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
12147c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
12247c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
12347c6ae99SBarry Smith     gxm = xm;
12447c6ae99SBarry Smith     gym = ym;
12547c6ae99SBarry Smith     gzm = zm;
12647c6ae99SBarry Smith     gxs = xs;
12747c6ae99SBarry Smith     gys = ys;
12847c6ae99SBarry Smith     gzs = zs;
12930729d88SBarry 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);
13047c6ae99SBarry Smith 
13147c6ae99SBarry Smith   if (dim == 1) {
13247c6ae99SBarry Smith     ierr = VecRestoreArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
13347c6ae99SBarry Smith   } else if (dim == 2) {
13447c6ae99SBarry Smith     ierr = VecRestoreArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
13547c6ae99SBarry Smith   } else if (dim == 3) {
13647c6ae99SBarry Smith     ierr = VecRestoreArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
13730729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
13847c6ae99SBarry Smith   PetscFunctionReturn(0);
13947c6ae99SBarry Smith }
14047c6ae99SBarry Smith 
14147c6ae99SBarry Smith #undef __FUNCT__
142aa219208SBarry Smith #define __FUNCT__ "DMDAVecGetArrayDOF"
14347c6ae99SBarry Smith /*@C
144aa219208SBarry Smith    DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
14547c6ae99SBarry Smith       the underlying vector and is indexed using the global dimensions.
14647c6ae99SBarry Smith 
147*2ddbf755SMatthew G. Knepley    Logically collective
14847c6ae99SBarry Smith 
14947c6ae99SBarry Smith    Input Parameter:
15047c6ae99SBarry Smith +  da - the distributed array
15147c6ae99SBarry Smith -  vec - the vector, either a vector the same size as one obtained with
152564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
15347c6ae99SBarry Smith 
15447c6ae99SBarry Smith    Output Parameter:
15547c6ae99SBarry Smith .  array - the array
15647c6ae99SBarry Smith 
15747c6ae99SBarry Smith    Notes:
158aa219208SBarry Smith     Call DMDAVecRestoreArrayDOF() once you have finished accessing the vector entries.
15947c6ae99SBarry Smith 
16047c6ae99SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
16147c6ae99SBarry Smith 
1621b82215eSBarry Smith     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayF90() and declare your array with one higher dimension,
1631b82215eSBarry Smith     see src/dm/examples/tutorials/ex11f90.F
1641b82215eSBarry Smith 
16547c6ae99SBarry Smith   Level: intermediate
16647c6ae99SBarry Smith 
16747c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
16847c6ae99SBarry Smith 
169aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF()
17047c6ae99SBarry Smith @*/
1717087cfbeSBarry Smith PetscErrorCode  DMDAVecGetArrayDOF(DM da,Vec vec,void *array)
17247c6ae99SBarry Smith {
17347c6ae99SBarry Smith   PetscErrorCode ierr;
17447c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
17547c6ae99SBarry Smith 
17647c6ae99SBarry Smith   PetscFunctionBegin;
177aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
178aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
1791321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
18047c6ae99SBarry Smith 
18147c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
18247c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
18347c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
18447c6ae99SBarry Smith     gxm = xm;
18547c6ae99SBarry Smith     gym = ym;
18647c6ae99SBarry Smith     gzm = zm;
18747c6ae99SBarry Smith     gxs = xs;
18847c6ae99SBarry Smith     gys = ys;
18947c6ae99SBarry Smith     gzs = zs;
19030729d88SBarry 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);
19147c6ae99SBarry Smith 
19247c6ae99SBarry Smith   if (dim == 1) {
19347c6ae99SBarry Smith     ierr = VecGetArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
19447c6ae99SBarry Smith   } else if (dim == 2) {
19547c6ae99SBarry Smith     ierr = VecGetArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
19647c6ae99SBarry Smith   } else if (dim == 3) {
19747c6ae99SBarry Smith     ierr = VecGetArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
19830729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
19947c6ae99SBarry Smith   PetscFunctionReturn(0);
20047c6ae99SBarry Smith }
20147c6ae99SBarry Smith 
20247c6ae99SBarry Smith #undef __FUNCT__
203aa219208SBarry Smith #define __FUNCT__ "DMDAVecRestoreArrayDOF"
20447c6ae99SBarry Smith /*@
205aa219208SBarry Smith    DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with DMDAVecGetArrayDOF()
20647c6ae99SBarry Smith 
207*2ddbf755SMatthew G. Knepley    Logically collective
20847c6ae99SBarry Smith 
20947c6ae99SBarry Smith    Input Parameter:
21047c6ae99SBarry Smith +  da - the distributed array
21147c6ae99SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
212564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
21347c6ae99SBarry Smith -  array - the array
21447c6ae99SBarry Smith 
21547c6ae99SBarry Smith   Level: intermediate
21647c6ae99SBarry Smith 
21747c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
21847c6ae99SBarry Smith 
219aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF()
22047c6ae99SBarry Smith @*/
2217087cfbeSBarry Smith PetscErrorCode  DMDAVecRestoreArrayDOF(DM da,Vec vec,void *array)
22247c6ae99SBarry Smith {
22347c6ae99SBarry Smith   PetscErrorCode ierr;
22447c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
22547c6ae99SBarry Smith 
22647c6ae99SBarry Smith   PetscFunctionBegin;
227aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
228aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
2291321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
23047c6ae99SBarry Smith 
23147c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
23247c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
23347c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
23447c6ae99SBarry Smith     gxm = xm;
23547c6ae99SBarry Smith     gym = ym;
23647c6ae99SBarry Smith     gzm = zm;
23747c6ae99SBarry Smith     gxs = xs;
23847c6ae99SBarry Smith     gys = ys;
23947c6ae99SBarry Smith     gzs = zs;
24047c6ae99SBarry Smith   }
24147c6ae99SBarry Smith 
24247c6ae99SBarry Smith   if (dim == 1) {
24347c6ae99SBarry Smith     ierr = VecRestoreArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
24447c6ae99SBarry Smith   } else if (dim == 2) {
24547c6ae99SBarry Smith     ierr = VecRestoreArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
24647c6ae99SBarry Smith   } else if (dim == 3) {
24747c6ae99SBarry Smith     ierr = VecRestoreArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
24830729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
24947c6ae99SBarry Smith   PetscFunctionReturn(0);
25047c6ae99SBarry Smith }
25147c6ae99SBarry Smith 
2525edff71fSBarry Smith #undef __FUNCT__
2535edff71fSBarry Smith #define __FUNCT__ "DMDAVecGetArrayRead"
2545edff71fSBarry Smith /*@C
2555edff71fSBarry Smith    DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
2565edff71fSBarry Smith       the underlying vector and is indexed using the global dimensions.
2575edff71fSBarry Smith 
258*2ddbf755SMatthew G. Knepley    Not collective
2595edff71fSBarry Smith 
2605edff71fSBarry Smith    Input Parameter:
2615edff71fSBarry Smith +  da - the distributed array
2625edff71fSBarry Smith -  vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
2635edff71fSBarry Smith 
2645edff71fSBarry Smith    Output Parameter:
2655edff71fSBarry Smith .  array - the array
2665edff71fSBarry Smith 
2675edff71fSBarry Smith    Notes:
2685edff71fSBarry Smith     Call DMDAVecRestoreArrayRead() once you have finished accessing the vector entries.
2695edff71fSBarry Smith 
2705edff71fSBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
2715edff71fSBarry Smith 
2725edff71fSBarry Smith     If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
2735edff71fSBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
2745edff71fSBarry Smith 
2755edff71fSBarry Smith     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
2765edff71fSBarry Smith 
2775edff71fSBarry Smith   Fortran Notes: From Fortran use DMDAVecGetArrayF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
2785edff71fSBarry 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
2795edff71fSBarry 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
2805edff71fSBarry Smith        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
281af0996ceSBarry Smith        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
2825edff71fSBarry Smith 
2835edff71fSBarry Smith   Due to bugs in the compiler DMDAVecGetArrayF90() does not work with gfortran versions before 4.5
2845edff71fSBarry Smith 
2855edff71fSBarry Smith   Level: intermediate
2865edff71fSBarry Smith 
2875edff71fSBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
2885edff71fSBarry Smith 
2895edff71fSBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF()
2905edff71fSBarry Smith           DMDAVecGetArrayDOF()
2915edff71fSBarry Smith @*/
2925edff71fSBarry Smith PetscErrorCode  DMDAVecGetArrayRead(DM da,Vec vec,void *array)
2935edff71fSBarry Smith {
2945edff71fSBarry Smith   PetscErrorCode ierr;
2955edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
2965edff71fSBarry Smith 
2975edff71fSBarry Smith   PetscFunctionBegin;
2985edff71fSBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
2995edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
3005edff71fSBarry Smith   PetscValidPointer(array, 3);
3015edff71fSBarry Smith   if (da->defaultSection) {
3025edff71fSBarry Smith     ierr = VecGetArrayRead(vec,(const PetscScalar**)array);CHKERRQ(ierr);
3035edff71fSBarry Smith     PetscFunctionReturn(0);
3045edff71fSBarry Smith   }
3055edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
3065edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
3075edff71fSBarry Smith   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
3085edff71fSBarry Smith 
3095edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
3105edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
3115edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
3125edff71fSBarry Smith     gxm = xm;
3135edff71fSBarry Smith     gym = ym;
3145edff71fSBarry Smith     gzm = zm;
3155edff71fSBarry Smith     gxs = xs;
3165edff71fSBarry Smith     gys = ys;
3175edff71fSBarry Smith     gzs = zs;
3185edff71fSBarry 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);
3195edff71fSBarry Smith 
3205edff71fSBarry Smith   if (dim == 1) {
3215edff71fSBarry Smith     ierr = VecGetArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
3225edff71fSBarry Smith   } else if (dim == 2) {
3235edff71fSBarry Smith     ierr = VecGetArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
3245edff71fSBarry Smith   } else if (dim == 3) {
3255edff71fSBarry Smith     ierr = VecGetArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
3265edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
3275edff71fSBarry Smith   PetscFunctionReturn(0);
3285edff71fSBarry Smith }
3295edff71fSBarry Smith 
3305edff71fSBarry Smith #undef __FUNCT__
3315edff71fSBarry Smith #define __FUNCT__ "DMDAVecRestoreArrayRead"
3325edff71fSBarry Smith /*@
3335edff71fSBarry Smith    DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with DMDAVecGetArrayRead()
3345edff71fSBarry Smith 
335*2ddbf755SMatthew G. Knepley    Not collective
3365edff71fSBarry Smith 
3375edff71fSBarry Smith    Input Parameter:
3385edff71fSBarry Smith +  da - the distributed array
3395edff71fSBarry Smith .  vec - the vector, either a vector the same size as one obtained with
3405edff71fSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
3415edff71fSBarry Smith -  array - the array, non-NULL pointer is zeroed
3425edff71fSBarry Smith 
3435edff71fSBarry Smith   Level: intermediate
3445edff71fSBarry Smith 
3455edff71fSBarry Smith   Fortran Notes: From Fortran use DMDAVecRestoreArrayF90()
3465edff71fSBarry Smith 
3475edff71fSBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
3485edff71fSBarry Smith 
3495edff71fSBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray()
3505edff71fSBarry Smith @*/
3515edff71fSBarry Smith PetscErrorCode  DMDAVecRestoreArrayRead(DM da,Vec vec,void *array)
3525edff71fSBarry Smith {
3535edff71fSBarry Smith   PetscErrorCode ierr;
3545edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
3555edff71fSBarry Smith 
3565edff71fSBarry Smith   PetscFunctionBegin;
3575edff71fSBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
3585edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
3595edff71fSBarry Smith   PetscValidPointer(array, 3);
3605edff71fSBarry Smith   if (da->defaultSection) {
3615edff71fSBarry Smith     ierr = VecRestoreArrayRead(vec,(const PetscScalar**)array);CHKERRQ(ierr);
3625edff71fSBarry Smith     PetscFunctionReturn(0);
3635edff71fSBarry Smith   }
3645edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
3655edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
3665edff71fSBarry Smith   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
3675edff71fSBarry Smith 
3685edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
3695edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
3705edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
3715edff71fSBarry Smith     gxm = xm;
3725edff71fSBarry Smith     gym = ym;
3735edff71fSBarry Smith     gzm = zm;
3745edff71fSBarry Smith     gxs = xs;
3755edff71fSBarry Smith     gys = ys;
3765edff71fSBarry Smith     gzs = zs;
3775edff71fSBarry 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);
3785edff71fSBarry Smith 
3795edff71fSBarry Smith   if (dim == 1) {
3805edff71fSBarry Smith     ierr = VecRestoreArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
3815edff71fSBarry Smith   } else if (dim == 2) {
3825edff71fSBarry Smith     ierr = VecRestoreArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
3835edff71fSBarry Smith   } else if (dim == 3) {
3845edff71fSBarry Smith     ierr = VecRestoreArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
3855edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
3865edff71fSBarry Smith   PetscFunctionReturn(0);
3875edff71fSBarry Smith }
3885edff71fSBarry Smith 
3895edff71fSBarry Smith #undef __FUNCT__
3905edff71fSBarry Smith #define __FUNCT__ "DMDAVecGetArrayDOFRead"
3915edff71fSBarry Smith /*@C
3925edff71fSBarry Smith    DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
3935edff71fSBarry Smith       the underlying vector and is indexed using the global dimensions.
3945edff71fSBarry Smith 
3955edff71fSBarry Smith    Not Collective
3965edff71fSBarry Smith 
3975edff71fSBarry Smith    Input Parameter:
3985edff71fSBarry Smith +  da - the distributed array
3995edff71fSBarry Smith -  vec - the vector, either a vector the same size as one obtained with
4005edff71fSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
4015edff71fSBarry Smith 
4025edff71fSBarry Smith    Output Parameter:
4035edff71fSBarry Smith .  array - the array
4045edff71fSBarry Smith 
4055edff71fSBarry Smith    Notes:
4065edff71fSBarry Smith     Call DMDAVecRestoreArrayDOFRead() once you have finished accessing the vector entries.
4075edff71fSBarry Smith 
4085edff71fSBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
4095edff71fSBarry Smith 
4105edff71fSBarry Smith     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayF90() and declare your array with one higher dimension,
4115edff71fSBarry Smith     see src/dm/examples/tutorials/ex11f90.F
4125edff71fSBarry Smith 
4135edff71fSBarry Smith   Level: intermediate
4145edff71fSBarry Smith 
4155edff71fSBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
4165edff71fSBarry Smith 
4175edff71fSBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF()
4185edff71fSBarry Smith @*/
4195edff71fSBarry Smith PetscErrorCode  DMDAVecGetArrayDOFRead(DM da,Vec vec,void *array)
4205edff71fSBarry Smith {
4215edff71fSBarry Smith   PetscErrorCode ierr;
4225edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
4235edff71fSBarry Smith 
4245edff71fSBarry Smith   PetscFunctionBegin;
4255edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
4265edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
4275edff71fSBarry Smith   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
4285edff71fSBarry Smith 
4295edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
4305edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
4315edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
4325edff71fSBarry Smith     gxm = xm;
4335edff71fSBarry Smith     gym = ym;
4345edff71fSBarry Smith     gzm = zm;
4355edff71fSBarry Smith     gxs = xs;
4365edff71fSBarry Smith     gys = ys;
4375edff71fSBarry Smith     gzs = zs;
4385edff71fSBarry 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);
4395edff71fSBarry Smith 
4405edff71fSBarry Smith   if (dim == 1) {
4415edff71fSBarry Smith     ierr = VecGetArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
4425edff71fSBarry Smith   } else if (dim == 2) {
4435edff71fSBarry Smith     ierr = VecGetArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
4445edff71fSBarry Smith   } else if (dim == 3) {
4455edff71fSBarry Smith     ierr = VecGetArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
4465edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
4475edff71fSBarry Smith   PetscFunctionReturn(0);
4485edff71fSBarry Smith }
4495edff71fSBarry Smith 
4505edff71fSBarry Smith #undef __FUNCT__
4515edff71fSBarry Smith #define __FUNCT__ "DMDAVecRestoreArrayDOFRead"
4525edff71fSBarry Smith /*@
4535edff71fSBarry Smith    DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with DMDAVecGetArrayDOFRead()
4545edff71fSBarry Smith 
4555edff71fSBarry Smith    Not Collective
4565edff71fSBarry Smith 
4575edff71fSBarry Smith    Input Parameter:
4585edff71fSBarry Smith +  da - the distributed array
4595edff71fSBarry Smith .  vec - the vector, either a vector the same size as one obtained with
4605edff71fSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
4615edff71fSBarry Smith -  array - the array
4625edff71fSBarry Smith 
4635edff71fSBarry Smith   Level: intermediate
4645edff71fSBarry Smith 
4655edff71fSBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
4665edff71fSBarry Smith 
4675edff71fSBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF()
4685edff71fSBarry Smith @*/
4695edff71fSBarry Smith PetscErrorCode  DMDAVecRestoreArrayDOFRead(DM da,Vec vec,void *array)
4705edff71fSBarry Smith {
4715edff71fSBarry Smith   PetscErrorCode ierr;
4725edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
4735edff71fSBarry Smith 
4745edff71fSBarry Smith   PetscFunctionBegin;
4755edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
4765edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
4775edff71fSBarry Smith   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
4785edff71fSBarry Smith 
4795edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
4805edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
4815edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
4825edff71fSBarry Smith     gxm = xm;
4835edff71fSBarry Smith     gym = ym;
4845edff71fSBarry Smith     gzm = zm;
4855edff71fSBarry Smith     gxs = xs;
4865edff71fSBarry Smith     gys = ys;
4875edff71fSBarry Smith     gzs = zs;
4885edff71fSBarry Smith   }
4895edff71fSBarry Smith 
4905edff71fSBarry Smith   if (dim == 1) {
4915edff71fSBarry Smith     ierr = VecRestoreArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
4925edff71fSBarry Smith   } else if (dim == 2) {
4935edff71fSBarry Smith     ierr = VecRestoreArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
4945edff71fSBarry Smith   } else if (dim == 3) {
4955edff71fSBarry Smith     ierr = VecRestoreArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
4965edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
4975edff71fSBarry Smith   PetscFunctionReturn(0);
4985edff71fSBarry Smith }
4995edff71fSBarry Smith 
50047c6ae99SBarry Smith 
50147c6ae99SBarry Smith 
50247c6ae99SBarry Smith 
50347c6ae99SBarry Smith 
50447c6ae99SBarry Smith 
50547c6ae99SBarry Smith 
50647c6ae99SBarry Smith 
50747c6ae99SBarry Smith 
50847c6ae99SBarry Smith 
50947c6ae99SBarry Smith 
51047c6ae99SBarry Smith 
51147c6ae99SBarry Smith 
512