1 2 #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 3 4 /*@C 5 DMDAVecGetArray - Returns a multiple dimension array that shares data with 6 the underlying vector and is indexed using the global dimensions. 7 8 Logically collective on Vec 9 10 Input Parameter: 11 + da - the distributed array 12 - vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector() 13 14 Output Parameter: 15 . array - the array 16 17 Notes: 18 Call DMDAVecRestoreArray() once you have finished accessing the vector entries. 19 20 In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]! 21 22 If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is 23 a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the 24 25 appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations. 26 27 Fortran Notes: 28 From Fortran use DMDAVecGetArrayF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate 29 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 30 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 31 array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from 32 DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine. 33 34 Due to bugs in the compiler DMDAVecGetArrayF90() does not work with gfortran versions before 4.5 35 36 Level: intermediate 37 38 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF() 39 DMDAVecGetArrayDOF() 40 @*/ 41 PetscErrorCode DMDAVecGetArray(DM da,Vec vec,void *array) 42 { 43 PetscErrorCode ierr; 44 PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 45 46 PetscFunctionBegin; 47 PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA); 48 PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 49 PetscValidPointer(array, 3); 50 if (da->defaultSection) { 51 ierr = VecGetArray(vec,(PetscScalar**)array);CHKERRQ(ierr); 52 PetscFunctionReturn(0); 53 } 54 ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 55 ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr); 56 ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 57 58 /* Handle case where user passes in global vector as opposed to local */ 59 ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr); 60 if (N == xm*ym*zm*dof) { 61 gxm = xm; 62 gym = ym; 63 gzm = zm; 64 gxs = xs; 65 gys = ys; 66 gzs = zs; 67 } 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); 68 69 if (dim == 1) { 70 ierr = VecGetArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr); 71 } else if (dim == 2) { 72 ierr = VecGetArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr); 73 } else if (dim == 3) { 74 ierr = VecGetArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr); 75 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim); 76 PetscFunctionReturn(0); 77 } 78 79 /*@ 80 DMDAVecRestoreArray - Restores a multiple dimension array obtained with DMDAVecGetArray() 81 82 Logically collective on Vec 83 84 Input Parameter: 85 + da - the distributed array 86 . vec - the vector, either a vector the same size as one obtained with 87 DMCreateGlobalVector() or DMCreateLocalVector() 88 - array - the array, non-NULL pointer is zeroed 89 90 Level: intermediate 91 92 Fortran Notes: 93 From Fortran use DMDAVecRestoreArrayF90() 94 95 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray() 96 @*/ 97 PetscErrorCode DMDAVecRestoreArray(DM da,Vec vec,void *array) 98 { 99 PetscErrorCode ierr; 100 PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 101 102 PetscFunctionBegin; 103 PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA); 104 PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 105 PetscValidPointer(array, 3); 106 if (da->defaultSection) { 107 ierr = VecRestoreArray(vec,(PetscScalar**)array);CHKERRQ(ierr); 108 PetscFunctionReturn(0); 109 } 110 ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 111 ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr); 112 ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 113 114 /* Handle case where user passes in global vector as opposed to local */ 115 ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr); 116 if (N == xm*ym*zm*dof) { 117 gxm = xm; 118 gym = ym; 119 gzm = zm; 120 gxs = xs; 121 gys = ys; 122 gzs = zs; 123 } 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); 124 125 if (dim == 1) { 126 ierr = VecRestoreArray1d(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr); 127 } else if (dim == 2) { 128 ierr = VecRestoreArray2d(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr); 129 } else if (dim == 3) { 130 ierr = VecRestoreArray3d(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr); 131 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim); 132 PetscFunctionReturn(0); 133 } 134 135 /*@C 136 DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with 137 the underlying vector and is indexed using the global dimensions. 138 139 Logically collective 140 141 Input Parameter: 142 + da - the distributed array 143 - vec - the vector, either a vector the same size as one obtained with 144 DMCreateGlobalVector() or DMCreateLocalVector() 145 146 Output Parameter: 147 . array - the array 148 149 Notes: 150 Call DMDAVecRestoreArrayDOF() once you have finished accessing the vector entries. 151 152 In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]! 153 154 In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use DMDAVecRestoreArrayF90() and declare your array with one higher dimension, 155 see src/dm/examples/tutorials/ex11f90.F 156 157 Level: intermediate 158 159 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF() 160 @*/ 161 PetscErrorCode DMDAVecGetArrayDOF(DM da,Vec vec,void *array) 162 { 163 PetscErrorCode ierr; 164 PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 165 166 PetscFunctionBegin; 167 ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 168 ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr); 169 ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 170 171 /* Handle case where user passes in global vector as opposed to local */ 172 ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr); 173 if (N == xm*ym*zm*dof) { 174 gxm = xm; 175 gym = ym; 176 gzm = zm; 177 gxs = xs; 178 gys = ys; 179 gzs = zs; 180 } 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); 181 182 if (dim == 1) { 183 ierr = VecGetArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr); 184 } else if (dim == 2) { 185 ierr = VecGetArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr); 186 } else if (dim == 3) { 187 ierr = VecGetArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr); 188 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim); 189 PetscFunctionReturn(0); 190 } 191 192 /*@ 193 DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with DMDAVecGetArrayDOF() 194 195 Logically collective 196 197 Input Parameter: 198 + da - the distributed array 199 . vec - the vector, either a vector the same size as one obtained with 200 DMCreateGlobalVector() or DMCreateLocalVector() 201 - array - the array 202 203 Level: intermediate 204 205 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF() 206 @*/ 207 PetscErrorCode DMDAVecRestoreArrayDOF(DM da,Vec vec,void *array) 208 { 209 PetscErrorCode ierr; 210 PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 211 212 PetscFunctionBegin; 213 ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 214 ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr); 215 ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 216 217 /* Handle case where user passes in global vector as opposed to local */ 218 ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr); 219 if (N == xm*ym*zm*dof) { 220 gxm = xm; 221 gym = ym; 222 gzm = zm; 223 gxs = xs; 224 gys = ys; 225 gzs = zs; 226 } 227 228 if (dim == 1) { 229 ierr = VecRestoreArray2d(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr); 230 } else if (dim == 2) { 231 ierr = VecRestoreArray3d(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr); 232 } else if (dim == 3) { 233 ierr = VecRestoreArray4d(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr); 234 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim); 235 PetscFunctionReturn(0); 236 } 237 238 /*@C 239 DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with 240 the underlying vector and is indexed using the global dimensions. 241 242 Not collective 243 244 Input Parameter: 245 + da - the distributed array 246 - vec - the vector, either a vector the same size as one obtained with DMCreateGlobalVector() or DMCreateLocalVector() 247 248 Output Parameter: 249 . array - the array 250 251 Notes: 252 Call DMDAVecRestoreArrayRead() once you have finished accessing the vector entries. 253 254 In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]! 255 256 If vec is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is 257 a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the 258 259 appropriate DMGlobalToLocalBegin() and DMGlobalToLocalEnd() to have correct values in the ghost locations. 260 261 Fortran Notes: 262 From Fortran use DMDAVecGetArrayReadF90() and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate 263 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 264 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 265 array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1) where the values are obtained from 266 DMDAGetCorners() for a global array or DMDAGetGhostCorners() for a local array. Include petsc/finclude/petscdmda.h90 to access this routine. 267 268 Due to bugs in the compiler DMDAVecGetArrayReadF90() does not work with gfortran versions before 4.5 269 270 Level: intermediate 271 272 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecRestoreArrayDOF() 273 DMDAVecGetArrayDOF() 274 @*/ 275 PetscErrorCode DMDAVecGetArrayRead(DM da,Vec vec,void *array) 276 { 277 PetscErrorCode ierr; 278 PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 279 280 PetscFunctionBegin; 281 PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA); 282 PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 283 PetscValidPointer(array, 3); 284 if (da->defaultSection) { 285 ierr = VecGetArrayRead(vec,(const PetscScalar**)array);CHKERRQ(ierr); 286 PetscFunctionReturn(0); 287 } 288 ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 289 ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr); 290 ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 291 292 /* Handle case where user passes in global vector as opposed to local */ 293 ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr); 294 if (N == xm*ym*zm*dof) { 295 gxm = xm; 296 gym = ym; 297 gzm = zm; 298 gxs = xs; 299 gys = ys; 300 gzs = zs; 301 } 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); 302 303 if (dim == 1) { 304 ierr = VecGetArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr); 305 } else if (dim == 2) { 306 ierr = VecGetArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr); 307 } else if (dim == 3) { 308 ierr = VecGetArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr); 309 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim); 310 PetscFunctionReturn(0); 311 } 312 313 /*@ 314 DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with DMDAVecGetArrayRead() 315 316 Not collective 317 318 Input Parameter: 319 + da - the distributed array 320 . vec - the vector, either a vector the same size as one obtained with 321 DMCreateGlobalVector() or DMCreateLocalVector() 322 - array - the array, non-NULL pointer is zeroed 323 324 Level: intermediate 325 326 Fortran Notes: 327 From Fortran use DMDAVecRestoreArrayReadF90() 328 329 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray() 330 @*/ 331 PetscErrorCode DMDAVecRestoreArrayRead(DM da,Vec vec,void *array) 332 { 333 PetscErrorCode ierr; 334 PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 335 336 PetscFunctionBegin; 337 PetscValidHeaderSpecificType(da, DM_CLASSID, 1,DMDA); 338 PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 339 PetscValidPointer(array, 3); 340 if (da->defaultSection) { 341 ierr = VecRestoreArrayRead(vec,(const PetscScalar**)array);CHKERRQ(ierr); 342 PetscFunctionReturn(0); 343 } 344 ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 345 ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr); 346 ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 347 348 /* Handle case where user passes in global vector as opposed to local */ 349 ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr); 350 if (N == xm*ym*zm*dof) { 351 gxm = xm; 352 gym = ym; 353 gzm = zm; 354 gxs = xs; 355 gys = ys; 356 gzs = zs; 357 } 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); 358 359 if (dim == 1) { 360 ierr = VecRestoreArray1dRead(vec,gxm*dof,gxs*dof,(PetscScalar**)array);CHKERRQ(ierr); 361 } else if (dim == 2) { 362 ierr = VecRestoreArray2dRead(vec,gym,gxm*dof,gys,gxs*dof,(PetscScalar***)array);CHKERRQ(ierr); 363 } else if (dim == 3) { 364 ierr = VecRestoreArray3dRead(vec,gzm,gym,gxm*dof,gzs,gys,gxs*dof,(PetscScalar****)array);CHKERRQ(ierr); 365 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim); 366 PetscFunctionReturn(0); 367 } 368 369 /*@C 370 DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with 371 the underlying vector and is indexed using the global dimensions. 372 373 Not Collective 374 375 Input Parameter: 376 + da - the distributed array 377 - vec - the vector, either a vector the same size as one obtained with 378 DMCreateGlobalVector() or DMCreateLocalVector() 379 380 Output Parameter: 381 . array - the array 382 383 Notes: 384 Call DMDAVecRestoreArrayDOFRead() once you have finished accessing the vector entries. 385 386 In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]! 387 388 In Fortran 90 you do not need a version of DMDAVecRestoreArrayDOF() just use DMDAVecRestoreArrayReadF90() and declare your array with one higher dimension, 389 see src/dm/examples/tutorials/ex11f90.F 390 391 Level: intermediate 392 393 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecRestoreArray(), DMDAVecGetArray(), DMDAVecRestoreArrayDOF() 394 @*/ 395 PetscErrorCode DMDAVecGetArrayDOFRead(DM da,Vec vec,void *array) 396 { 397 PetscErrorCode ierr; 398 PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 399 400 PetscFunctionBegin; 401 ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 402 ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr); 403 ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 404 405 /* Handle case where user passes in global vector as opposed to local */ 406 ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr); 407 if (N == xm*ym*zm*dof) { 408 gxm = xm; 409 gym = ym; 410 gzm = zm; 411 gxs = xs; 412 gys = ys; 413 gzs = zs; 414 } 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); 415 416 if (dim == 1) { 417 ierr = VecGetArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr); 418 } else if (dim == 2) { 419 ierr = VecGetArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr); 420 } else if (dim == 3) { 421 ierr = VecGetArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr); 422 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim); 423 PetscFunctionReturn(0); 424 } 425 426 /*@ 427 DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with DMDAVecGetArrayDOFRead() 428 429 Not Collective 430 431 Input Parameter: 432 + da - the distributed array 433 . vec - the vector, either a vector the same size as one obtained with 434 DMCreateGlobalVector() or DMCreateLocalVector() 435 - array - the array 436 437 Level: intermediate 438 439 .seealso: DMDAGetGhostCorners(), DMDAGetCorners(), VecGetArray(), VecRestoreArray(), DMDAVecGetArray(), DMDAVecGetArrayDOF(), DMDAVecRestoreArrayDOF() 440 @*/ 441 PetscErrorCode DMDAVecRestoreArrayDOFRead(DM da,Vec vec,void *array) 442 { 443 PetscErrorCode ierr; 444 PetscInt xs,ys,zs,xm,ym,zm,gxs,gys,gzs,gxm,gym,gzm,N,dim,dof; 445 446 PetscFunctionBegin; 447 ierr = DMDAGetCorners(da,&xs,&ys,&zs,&xm,&ym,&zm);CHKERRQ(ierr); 448 ierr = DMDAGetGhostCorners(da,&gxs,&gys,&gzs,&gxm,&gym,&gzm);CHKERRQ(ierr); 449 ierr = DMDAGetInfo(da,&dim,0,0,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr); 450 451 /* Handle case where user passes in global vector as opposed to local */ 452 ierr = VecGetLocalSize(vec,&N);CHKERRQ(ierr); 453 if (N == xm*ym*zm*dof) { 454 gxm = xm; 455 gym = ym; 456 gzm = zm; 457 gxs = xs; 458 gys = ys; 459 gzs = zs; 460 } 461 462 if (dim == 1) { 463 ierr = VecRestoreArray2dRead(vec,gxm,dof,gxs,0,(PetscScalar***)array);CHKERRQ(ierr); 464 } else if (dim == 2) { 465 ierr = VecRestoreArray3dRead(vec,gym,gxm,dof,gys,gxs,0,(PetscScalar****)array);CHKERRQ(ierr); 466 } else if (dim == 3) { 467 ierr = VecRestoreArray4dRead(vec,gzm,gym,gxm,dof,gzs,gys,gxs,0,(PetscScalar*****)array);CHKERRQ(ierr); 468 } else SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_CORRUPT,"DMDA dimension not 1, 2, or 3, it is %D\n",dim); 469 PetscFunctionReturn(0); 470 } 471 472 473 474 475 476 477 478 479 480 481 482 483 484