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