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 8d083f849SBarry Smith Logically collective on da 947c6ae99SBarry Smith 10d8d19677SJose E. Roman Input Parameters: 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() 393bd220d7SPatrick Sanan DMDAVecGetArrayDOF(), DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), 403bd220d7SPatrick Sanan DMStagVecGetArray() 4147c6ae99SBarry Smith @*/ 427087cfbeSBarry Smith PetscErrorCode DMDAVecGetArray(DM da,Vec vec,void *array) 4347c6ae99SBarry Smith { 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); 509566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 519566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm)); 529566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL)); 5347c6ae99SBarry Smith 5447c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 559566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec,&N)); 5647c6ae99SBarry Smith if (N == xm*ym*zm*dof) { 5747c6ae99SBarry Smith gxm = xm; 5847c6ae99SBarry Smith gym = ym; 5947c6ae99SBarry Smith gzm = zm; 6047c6ae99SBarry Smith gxs = xs; 6147c6ae99SBarry Smith gys = ys; 6247c6ae99SBarry Smith gzs = zs; 63*63a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm*gym*gzm*dof,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT,N,xm*ym*zm*dof,gxm*gym*gzm*dof); 6447c6ae99SBarry Smith 6547c6ae99SBarry Smith if (dim == 1) { 669566063dSJacob Faibussowitsch PetscCall(VecGetArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array)); 6747c6ae99SBarry Smith } else if (dim == 2) { 689566063dSJacob Faibussowitsch PetscCall(VecGetArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array)); 6947c6ae99SBarry Smith } else if (dim == 3) { 709566063dSJacob Faibussowitsch PetscCall(VecGetArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array)); 71*63a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT,dim); 7247c6ae99SBarry Smith PetscFunctionReturn(0); 7347c6ae99SBarry Smith } 7447c6ae99SBarry Smith 7547c6ae99SBarry Smith /*@ 76aa219208SBarry Smith DMDAVecRestoreArray - Restores a multiple dimension array obtained with DMDAVecGetArray() 7747c6ae99SBarry Smith 78d083f849SBarry Smith Logically collective on da 7947c6ae99SBarry Smith 80d8d19677SJose E. Roman Input Parameters: 8147c6ae99SBarry Smith + da - the distributed array 8247c6ae99SBarry Smith . vec - the vector, either a vector the same size as one obtained with 83564755cdSBarry Smith DMCreateGlobalVector() or DMCreateLocalVector() 84e5c84f05SJed Brown - array - the array, non-NULL pointer is zeroed 8547c6ae99SBarry Smith 8647c6ae99SBarry Smith Level: intermediate 8747c6ae99SBarry Smith 8895452b02SPatrick Sanan Fortran Notes: 89fdc842d1SBarry Smith From Fortran use DMDAVecRestoreArayF90() 9047c6ae99SBarry Smith 91fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), 923bd220d7SPatrick Sanan DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), 933bd220d7SPatrick Sanan DMDStagVecRestoreArray() 9447c6ae99SBarry Smith @*/ 957087cfbeSBarry Smith PetscErrorCode DMDAVecRestoreArray(DM da,Vec vec,void *array) 9647c6ae99SBarry Smith { 9747c6ae99SBarry Smith PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 9847c6ae99SBarry Smith 9947c6ae99SBarry Smith PetscFunctionBegin; 100a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA); 1016db82c94SMatthew G Knepley PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 1026db82c94SMatthew G Knepley PetscValidPointer(array, 3); 1039566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 1049566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm)); 1059566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL)); 10647c6ae99SBarry Smith 10747c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 1089566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec,&N)); 10947c6ae99SBarry Smith if (N == xm*ym*zm*dof) { 11047c6ae99SBarry Smith gxm = xm; 11147c6ae99SBarry Smith gym = ym; 11247c6ae99SBarry Smith gzm = zm; 11347c6ae99SBarry Smith gxs = xs; 11447c6ae99SBarry Smith gys = ys; 11547c6ae99SBarry Smith gzs = zs; 116*63a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm*gym*gzm*dof,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT,N,xm*ym*zm*dof,gxm*gym*gzm*dof); 11747c6ae99SBarry Smith 11847c6ae99SBarry Smith if (dim == 1) { 1199566063dSJacob Faibussowitsch PetscCall(VecRestoreArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array)); 12047c6ae99SBarry Smith } else if (dim == 2) { 1219566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array)); 12247c6ae99SBarry Smith } else if (dim == 3) { 1239566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array)); 124*63a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT,dim); 12547c6ae99SBarry Smith PetscFunctionReturn(0); 12647c6ae99SBarry Smith } 12747c6ae99SBarry Smith 12847c6ae99SBarry Smith /*@C 129fdc842d1SBarry Smith DMDAVecGetArrayWrite - Returns a multiple dimension array that shares data with 130fdc842d1SBarry Smith the underlying vector and is indexed using the global dimensions. 131fdc842d1SBarry Smith 132fdc842d1SBarry Smith Logically collective on Vec 133fdc842d1SBarry Smith 134d8d19677SJose E. Roman Input Parameters: 135fdc842d1SBarry Smith + da - the distributed array 136fdc842d1SBarry Smith - vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector() 137fdc842d1SBarry Smith 138fdc842d1SBarry Smith Output Parameter: 139fdc842d1SBarry Smith . array - the array 140fdc842d1SBarry Smith 141fdc842d1SBarry Smith Notes: 142fdc842d1SBarry Smith Call DMDAVecRestoreArray() once you have finished accessing the vector entries. 143fdc842d1SBarry Smith 144fdc842d1SBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]! 145fdc842d1SBarry Smith 146fdc842d1SBarry Smith If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is 147fdc842d1SBarry Smith a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the 148fdc842d1SBarry Smith 149fdc842d1SBarry Smith appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations. 150fdc842d1SBarry Smith 151fdc842d1SBarry Smith Fortran Notes: 152fdc842d1SBarry Smith From Fortran use DMDAVecGetArrayF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate 153fdc842d1SBarry 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 154fdc842d1SBarry 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 155fdc842d1SBarry Smith array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from 156fdc842d1SBarry Smith DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine. 157fdc842d1SBarry Smith 158fdc842d1SBarry Smith Due to bugs in the compiler DMDAVecGetArrayF90() does not work with gfortran versions before 4.5 159fdc842d1SBarry Smith 160fdc842d1SBarry Smith Level: intermediate 161fdc842d1SBarry Smith 162fdc842d1SBarry Smith Developer Notes: This has code duplication with DMDAVecGetArray() and DMDAVecGetArrayRead() 163fdc842d1SBarry Smith 164fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArrayWrite(), DMDAVecRestoreArrayDOF() 165fdc842d1SBarry Smith DMDAVecGetArrayDOF(), DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead() 166fdc842d1SBarry Smith @*/ 167fdc842d1SBarry Smith PetscErrorCode DMDAVecGetArrayWrite(DM da,Vec vec,void *array) 168fdc842d1SBarry Smith { 169fdc842d1SBarry Smith PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 170fdc842d1SBarry Smith 171fdc842d1SBarry Smith PetscFunctionBegin; 172fdc842d1SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA); 173fdc842d1SBarry Smith PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 174fdc842d1SBarry Smith PetscValidPointer(array, 3); 1751bb6d2a8SBarry Smith if (da->localSection) { 1769566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(vec,(PetscScalar**)array)); 177fdc842d1SBarry Smith PetscFunctionReturn(0); 178fdc842d1SBarry Smith } 1799566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 1809566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm)); 1819566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL)); 182fdc842d1SBarry Smith 183fdc842d1SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 1849566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec,&N)); 185fdc842d1SBarry Smith if (N == xm*ym*zm*dof) { 186fdc842d1SBarry Smith gxm = xm; 187fdc842d1SBarry Smith gym = ym; 188fdc842d1SBarry Smith gzm = zm; 189fdc842d1SBarry Smith gxs = xs; 190fdc842d1SBarry Smith gys = ys; 191fdc842d1SBarry Smith gzs = zs; 192*63a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm*gym*gzm*dof,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT,N,xm*ym*zm*dof,gxm*gym*gzm*dof); 193fdc842d1SBarry Smith 194fdc842d1SBarry Smith if (dim == 1) { 1959566063dSJacob Faibussowitsch PetscCall(VecGetArray1dWrite(vec,gxm*dof,gxs*dof,(PetscScalar**)array)); 196fdc842d1SBarry Smith } else if (dim == 2) { 1979566063dSJacob Faibussowitsch PetscCall(VecGetArray2dWrite(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array)); 198fdc842d1SBarry Smith } else if (dim == 3) { 1999566063dSJacob Faibussowitsch PetscCall(VecGetArray3dWrite(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array)); 200*63a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT,dim); 201fdc842d1SBarry Smith PetscFunctionReturn(0); 202fdc842d1SBarry Smith } 203fdc842d1SBarry Smith 204fdc842d1SBarry Smith /*@ 205fdc842d1SBarry Smith DMDAVecRestoreArrayWrite - Restores a multiple dimension array obtained with DMDAVecGetArrayWrite() 206fdc842d1SBarry Smith 207fdc842d1SBarry Smith Logically collective on Vec 208fdc842d1SBarry Smith 209d8d19677SJose E. Roman Input Parameters: 210fdc842d1SBarry Smith + da - the distributed array 211fdc842d1SBarry Smith . vec - the vector, either a vector the same size as one obtained with 212fdc842d1SBarry Smith DMCreateGlobalVector() or DMCreateLocalVector() 213fdc842d1SBarry Smith - array - the array, non-NULL pointer is zeroed 214fdc842d1SBarry Smith 215fdc842d1SBarry Smith Level: intermediate 216fdc842d1SBarry Smith 217fdc842d1SBarry Smith Fortran Notes: 218fdc842d1SBarry Smith From Fortran use DMDAVecRestoreArayF90() 219fdc842d1SBarry Smith 220fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArrayWrite(), 221fdc842d1SBarry Smith DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead() 222fdc842d1SBarry Smith @*/ 223fdc842d1SBarry Smith PetscErrorCode DMDAVecRestoreArrayWrite(DM da,Vec vec,void *array) 224fdc842d1SBarry Smith { 225fdc842d1SBarry Smith PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 226fdc842d1SBarry Smith 227fdc842d1SBarry Smith PetscFunctionBegin; 228fdc842d1SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA); 229fdc842d1SBarry Smith PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 230fdc842d1SBarry Smith PetscValidPointer(array, 3); 2311bb6d2a8SBarry Smith if (da->localSection) { 2329566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(vec,(PetscScalar**)array)); 233fdc842d1SBarry Smith PetscFunctionReturn(0); 234fdc842d1SBarry Smith } 2359566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 2369566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm)); 2379566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL)); 238fdc842d1SBarry Smith 239fdc842d1SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 2409566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec,&N)); 241fdc842d1SBarry Smith if (N == xm*ym*zm*dof) { 242fdc842d1SBarry Smith gxm = xm; 243fdc842d1SBarry Smith gym = ym; 244fdc842d1SBarry Smith gzm = zm; 245fdc842d1SBarry Smith gxs = xs; 246fdc842d1SBarry Smith gys = ys; 247fdc842d1SBarry Smith gzs = zs; 248*63a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm*gym*gzm*dof,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT,N,xm*ym*zm*dof,gxm*gym*gzm*dof); 249fdc842d1SBarry Smith 250fdc842d1SBarry Smith if (dim == 1) { 2519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray1dWrite(vec,gxm*dof,gxs*dof,(PetscScalar**)array)); 252fdc842d1SBarry Smith } else if (dim == 2) { 2539566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2dWrite(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array)); 254fdc842d1SBarry Smith } else if (dim == 3) { 2559566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3dWrite(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array)); 256*63a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT,dim); 257fdc842d1SBarry Smith PetscFunctionReturn(0); 258fdc842d1SBarry Smith } 259fdc842d1SBarry Smith 260fdc842d1SBarry Smith /*@C 261aa219208SBarry Smith DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with 26247c6ae99SBarry Smith the underlying vector and is indexed using the global dimensions. 26347c6ae99SBarry Smith 2642ddbf755SMatthew G. Knepley Logically collective 26547c6ae99SBarry Smith 266d8d19677SJose E. Roman Input Parameters: 26747c6ae99SBarry Smith + da - the distributed array 26847c6ae99SBarry Smith - vec - the vector, either a vector the same size as one obtained with 269564755cdSBarry Smith DMCreateGlobalVector() or DMCreateLocalVector() 27047c6ae99SBarry Smith 27147c6ae99SBarry Smith Output Parameter: 27247c6ae99SBarry Smith . array - the array 27347c6ae99SBarry Smith 27447c6ae99SBarry Smith Notes: 275aa219208SBarry Smith Call DMDAVecRestoreArrayDOF() once you have finished accessing the vector entries. 27647c6ae99SBarry Smith 27747c6ae99SBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]! 27847c6ae99SBarry Smith 2791b82215eSBarry Smith In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use DMDAVecRestoreArrayF90() and declare your array with one higher dimension, 280c4762a1bSJed Brown see src/dm/tutorials/ex11f90.F 2811b82215eSBarry Smith 28247c6ae99SBarry Smith Level: intermediate 28347c6ae99SBarry Smith 284fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF(), 285856ac3e9SJed Brown DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), DMStagVecGetArrayDOF(), DMDAVecGetArrayDOFRead() 28647c6ae99SBarry Smith @*/ 2877087cfbeSBarry Smith PetscErrorCode DMDAVecGetArrayDOF(DM da,Vec vec,void *array) 28847c6ae99SBarry Smith { 28947c6ae99SBarry Smith PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 29047c6ae99SBarry Smith 29147c6ae99SBarry Smith PetscFunctionBegin; 2929566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 2939566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm)); 2949566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL)); 29547c6ae99SBarry Smith 29647c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 2979566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec,&N)); 29847c6ae99SBarry Smith if (N == xm*ym*zm*dof) { 29947c6ae99SBarry Smith gxm = xm; 30047c6ae99SBarry Smith gym = ym; 30147c6ae99SBarry Smith gzm = zm; 30247c6ae99SBarry Smith gxs = xs; 30347c6ae99SBarry Smith gys = ys; 30447c6ae99SBarry Smith gzs = zs; 305*63a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm*gym*gzm*dof,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT,N,xm*ym*zm*dof,gxm*gym*gzm*dof); 30647c6ae99SBarry Smith 30747c6ae99SBarry Smith if (dim == 1) { 3089566063dSJacob Faibussowitsch PetscCall(VecGetArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array)); 30947c6ae99SBarry Smith } else if (dim == 2) { 3109566063dSJacob Faibussowitsch PetscCall(VecGetArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array)); 31147c6ae99SBarry Smith } else if (dim == 3) { 3129566063dSJacob Faibussowitsch PetscCall(VecGetArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array)); 313*63a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT,dim); 31447c6ae99SBarry Smith PetscFunctionReturn(0); 31547c6ae99SBarry Smith } 31647c6ae99SBarry Smith 31747c6ae99SBarry Smith /*@ 318aa219208SBarry Smith DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with DMDAVecGetArrayDOF() 31947c6ae99SBarry Smith 3202ddbf755SMatthew G. Knepley Logically collective 32147c6ae99SBarry Smith 322d8d19677SJose E. Roman Input Parameters: 32347c6ae99SBarry Smith + da - the distributed array 32447c6ae99SBarry Smith . vec - the vector, either a vector the same size as one obtained with 325564755cdSBarry Smith DMCreateGlobalVector() or DMCreateLocalVector() 32647c6ae99SBarry Smith - array - the array 32747c6ae99SBarry Smith 32847c6ae99SBarry Smith Level: intermediate 32947c6ae99SBarry Smith 330fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF(), 331856ac3e9SJed Brown DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), DMStagVecRestoreArrayDOF(), DMDAVecRestoreArrayDOFRead() 33247c6ae99SBarry Smith @*/ 3337087cfbeSBarry Smith PetscErrorCode DMDAVecRestoreArrayDOF(DM da,Vec vec,void *array) 33447c6ae99SBarry Smith { 33547c6ae99SBarry Smith PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 33647c6ae99SBarry Smith 33747c6ae99SBarry Smith PetscFunctionBegin; 3389566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 3399566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm)); 3409566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL)); 34147c6ae99SBarry Smith 34247c6ae99SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 3439566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec,&N)); 34447c6ae99SBarry Smith if (N == xm*ym*zm*dof) { 34547c6ae99SBarry Smith gxm = xm; 34647c6ae99SBarry Smith gym = ym; 34747c6ae99SBarry Smith gzm = zm; 34847c6ae99SBarry Smith gxs = xs; 34947c6ae99SBarry Smith gys = ys; 35047c6ae99SBarry Smith gzs = zs; 35147c6ae99SBarry Smith } 35247c6ae99SBarry Smith 35347c6ae99SBarry Smith if (dim == 1) { 3549566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array)); 35547c6ae99SBarry Smith } else if (dim == 2) { 3569566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array)); 35747c6ae99SBarry Smith } else if (dim == 3) { 3589566063dSJacob Faibussowitsch PetscCall(VecRestoreArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array)); 359*63a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT,dim); 36047c6ae99SBarry Smith PetscFunctionReturn(0); 36147c6ae99SBarry Smith } 36247c6ae99SBarry Smith 3635edff71fSBarry Smith /*@C 3645edff71fSBarry Smith DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with 3655edff71fSBarry Smith the underlying vector and is indexed using the global dimensions. 3665edff71fSBarry Smith 3672ddbf755SMatthew G. Knepley Not collective 3685edff71fSBarry Smith 369d8d19677SJose E. Roman Input Parameters: 3705edff71fSBarry Smith + da - the distributed array 3715edff71fSBarry Smith - vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector() 3725edff71fSBarry Smith 3735edff71fSBarry Smith Output Parameter: 3745edff71fSBarry Smith . array - the array 3755edff71fSBarry Smith 3765edff71fSBarry Smith Notes: 3775edff71fSBarry Smith Call DMDAVecRestoreArrayRead() once you have finished accessing the vector entries. 3785edff71fSBarry Smith 3795edff71fSBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]! 3805edff71fSBarry Smith 3815edff71fSBarry Smith If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is 3825edff71fSBarry Smith a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the 3835edff71fSBarry Smith 3845edff71fSBarry Smith appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations. 3855edff71fSBarry Smith 38695452b02SPatrick Sanan Fortran Notes: 38795452b02SPatrick Sanan From Fortran use DMDAVecGetArrayReadF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate 3885edff71fSBarry 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 3895edff71fSBarry 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 3905edff71fSBarry Smith array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from 391af0996ceSBarry Smith DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine. 3925edff71fSBarry Smith 3934ebdaf59SBarry Smith Due to bugs in the compiler DMDAVecGetArrayReadF90() does not work with gfortran versions before 4.5 3945edff71fSBarry Smith 3955edff71fSBarry Smith Level: intermediate 3965edff71fSBarry Smith 397fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArrayRead(), DMDAVecRestoreArrayDOF() 3983bd220d7SPatrick Sanan DMDAVecGetArrayDOF(), DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), 3993bd220d7SPatrick Sanan DMStagVecGetArrayRead() 4005edff71fSBarry Smith @*/ 4015edff71fSBarry Smith PetscErrorCode DMDAVecGetArrayRead(DM da,Vec vec,void *array) 4025edff71fSBarry Smith { 4035edff71fSBarry Smith PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 4045edff71fSBarry Smith 4055edff71fSBarry Smith PetscFunctionBegin; 406a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA); 4075edff71fSBarry Smith PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 4085edff71fSBarry Smith PetscValidPointer(array, 3); 4099566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 4109566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm)); 4119566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL)); 4125edff71fSBarry Smith 4135edff71fSBarry Smith /* Handle case where user passes in global vector as opposed to local */ 4149566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec,&N)); 4155edff71fSBarry Smith if (N == xm*ym*zm*dof) { 4165edff71fSBarry Smith gxm = xm; 4175edff71fSBarry Smith gym = ym; 4185edff71fSBarry Smith gzm = zm; 4195edff71fSBarry Smith gxs = xs; 4205edff71fSBarry Smith gys = ys; 4215edff71fSBarry Smith gzs = zs; 422*63a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm*gym*gzm*dof,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT,N,xm*ym*zm*dof,gxm*gym*gzm*dof); 4235edff71fSBarry Smith 4245edff71fSBarry Smith if (dim == 1) { 4259566063dSJacob Faibussowitsch PetscCall(VecGetArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array)); 4265edff71fSBarry Smith } else if (dim == 2) { 4279566063dSJacob Faibussowitsch PetscCall(VecGetArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array)); 4285edff71fSBarry Smith } else if (dim == 3) { 4299566063dSJacob Faibussowitsch PetscCall(VecGetArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array)); 430*63a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT,dim); 4315edff71fSBarry Smith PetscFunctionReturn(0); 4325edff71fSBarry Smith } 4335edff71fSBarry Smith 4345edff71fSBarry Smith /*@ 4355edff71fSBarry Smith DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with DMDAVecGetArrayRead() 4365edff71fSBarry Smith 4372ddbf755SMatthew G. Knepley Not collective 4385edff71fSBarry Smith 439d8d19677SJose E. Roman Input Parameters: 4405edff71fSBarry Smith + da - the distributed array 4415edff71fSBarry Smith . vec - the vector, either a vector the same size as one obtained with 4425edff71fSBarry Smith DMCreateGlobalVector() or DMCreateLocalVector() 4435edff71fSBarry Smith - array - the array, non-NULL pointer is zeroed 4445edff71fSBarry Smith 4455edff71fSBarry Smith Level: intermediate 4465edff71fSBarry Smith 44795452b02SPatrick Sanan Fortran Notes: 44895452b02SPatrick Sanan From Fortran use DMDAVecRestoreArrayReadF90() 4495edff71fSBarry Smith 450fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArrayRead(), 4513bd220d7SPatrick Sanan DMDAVecGetArray(), DMDAVecRestoreArray(), DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), 4523bd220d7SPatrick Sanan DMStagVecRestoreArrayRead() 4535edff71fSBarry Smith @*/ 4545edff71fSBarry Smith PetscErrorCode DMDAVecRestoreArrayRead(DM da,Vec vec,void *array) 4555edff71fSBarry Smith { 4565edff71fSBarry Smith PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 4575edff71fSBarry Smith 4585edff71fSBarry Smith PetscFunctionBegin; 459a9a02de4SBarry Smith PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA); 4605edff71fSBarry Smith PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 4615edff71fSBarry Smith PetscValidPointer(array, 3); 4629566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 4639566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm)); 4649566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL)); 4655edff71fSBarry Smith 4665edff71fSBarry Smith /* Handle case where user passes in global vector as opposed to local */ 4679566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec,&N)); 4685edff71fSBarry Smith if (N == xm*ym*zm*dof) { 4695edff71fSBarry Smith gxm = xm; 4705edff71fSBarry Smith gym = ym; 4715edff71fSBarry Smith gzm = zm; 4725edff71fSBarry Smith gxs = xs; 4735edff71fSBarry Smith gys = ys; 4745edff71fSBarry Smith gzs = zs; 475*63a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm*gym*gzm*dof,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT,N,xm*ym*zm*dof,gxm*gym*gzm*dof); 4765edff71fSBarry Smith 4775edff71fSBarry Smith if (dim == 1) { 4789566063dSJacob Faibussowitsch PetscCall(VecRestoreArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array)); 4795edff71fSBarry Smith } else if (dim == 2) { 4809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array)); 4815edff71fSBarry Smith } else if (dim == 3) { 4829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array)); 483*63a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT,dim); 4845edff71fSBarry Smith PetscFunctionReturn(0); 4855edff71fSBarry Smith } 4865edff71fSBarry Smith 4875edff71fSBarry Smith /*@C 4885edff71fSBarry Smith DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with 4895edff71fSBarry Smith the underlying vector and is indexed using the global dimensions. 4905edff71fSBarry Smith 4915edff71fSBarry Smith Not Collective 4925edff71fSBarry Smith 493d8d19677SJose E. Roman Input Parameters: 4945edff71fSBarry Smith + da - the distributed array 4955edff71fSBarry Smith - vec - the vector, either a vector the same size as one obtained with 4965edff71fSBarry Smith DMCreateGlobalVector() or DMCreateLocalVector() 4975edff71fSBarry Smith 4985edff71fSBarry Smith Output Parameter: 4995edff71fSBarry Smith . array - the array 5005edff71fSBarry Smith 5015edff71fSBarry Smith Notes: 5025edff71fSBarry Smith Call DMDAVecRestoreArrayDOFRead() once you have finished accessing the vector entries. 5035edff71fSBarry Smith 5045edff71fSBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]! 5055edff71fSBarry Smith 5064ebdaf59SBarry Smith In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use DMDAVecRestoreArrayReadF90() and declare your array with one higher dimension, 507c4762a1bSJed Brown see src/dm/tutorials/ex11f90.F 5085edff71fSBarry Smith 5095edff71fSBarry Smith Level: intermediate 5105edff71fSBarry Smith 511856ac3e9SJed Brown .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), 5123bd220d7SPatrick Sanan DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), DMStagVecGetArrayDOFRead() 5135edff71fSBarry Smith @*/ 5145edff71fSBarry Smith PetscErrorCode DMDAVecGetArrayDOFRead(DM da,Vec vec,void *array) 5155edff71fSBarry Smith { 5165edff71fSBarry Smith PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 5175edff71fSBarry Smith 5185edff71fSBarry Smith PetscFunctionBegin; 5199566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 5209566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm)); 5219566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL)); 5225edff71fSBarry Smith 5235edff71fSBarry Smith /* Handle case where user passes in global vector as opposed to local */ 5249566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec,&N)); 5255edff71fSBarry Smith if (N == xm*ym*zm*dof) { 5265edff71fSBarry Smith gxm = xm; 5275edff71fSBarry Smith gym = ym; 5285edff71fSBarry Smith gzm = zm; 5295edff71fSBarry Smith gxs = xs; 5305edff71fSBarry Smith gys = ys; 5315edff71fSBarry Smith gzs = zs; 532*63a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm*gym*gzm*dof,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT,N,xm*ym*zm*dof,gxm*gym*gzm*dof); 5335edff71fSBarry Smith 5345edff71fSBarry Smith if (dim == 1) { 5359566063dSJacob Faibussowitsch PetscCall(VecGetArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array)); 5365edff71fSBarry Smith } else if (dim == 2) { 5379566063dSJacob Faibussowitsch PetscCall(VecGetArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array)); 5385edff71fSBarry Smith } else if (dim == 3) { 5399566063dSJacob Faibussowitsch PetscCall(VecGetArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array)); 540*63a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT,dim); 5415edff71fSBarry Smith PetscFunctionReturn(0); 5425edff71fSBarry Smith } 5435edff71fSBarry Smith 5445edff71fSBarry Smith /*@ 5455edff71fSBarry Smith DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with DMDAVecGetArrayDOFRead() 5465edff71fSBarry Smith 5475edff71fSBarry Smith Not Collective 5485edff71fSBarry Smith 549d8d19677SJose E. Roman Input Parameters: 5505edff71fSBarry Smith + da - the distributed array 5515edff71fSBarry Smith . vec - the vector, either a vector the same size as one obtained with 5525edff71fSBarry Smith DMCreateGlobalVector() or DMCreateLocalVector() 5535edff71fSBarry Smith - array - the array 5545edff71fSBarry Smith 5555edff71fSBarry Smith Level: intermediate 5565edff71fSBarry Smith 557fdc842d1SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF(), 5583bd220d7SPatrick Sanan DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayRead(), DMDAVecRestoreArrayRead(), DMStagVecRestoreArrayDOFRead() 5595edff71fSBarry Smith @*/ 5605edff71fSBarry Smith PetscErrorCode DMDAVecRestoreArrayDOFRead(DM da,Vec vec,void *array) 5615edff71fSBarry Smith { 5625edff71fSBarry Smith PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 5635edff71fSBarry Smith 5645edff71fSBarry Smith PetscFunctionBegin; 5659566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 5669566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm)); 5679566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL)); 5685edff71fSBarry Smith 5695edff71fSBarry Smith /* Handle case where user passes in global vector as opposed to local */ 5709566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec,&N)); 5715edff71fSBarry Smith if (N == xm*ym*zm*dof) { 5725edff71fSBarry Smith gxm = xm; 5735edff71fSBarry Smith gym = ym; 5745edff71fSBarry Smith gzm = zm; 5755edff71fSBarry Smith gxs = xs; 5765edff71fSBarry Smith gys = ys; 5775edff71fSBarry Smith gzs = zs; 5785edff71fSBarry Smith } 5795edff71fSBarry Smith 5805edff71fSBarry Smith if (dim == 1) { 5819566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array)); 5825edff71fSBarry Smith } else if (dim == 2) { 5839566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array)); 5845edff71fSBarry Smith } else if (dim == 3) { 5859566063dSJacob Faibussowitsch PetscCall(VecRestoreArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array)); 586*63a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT,dim); 5875edff71fSBarry Smith PetscFunctionReturn(0); 5885edff71fSBarry Smith } 5895edff71fSBarry Smith 5901e5d2365SBarry Smith /*@C 5911e5d2365SBarry Smith DMDAVecGetArrayDOFWrite - Returns a multiple dimension array that shares data with 5921e5d2365SBarry Smith the underlying vector and is indexed using the global dimensions. 5931e5d2365SBarry Smith 5941e5d2365SBarry Smith Not Collective 5951e5d2365SBarry Smith 596d8d19677SJose E. Roman Input Parameters: 5971e5d2365SBarry Smith + da - the distributed array 5981e5d2365SBarry Smith - vec - the vector, either a vector the same size as one obtained with 5991e5d2365SBarry Smith DMCreateGlobalVector() or DMCreateLocalVector() 6001e5d2365SBarry Smith 6011e5d2365SBarry Smith Output Parameter: 6021e5d2365SBarry Smith . array - the array 6031e5d2365SBarry Smith 6041e5d2365SBarry Smith Notes: 6051e5d2365SBarry Smith Call DMDAVecRestoreArrayDOFWrite() once you have finished accessing the vector entries. 6061e5d2365SBarry Smith 6071e5d2365SBarry Smith In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]! 6081e5d2365SBarry Smith 6091e5d2365SBarry Smith In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use DMDAVecRestoreArrayWriteF90() and declare your array with one higher dimension, 6101e5d2365SBarry Smith see src/dm/tutorials/ex11f90.F 6111e5d2365SBarry Smith 6121e5d2365SBarry Smith Level: intermediate 6131e5d2365SBarry Smith 6141e5d2365SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), 6151e5d2365SBarry Smith DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMStagVecGetArrayDOFWrite(), 6161e5d2365SBarry Smith DMStagVecRestoreArrayDOFRead(), DMStagVecRestoreArrayDOFRead() 6171e5d2365SBarry Smith @*/ 6181e5d2365SBarry Smith PetscErrorCode DMDAVecGetArrayDOFWrite(DM da,Vec vec,void *array) 6191e5d2365SBarry Smith { 6201e5d2365SBarry Smith PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 6211e5d2365SBarry Smith 6221e5d2365SBarry Smith PetscFunctionBegin; 6239566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 6249566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm)); 6259566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL)); 6261e5d2365SBarry Smith 6271e5d2365SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 6289566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec,&N)); 6291e5d2365SBarry Smith if (N == xm*ym*zm*dof) { 6301e5d2365SBarry Smith gxm = xm; 6311e5d2365SBarry Smith gym = ym; 6321e5d2365SBarry Smith gzm = zm; 6331e5d2365SBarry Smith gxs = xs; 6341e5d2365SBarry Smith gys = ys; 6351e5d2365SBarry Smith gzs = zs; 636*63a3b9bcSJacob Faibussowitsch } else PetscCheck(N == gxm*gym*gzm*dof,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Vector local size %" PetscInt_FMT " is not compatible with DMDA local sizes %" PetscInt_FMT " %" PetscInt_FMT,N,xm*ym*zm*dof,gxm*gym*gzm*dof); 6371e5d2365SBarry Smith 6381e5d2365SBarry Smith if (dim == 1) { 6399566063dSJacob Faibussowitsch PetscCall(VecGetArray2dWrite(vec,gxm,dof,gxs,0,(PetscScalar***)array)); 6401e5d2365SBarry Smith } else if (dim == 2) { 6419566063dSJacob Faibussowitsch PetscCall(VecGetArray3dWrite(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array)); 6421e5d2365SBarry Smith } else if (dim == 3) { 6439566063dSJacob Faibussowitsch PetscCall(VecGetArray4dWrite(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array)); 644*63a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT,dim); 6451e5d2365SBarry Smith PetscFunctionReturn(0); 6461e5d2365SBarry Smith } 6471e5d2365SBarry Smith 6481e5d2365SBarry Smith /*@ 6491e5d2365SBarry Smith DMDAVecRestoreArrayDOFWrite - Restores a multiple dimension array obtained with DMDAVecGetArrayDOFWrite() 6501e5d2365SBarry Smith 6511e5d2365SBarry Smith Not Collective 6521e5d2365SBarry Smith 653d8d19677SJose E. Roman Input Parameters: 6541e5d2365SBarry Smith + da - the distributed array 6551e5d2365SBarry Smith . vec - the vector, either a vector the same size as one obtained with 6561e5d2365SBarry Smith DMCreateGlobalVector() or DMCreateLocalVector() 6571e5d2365SBarry Smith - array - the array 6581e5d2365SBarry Smith 6591e5d2365SBarry Smith Level: intermediate 6601e5d2365SBarry Smith 6611e5d2365SBarry Smith .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF(), 6621e5d2365SBarry Smith DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMDAVecGetArrayWrite(), DMDAVecRestoreArrayWrite(), DMStagVecRestoreArrayDOFWrite(), 6631e5d2365SBarry Smith DMStagVecRestoreArrayDOFRead(), DMStagVecRestoreArrayDOFRead() 6641e5d2365SBarry Smith @*/ 6651e5d2365SBarry Smith PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM da,Vec vec,void *array) 6661e5d2365SBarry Smith { 6671e5d2365SBarry Smith PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 6681e5d2365SBarry Smith 6691e5d2365SBarry Smith PetscFunctionBegin; 6709566063dSJacob Faibussowitsch PetscCall(DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm)); 6719566063dSJacob Faibussowitsch PetscCall(DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm)); 6729566063dSJacob Faibussowitsch PetscCall(DMDAGetInfo(da,&dim,NULL,NULL,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL)); 6731e5d2365SBarry Smith 6741e5d2365SBarry Smith /* Handle case where user passes in global vector as opposed to local */ 6759566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(vec,&N)); 6761e5d2365SBarry Smith if (N == xm*ym*zm*dof) { 6771e5d2365SBarry Smith gxm = xm; 6781e5d2365SBarry Smith gym = ym; 6791e5d2365SBarry Smith gzm = zm; 6801e5d2365SBarry Smith gxs = xs; 6811e5d2365SBarry Smith gys = ys; 6821e5d2365SBarry Smith gzs = zs; 6831e5d2365SBarry Smith } 6841e5d2365SBarry Smith 6851e5d2365SBarry Smith if (dim == 1) { 6869566063dSJacob Faibussowitsch PetscCall(VecRestoreArray2dWrite(vec,gxm,dof,gxs,0,(PetscScalar***)array)); 6871e5d2365SBarry Smith } else if (dim == 2) { 6889566063dSJacob Faibussowitsch PetscCall(VecRestoreArray3dWrite(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array)); 6891e5d2365SBarry Smith } else if (dim == 3) { 6909566063dSJacob Faibussowitsch PetscCall(VecRestoreArray4dWrite(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array)); 691*63a3b9bcSJacob Faibussowitsch } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT,dim); 6921e5d2365SBarry Smith PetscFunctionReturn(0); 6931e5d2365SBarry Smith } 694