xref: /petsc/src/dm/impls/da/dagetarray.c (revision 95452b02e12c0ee11232c7ff2b24b568a8e07e43)
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 
82ddbf755SMatthew G. Knepley    Logically collective on Vec
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 
27*95452b02SPatrick Sanan   Fortran Notes:
28*95452b02SPatrick 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 
3847c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
3947c6ae99SBarry Smith 
40aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF()
416db82c94SMatthew G Knepley           DMDAVecGetArrayDOF()
4247c6ae99SBarry Smith @*/
437087cfbeSBarry Smith PetscErrorCode  DMDAVecGetArray(DM da,Vec vec,void *array)
4447c6ae99SBarry Smith {
4547c6ae99SBarry Smith   PetscErrorCode ierr;
4647c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
4747c6ae99SBarry Smith 
4847c6ae99SBarry Smith   PetscFunctionBegin;
496db82c94SMatthew G Knepley   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
506db82c94SMatthew G Knepley   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
516db82c94SMatthew G Knepley   PetscValidPointer(array, 3);
52aa1993deSMatthew G Knepley   if (da->defaultSection) {
53aa1993deSMatthew G Knepley     ierr = VecGetArray(vec,(PetscScalar**)array);CHKERRQ(ierr);
54aa1993deSMatthew G Knepley     PetscFunctionReturn(0);
55aa1993deSMatthew G Knepley   }
56aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
57aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
581321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
5947c6ae99SBarry Smith 
6047c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
6147c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
6247c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
6347c6ae99SBarry Smith     gxm = xm;
6447c6ae99SBarry Smith     gym = ym;
6547c6ae99SBarry Smith     gzm = zm;
6647c6ae99SBarry Smith     gxs = xs;
6747c6ae99SBarry Smith     gys = ys;
6847c6ae99SBarry Smith     gzs = zs;
6930729d88SBarry 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);
7047c6ae99SBarry Smith 
7147c6ae99SBarry Smith   if (dim == 1) {
7247c6ae99SBarry Smith     ierr = VecGetArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
7347c6ae99SBarry Smith   } else if (dim == 2) {
7447c6ae99SBarry Smith     ierr = VecGetArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
7547c6ae99SBarry Smith   } else if (dim == 3) {
7647c6ae99SBarry Smith     ierr = VecGetArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
7730729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
7847c6ae99SBarry Smith   PetscFunctionReturn(0);
7947c6ae99SBarry Smith }
8047c6ae99SBarry Smith 
8147c6ae99SBarry Smith /*@
82aa219208SBarry Smith    DMDAVecRestoreArray - Restores a multiple dimension array obtained with DMDAVecGetArray()
8347c6ae99SBarry Smith 
842ddbf755SMatthew G. Knepley    Logically collective on Vec
8547c6ae99SBarry Smith 
8647c6ae99SBarry Smith    Input Parameter:
8747c6ae99SBarry Smith +  da - the distributed array
8847c6ae99SBarry Smith .  vec - the vector, either a vector the same size as one obtained with
89564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
90e5c84f05SJed Brown -  array - the array, non-NULL pointer is zeroed
9147c6ae99SBarry Smith 
9247c6ae99SBarry Smith   Level: intermediate
9347c6ae99SBarry Smith 
94*95452b02SPatrick Sanan   Fortran Notes:
95*95452b02SPatrick Sanan     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);
110aa1993deSMatthew G Knepley   if (da->defaultSection) {
111aa1993deSMatthew G Knepley     ierr = VecRestoreArray(vec,(PetscScalar**)array);CHKERRQ(ierr);
112aa1993deSMatthew G Knepley     PetscFunctionReturn(0);
113aa1993deSMatthew G Knepley   }
114aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
115aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
1161321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
11747c6ae99SBarry Smith 
11847c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
11947c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
12047c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
12147c6ae99SBarry Smith     gxm = xm;
12247c6ae99SBarry Smith     gym = ym;
12347c6ae99SBarry Smith     gzm = zm;
12447c6ae99SBarry Smith     gxs = xs;
12547c6ae99SBarry Smith     gys = ys;
12647c6ae99SBarry Smith     gzs = zs;
12730729d88SBarry 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);
12847c6ae99SBarry Smith 
12947c6ae99SBarry Smith   if (dim == 1) {
13047c6ae99SBarry Smith     ierr = VecRestoreArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
13147c6ae99SBarry Smith   } else if (dim == 2) {
13247c6ae99SBarry Smith     ierr = VecRestoreArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
13347c6ae99SBarry Smith   } else if (dim == 3) {
13447c6ae99SBarry Smith     ierr = VecRestoreArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
13530729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
13647c6ae99SBarry Smith   PetscFunctionReturn(0);
13747c6ae99SBarry Smith }
13847c6ae99SBarry Smith 
13947c6ae99SBarry Smith /*@C
140aa219208SBarry Smith    DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
14147c6ae99SBarry Smith       the underlying vector and is indexed using the global dimensions.
14247c6ae99SBarry Smith 
1432ddbf755SMatthew G. Knepley    Logically collective
14447c6ae99SBarry Smith 
14547c6ae99SBarry Smith    Input Parameter:
14647c6ae99SBarry Smith +  da - the distributed array
14747c6ae99SBarry Smith -  vec - the vector, either a vector the same size as one obtained with
148564755cdSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
14947c6ae99SBarry Smith 
15047c6ae99SBarry Smith    Output Parameter:
15147c6ae99SBarry Smith .  array - the array
15247c6ae99SBarry Smith 
15347c6ae99SBarry Smith    Notes:
154aa219208SBarry Smith     Call DMDAVecRestoreArrayDOF() once you have finished accessing the vector entries.
15547c6ae99SBarry Smith 
15647c6ae99SBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
15747c6ae99SBarry Smith 
1581b82215eSBarry Smith     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayF90() and declare your array with one higher dimension,
1591b82215eSBarry Smith     see src/dm/examples/tutorials/ex11f90.F
1601b82215eSBarry Smith 
16147c6ae99SBarry Smith   Level: intermediate
16247c6ae99SBarry Smith 
16347c6ae99SBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
16447c6ae99SBarry Smith 
165aa219208SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF()
16647c6ae99SBarry Smith @*/
1677087cfbeSBarry Smith PetscErrorCode  DMDAVecGetArrayDOF(DM da,Vec vec,void *array)
16847c6ae99SBarry Smith {
16947c6ae99SBarry Smith   PetscErrorCode ierr;
17047c6ae99SBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
17147c6ae99SBarry Smith 
17247c6ae99SBarry Smith   PetscFunctionBegin;
173aa219208SBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
174aa219208SBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
1751321219cSEthan Coon   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
17647c6ae99SBarry Smith 
17747c6ae99SBarry Smith   /* Handle case where user passes in global vector as opposed to local */
17847c6ae99SBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
17947c6ae99SBarry Smith   if (N == xm*ym*zm*dof) {
18047c6ae99SBarry Smith     gxm = xm;
18147c6ae99SBarry Smith     gym = ym;
18247c6ae99SBarry Smith     gzm = zm;
18347c6ae99SBarry Smith     gxs = xs;
18447c6ae99SBarry Smith     gys = ys;
18547c6ae99SBarry Smith     gzs = zs;
18630729d88SBarry 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);
18747c6ae99SBarry Smith 
18847c6ae99SBarry Smith   if (dim == 1) {
18947c6ae99SBarry Smith     ierr = VecGetArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
19047c6ae99SBarry Smith   } else if (dim == 2) {
19147c6ae99SBarry Smith     ierr = VecGetArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
19247c6ae99SBarry Smith   } else if (dim == 3) {
19347c6ae99SBarry Smith     ierr = VecGetArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
19430729d88SBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
19547c6ae99SBarry Smith   PetscFunctionReturn(0);
19647c6ae99SBarry Smith }
19747c6ae99SBarry Smith 
19847c6ae99SBarry Smith /*@
199aa219208SBarry Smith    DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with DMDAVecGetArrayDOF()
20047c6ae99SBarry Smith 
2012ddbf755SMatthew G. Knepley    Logically 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 
2465edff71fSBarry Smith /*@C
2475edff71fSBarry Smith    DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
2485edff71fSBarry Smith       the underlying vector and is indexed using the global dimensions.
2495edff71fSBarry Smith 
2502ddbf755SMatthew G. Knepley    Not collective
2515edff71fSBarry Smith 
2525edff71fSBarry Smith    Input Parameter:
2535edff71fSBarry Smith +  da - the distributed array
2545edff71fSBarry Smith -  vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector()
2555edff71fSBarry Smith 
2565edff71fSBarry Smith    Output Parameter:
2575edff71fSBarry Smith .  array - the array
2585edff71fSBarry Smith 
2595edff71fSBarry Smith    Notes:
2605edff71fSBarry Smith     Call DMDAVecRestoreArrayRead() once you have finished accessing the vector entries.
2615edff71fSBarry Smith 
2625edff71fSBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
2635edff71fSBarry Smith 
2645edff71fSBarry Smith     If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
2655edff71fSBarry Smith     a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
2665edff71fSBarry Smith 
2675edff71fSBarry Smith     appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations.
2685edff71fSBarry Smith 
269*95452b02SPatrick Sanan   Fortran Notes:
270*95452b02SPatrick Sanan     From Fortran use DMDAVecGetArrayReadF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
2715edff71fSBarry 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
2725edff71fSBarry 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
2735edff71fSBarry Smith        array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from
274af0996ceSBarry Smith        DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine.
2755edff71fSBarry Smith 
2764ebdaf59SBarry Smith   Due to bugs in the compiler DMDAVecGetArrayReadF90() does not work with gfortran versions before 4.5
2775edff71fSBarry Smith 
2785edff71fSBarry Smith   Level: intermediate
2795edff71fSBarry Smith 
2805edff71fSBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
2815edff71fSBarry Smith 
2825edff71fSBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF()
2835edff71fSBarry Smith           DMDAVecGetArrayDOF()
2845edff71fSBarry Smith @*/
2855edff71fSBarry Smith PetscErrorCode  DMDAVecGetArrayRead(DM da,Vec vec,void *array)
2865edff71fSBarry Smith {
2875edff71fSBarry Smith   PetscErrorCode ierr;
2885edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
2895edff71fSBarry Smith 
2905edff71fSBarry Smith   PetscFunctionBegin;
2915edff71fSBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
2925edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
2935edff71fSBarry Smith   PetscValidPointer(array, 3);
2945edff71fSBarry Smith   if (da->defaultSection) {
2955edff71fSBarry Smith     ierr = VecGetArrayRead(vec,(const PetscScalar**)array);CHKERRQ(ierr);
2965edff71fSBarry Smith     PetscFunctionReturn(0);
2975edff71fSBarry Smith   }
2985edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
2995edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
3005edff71fSBarry Smith   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
3015edff71fSBarry Smith 
3025edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
3035edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
3045edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
3055edff71fSBarry Smith     gxm = xm;
3065edff71fSBarry Smith     gym = ym;
3075edff71fSBarry Smith     gzm = zm;
3085edff71fSBarry Smith     gxs = xs;
3095edff71fSBarry Smith     gys = ys;
3105edff71fSBarry Smith     gzs = zs;
3115edff71fSBarry 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);
3125edff71fSBarry Smith 
3135edff71fSBarry Smith   if (dim == 1) {
3145edff71fSBarry Smith     ierr = VecGetArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
3155edff71fSBarry Smith   } else if (dim == 2) {
3165edff71fSBarry Smith     ierr = VecGetArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
3175edff71fSBarry Smith   } else if (dim == 3) {
3185edff71fSBarry Smith     ierr = VecGetArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
3195edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
3205edff71fSBarry Smith   PetscFunctionReturn(0);
3215edff71fSBarry Smith }
3225edff71fSBarry Smith 
3235edff71fSBarry Smith /*@
3245edff71fSBarry Smith    DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with DMDAVecGetArrayRead()
3255edff71fSBarry Smith 
3262ddbf755SMatthew G. Knepley    Not collective
3275edff71fSBarry Smith 
3285edff71fSBarry Smith    Input Parameter:
3295edff71fSBarry Smith +  da - the distributed array
3305edff71fSBarry Smith .  vec - the vector, either a vector the same size as one obtained with
3315edff71fSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
3325edff71fSBarry Smith -  array - the array, non-NULL pointer is zeroed
3335edff71fSBarry Smith 
3345edff71fSBarry Smith   Level: intermediate
3355edff71fSBarry Smith 
336*95452b02SPatrick Sanan   Fortran Notes:
337*95452b02SPatrick Sanan     From Fortran use DMDAVecRestoreArrayReadF90()
3385edff71fSBarry Smith 
3395edff71fSBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
3405edff71fSBarry Smith 
3415edff71fSBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray()
3425edff71fSBarry Smith @*/
3435edff71fSBarry Smith PetscErrorCode  DMDAVecRestoreArrayRead(DM da,Vec vec,void *array)
3445edff71fSBarry Smith {
3455edff71fSBarry Smith   PetscErrorCode ierr;
3465edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
3475edff71fSBarry Smith 
3485edff71fSBarry Smith   PetscFunctionBegin;
3495edff71fSBarry Smith   PetscValidHeaderSpecific(da, DM_CLASSID, 1);
3505edff71fSBarry Smith   PetscValidHeaderSpecific(vec, VEC_CLASSID, 2);
3515edff71fSBarry Smith   PetscValidPointer(array, 3);
3525edff71fSBarry Smith   if (da->defaultSection) {
3535edff71fSBarry Smith     ierr = VecRestoreArrayRead(vec,(const PetscScalar**)array);CHKERRQ(ierr);
3545edff71fSBarry Smith     PetscFunctionReturn(0);
3555edff71fSBarry Smith   }
3565edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
3575edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
3585edff71fSBarry Smith   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
3595edff71fSBarry Smith 
3605edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
3615edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
3625edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
3635edff71fSBarry Smith     gxm = xm;
3645edff71fSBarry Smith     gym = ym;
3655edff71fSBarry Smith     gzm = zm;
3665edff71fSBarry Smith     gxs = xs;
3675edff71fSBarry Smith     gys = ys;
3685edff71fSBarry Smith     gzs = zs;
3695edff71fSBarry 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);
3705edff71fSBarry Smith 
3715edff71fSBarry Smith   if (dim == 1) {
3725edff71fSBarry Smith     ierr = VecRestoreArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr);
3735edff71fSBarry Smith   } else if (dim == 2) {
3745edff71fSBarry Smith     ierr = VecRestoreArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr);
3755edff71fSBarry Smith   } else if (dim == 3) {
3765edff71fSBarry Smith     ierr = VecRestoreArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr);
3775edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
3785edff71fSBarry Smith   PetscFunctionReturn(0);
3795edff71fSBarry Smith }
3805edff71fSBarry Smith 
3815edff71fSBarry Smith /*@C
3825edff71fSBarry Smith    DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
3835edff71fSBarry Smith       the underlying vector and is indexed using the global dimensions.
3845edff71fSBarry Smith 
3855edff71fSBarry Smith    Not Collective
3865edff71fSBarry Smith 
3875edff71fSBarry Smith    Input Parameter:
3885edff71fSBarry Smith +  da - the distributed array
3895edff71fSBarry Smith -  vec - the vector, either a vector the same size as one obtained with
3905edff71fSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
3915edff71fSBarry Smith 
3925edff71fSBarry Smith    Output Parameter:
3935edff71fSBarry Smith .  array - the array
3945edff71fSBarry Smith 
3955edff71fSBarry Smith    Notes:
3965edff71fSBarry Smith     Call DMDAVecRestoreArrayDOFRead() once you have finished accessing the vector entries.
3975edff71fSBarry Smith 
3985edff71fSBarry Smith     In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
3995edff71fSBarry Smith 
4004ebdaf59SBarry Smith     In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use  DMDAVecRestoreArrayReadF90() and declare your array with one higher dimension,
4015edff71fSBarry Smith     see src/dm/examples/tutorials/ex11f90.F
4025edff71fSBarry Smith 
4035edff71fSBarry Smith   Level: intermediate
4045edff71fSBarry Smith 
4055edff71fSBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
4065edff71fSBarry Smith 
4075edff71fSBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF()
4085edff71fSBarry Smith @*/
4095edff71fSBarry Smith PetscErrorCode  DMDAVecGetArrayDOFRead(DM da,Vec vec,void *array)
4105edff71fSBarry Smith {
4115edff71fSBarry Smith   PetscErrorCode ierr;
4125edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
4135edff71fSBarry Smith 
4145edff71fSBarry Smith   PetscFunctionBegin;
4155edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
4165edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
4175edff71fSBarry Smith   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
4185edff71fSBarry Smith 
4195edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
4205edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
4215edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
4225edff71fSBarry Smith     gxm = xm;
4235edff71fSBarry Smith     gym = ym;
4245edff71fSBarry Smith     gzm = zm;
4255edff71fSBarry Smith     gxs = xs;
4265edff71fSBarry Smith     gys = ys;
4275edff71fSBarry Smith     gzs = zs;
4285edff71fSBarry 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);
4295edff71fSBarry Smith 
4305edff71fSBarry Smith   if (dim == 1) {
4315edff71fSBarry Smith     ierr = VecGetArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
4325edff71fSBarry Smith   } else if (dim == 2) {
4335edff71fSBarry Smith     ierr = VecGetArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
4345edff71fSBarry Smith   } else if (dim == 3) {
4355edff71fSBarry Smith     ierr = VecGetArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
4365edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
4375edff71fSBarry Smith   PetscFunctionReturn(0);
4385edff71fSBarry Smith }
4395edff71fSBarry Smith 
4405edff71fSBarry Smith /*@
4415edff71fSBarry Smith    DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with DMDAVecGetArrayDOFRead()
4425edff71fSBarry Smith 
4435edff71fSBarry Smith    Not Collective
4445edff71fSBarry Smith 
4455edff71fSBarry Smith    Input Parameter:
4465edff71fSBarry Smith +  da - the distributed array
4475edff71fSBarry Smith .  vec - the vector, either a vector the same size as one obtained with
4485edff71fSBarry Smith          DMCreateGlobalVector() or DMCreateLocalVector()
4495edff71fSBarry Smith -  array - the array
4505edff71fSBarry Smith 
4515edff71fSBarry Smith   Level: intermediate
4525edff71fSBarry Smith 
4535edff71fSBarry Smith .keywords: distributed array, get, corners, nodes, local indices, coordinates
4545edff71fSBarry Smith 
4555edff71fSBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF()
4565edff71fSBarry Smith @*/
4575edff71fSBarry Smith PetscErrorCode  DMDAVecRestoreArrayDOFRead(DM da,Vec vec,void *array)
4585edff71fSBarry Smith {
4595edff71fSBarry Smith   PetscErrorCode ierr;
4605edff71fSBarry Smith   PetscInt       xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof;
4615edff71fSBarry Smith 
4625edff71fSBarry Smith   PetscFunctionBegin;
4635edff71fSBarry Smith   ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr);
4645edff71fSBarry Smith   ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr);
4655edff71fSBarry Smith   ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
4665edff71fSBarry Smith 
4675edff71fSBarry Smith   /* Handle case where user passes in global vector as opposed to local */
4685edff71fSBarry Smith   ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr);
4695edff71fSBarry Smith   if (N == xm*ym*zm*dof) {
4705edff71fSBarry Smith     gxm = xm;
4715edff71fSBarry Smith     gym = ym;
4725edff71fSBarry Smith     gzm = zm;
4735edff71fSBarry Smith     gxs = xs;
4745edff71fSBarry Smith     gys = ys;
4755edff71fSBarry Smith     gzs = zs;
4765edff71fSBarry Smith   }
4775edff71fSBarry Smith 
4785edff71fSBarry Smith   if (dim == 1) {
4795edff71fSBarry Smith     ierr = VecRestoreArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr);
4805edff71fSBarry Smith   } else if (dim == 2) {
4815edff71fSBarry Smith     ierr = VecRestoreArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr);
4825edff71fSBarry Smith   } else if (dim == 3) {
4835edff71fSBarry Smith     ierr = VecRestoreArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr);
4845edff71fSBarry Smith   } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim);
4855edff71fSBarry Smith   PetscFunctionReturn(0);
4865edff71fSBarry Smith }
4875edff71fSBarry Smith 
48847c6ae99SBarry Smith 
48947c6ae99SBarry Smith 
49047c6ae99SBarry Smith 
49147c6ae99SBarry Smith 
49247c6ae99SBarry Smith 
49347c6ae99SBarry Smith 
49447c6ae99SBarry Smith 
49547c6ae99SBarry Smith 
49647c6ae99SBarry Smith 
49747c6ae99SBarry Smith 
49847c6ae99SBarry Smith 
49947c6ae99SBarry Smith 
500