1 #include <petsc/private/dmdaimpl.h> /*I "petscdmda.h" I*/ 2 3 /*MC 4 DMDAVecGetArrayF90 - see Fortran Notes at `DMDAVecGetArray()` 5 6 Level: intermediate 7 M*/ 8 9 /*@C 10 DMDAVecGetArray - Returns a multiple dimension array that shares data with 11 the underlying vector and is indexed using the global or local dimensions of a `DMDA`. 12 13 Logically Collective 14 15 Input Parameters: 16 + da - the `DMDA` 17 - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 18 19 Output Parameter: 20 . array - the array 21 22 Level: intermediate 23 24 Notes: 25 Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries. 26 27 In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]! 28 29 If `vec` is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is 30 a global vector then the ghost points are not accessible. Of course, with a local vector you will have had to do the 31 appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations. 32 33 The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from 34 `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 35 36 Fortran Notes: 37 Use `DMDAVecGetArray()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate 38 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 39 dimension of the `DMDA`. 40 41 The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise 42 `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 43 `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 44 45 .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()` 46 `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, 47 `DMStagVecGetArray()` 48 @*/ 49 PetscErrorCode DMDAVecGetArray(DM da, Vec vec, void *array) 50 { 51 PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 52 53 PetscFunctionBegin; 54 PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 55 PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 56 PetscAssertPointer(array, 3); 57 PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 58 PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 59 PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 60 61 /* Handle case where user passes in global vector as opposed to local */ 62 PetscCall(VecGetLocalSize(vec, &N)); 63 if (N == xm * ym * zm * dof) { 64 gxm = xm; 65 gym = ym; 66 gzm = zm; 67 gxs = xs; 68 gys = ys; 69 gzs = zs; 70 } 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); 71 72 if (dim == 1) { 73 PetscCall(VecGetArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 74 } else if (dim == 2) { 75 PetscCall(VecGetArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 76 } else if (dim == 3) { 77 PetscCall(VecGetArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 78 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 79 PetscFunctionReturn(PETSC_SUCCESS); 80 } 81 82 /*MC 83 DMDAVecRestoreArrayF90 - see Fortran Notes at `DMDAVecRestoreArray()` 84 85 Level: intermediate 86 M*/ 87 88 /*@ 89 DMDAVecRestoreArray - Restores a multiple dimension array obtained with `DMDAVecGetArray()` 90 91 Logically Collective 92 93 Input Parameters: 94 + da - the `DMDA` 95 . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 96 - array - the `array` pointer 97 98 Level: intermediate 99 100 .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, 101 `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, 102 `DMStagVecRestoreArray()` 103 @*/ 104 PetscErrorCode DMDAVecRestoreArray(DM da, Vec vec, void *array) 105 { 106 PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 107 108 PetscFunctionBegin; 109 PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 110 PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 111 PetscAssertPointer(array, 3); 112 PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 113 PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 114 PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 115 116 /* Handle case where user passes in global vector as opposed to local */ 117 PetscCall(VecGetLocalSize(vec, &N)); 118 if (N == xm * ym * zm * dof) { 119 gxm = xm; 120 gym = ym; 121 gzm = zm; 122 gxs = xs; 123 gys = ys; 124 gzs = zs; 125 } 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); 126 127 if (dim == 1) { 128 PetscCall(VecRestoreArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 129 } else if (dim == 2) { 130 PetscCall(VecRestoreArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 131 } else if (dim == 3) { 132 PetscCall(VecRestoreArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 133 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 134 PetscFunctionReturn(PETSC_SUCCESS); 135 } 136 137 /*MC 138 DMDAVecGetArrayWriteF90 - see Fortran Notes at `DMDAVecGetArrayWrite()` 139 140 Level: intermediate 141 M*/ 142 143 /*@C 144 DMDAVecGetArrayWrite - Returns a multiple dimension array that shares data with 145 the underlying vector and is indexed using the global or local dimensions of a `DMDA`. 146 147 Logically Collective 148 149 Input Parameters: 150 + da - the `DMDA` 151 - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 152 153 Output Parameter: 154 . array - the array 155 156 Level: intermediate 157 158 Notes: 159 Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries. 160 161 In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]! 162 163 if `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is 164 a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the 165 appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations. 166 167 The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from 168 `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 169 170 Fortran Notes: 171 Use `DMDAVecGetArrayWrite()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate 172 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 173 dimension of the `DMDA`. 174 175 The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise 176 `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 177 `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 178 179 Developer Note: 180 This has code duplication with `DMDAVecGetArray()` and `DMDAVecGetArrayRead()` 181 182 .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecRestoreArrayDOF()` 183 `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()` 184 @*/ 185 PetscErrorCode DMDAVecGetArrayWrite(DM da, Vec vec, void *array) 186 { 187 PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 188 189 PetscFunctionBegin; 190 PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 191 PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 192 PetscAssertPointer(array, 3); 193 if (da->localSection) { 194 PetscCall(VecGetArrayWrite(vec, (PetscScalar **)array)); 195 PetscFunctionReturn(PETSC_SUCCESS); 196 } 197 PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 198 PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 199 PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 200 201 /* Handle case where user passes in global vector as opposed to local */ 202 PetscCall(VecGetLocalSize(vec, &N)); 203 if (N == xm * ym * zm * dof) { 204 gxm = xm; 205 gym = ym; 206 gzm = zm; 207 gxs = xs; 208 gys = ys; 209 gzs = zs; 210 } 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); 211 212 if (dim == 1) { 213 PetscCall(VecGetArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 214 } else if (dim == 2) { 215 PetscCall(VecGetArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 216 } else if (dim == 3) { 217 PetscCall(VecGetArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 218 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 219 PetscFunctionReturn(PETSC_SUCCESS); 220 } 221 222 /*MC 223 DMDAVecRestoreArrayWriteF90 - see Fortran Notes at `DMDAVecRestoreArrayWrite()` 224 225 Level: intermediate 226 M*/ 227 228 /*@ 229 DMDAVecRestoreArrayWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayWrite()` 230 231 Logically Collective 232 233 Input Parameters: 234 + da - the `DMDA` 235 . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 236 - array - the `array` pointer 237 238 Level: intermediate 239 240 .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayWrite()`, 241 `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()` 242 @*/ 243 PetscErrorCode DMDAVecRestoreArrayWrite(DM da, Vec vec, void *array) 244 { 245 PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 246 247 PetscFunctionBegin; 248 PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 249 PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 250 PetscAssertPointer(array, 3); 251 if (da->localSection) { 252 PetscCall(VecRestoreArray(vec, (PetscScalar **)array)); 253 PetscFunctionReturn(PETSC_SUCCESS); 254 } 255 PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 256 PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 257 PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 258 259 /* Handle case where user passes in global vector as opposed to local */ 260 PetscCall(VecGetLocalSize(vec, &N)); 261 if (N == xm * ym * zm * dof) { 262 gxm = xm; 263 gym = ym; 264 gzm = zm; 265 gxs = xs; 266 gys = ys; 267 gzs = zs; 268 } 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); 269 270 if (dim == 1) { 271 PetscCall(VecRestoreArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 272 } else if (dim == 2) { 273 PetscCall(VecRestoreArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 274 } else if (dim == 3) { 275 PetscCall(VecRestoreArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 276 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 277 PetscFunctionReturn(PETSC_SUCCESS); 278 } 279 280 /*@C 281 DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with 282 the underlying vector and is indexed using the global or local dimensions of a `DMDA` 283 284 Logically Collective 285 286 Input Parameters: 287 + da - the `DMDA` 288 - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 289 290 Output Parameter: 291 . array - the `array` pointer 292 293 Level: intermediate 294 295 Notes: 296 Call `DMDAVecRestoreArrayDOF()` once you have finished accessing the vector entries. 297 298 In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF] 299 300 The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1][0:ndof-1]` where the values are obtained from 301 `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 302 303 Fortran Notes: 304 Use `DMDAVecGetArray()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate 305 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 306 dimension of the `DMDA`. 307 308 The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when ndof is 1) otherwise 309 `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 310 `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 311 312 .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecRestoreArrayDOF()`, 313 `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, `DMDAVecGetArrayDOFRead()` 314 @*/ 315 PetscErrorCode DMDAVecGetArrayDOF(DM da, Vec vec, void *array) 316 { 317 PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 318 319 PetscFunctionBegin; 320 PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 321 PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 322 PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 323 324 /* Handle case where user passes in global vector as opposed to local */ 325 PetscCall(VecGetLocalSize(vec, &N)); 326 if (N == xm * ym * zm * dof) { 327 gxm = xm; 328 gym = ym; 329 gzm = zm; 330 gxs = xs; 331 gys = ys; 332 gzs = zs; 333 } 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); 334 335 if (dim == 1) { 336 PetscCall(VecGetArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array)); 337 } else if (dim == 2) { 338 PetscCall(VecGetArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array)); 339 } else if (dim == 3) { 340 PetscCall(VecGetArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array)); 341 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 342 PetscFunctionReturn(PETSC_SUCCESS); 343 } 344 345 /*@ 346 DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOF()` 347 348 Logically Collective 349 350 Input Parameters: 351 + da - the `DMDA` 352 . vec - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 353 - array - the `array` point 354 355 Level: intermediate 356 357 .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, 358 `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()` 359 @*/ 360 PetscErrorCode DMDAVecRestoreArrayDOF(DM da, Vec vec, void *array) 361 { 362 PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 363 364 PetscFunctionBegin; 365 PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 366 PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 367 PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 368 369 /* Handle case where user passes in global vector as opposed to local */ 370 PetscCall(VecGetLocalSize(vec, &N)); 371 if (N == xm * ym * zm * dof) { 372 gxm = xm; 373 gym = ym; 374 gzm = zm; 375 gxs = xs; 376 gys = ys; 377 gzs = zs; 378 } 379 380 if (dim == 1) { 381 PetscCall(VecRestoreArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array)); 382 } else if (dim == 2) { 383 PetscCall(VecRestoreArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array)); 384 } else if (dim == 3) { 385 PetscCall(VecRestoreArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array)); 386 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 387 PetscFunctionReturn(PETSC_SUCCESS); 388 } 389 390 /*MC 391 DMDAVecGetArrayReadF90 - see Fortran Notes at `DMDAVecGetArrayRead()` 392 393 Level: intermediate 394 M*/ 395 396 /*@C 397 DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with 398 the underlying vector and is indexed using the global or local dimensions of a `DMDA`. 399 400 Not Collective 401 402 Input Parameters: 403 + da - the `DMDA` 404 - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 405 406 Output Parameter: 407 . array - the array 408 409 Level: intermediate 410 411 Notes: 412 Call `DMDAVecRestoreArrayRead()` once you have finished accessing the vector entries. 413 414 In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]! 415 416 If `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is 417 a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the 418 appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations. 419 420 The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from 421 `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 422 423 Fortran Notes: 424 Use `DMDAVecGetArrayRead()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate 425 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 426 dimension of the `DMDA`. 427 428 The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise 429 `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 430 `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 431 432 .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, 433 `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayRead()`, 434 `DMDAVecRestoreArrayDOF()`, `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, 435 `DMDAVecRestoreArray()`, `DMStagVecGetArrayRead()` 436 @*/ 437 PetscErrorCode DMDAVecGetArrayRead(DM da, Vec vec, void *array) 438 { 439 PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 440 441 PetscFunctionBegin; 442 PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 443 PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 444 PetscAssertPointer(array, 3); 445 PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 446 PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 447 PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 448 449 /* Handle case where user passes in global vector as opposed to local */ 450 PetscCall(VecGetLocalSize(vec, &N)); 451 if (N == xm * ym * zm * dof) { 452 gxm = xm; 453 gym = ym; 454 gzm = zm; 455 gxs = xs; 456 gys = ys; 457 gzs = zs; 458 } 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); 459 460 if (dim == 1) { 461 PetscCall(VecGetArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 462 } else if (dim == 2) { 463 PetscCall(VecGetArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 464 } else if (dim == 3) { 465 PetscCall(VecGetArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 466 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 467 PetscFunctionReturn(PETSC_SUCCESS); 468 } 469 470 /*MC 471 DMDAVecRestoreArrayReadF90 - see Fortran Notes at `DMDAVecRestoreArrayRead()` 472 473 Level: intermediate 474 M*/ 475 476 /*@ 477 DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayRead()` 478 479 Not Collective 480 481 Input Parameters: 482 + da - the `DMDA` 483 . vec - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 484 - array - the `array` pointer 485 486 Level: intermediate 487 488 .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayRead()`, 489 `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, 490 `DMStagVecRestoreArrayRead()` 491 @*/ 492 PetscErrorCode DMDAVecRestoreArrayRead(DM da, Vec vec, void *array) 493 { 494 PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 495 496 PetscFunctionBegin; 497 PetscValidHeaderSpecificType(da, DM_CLASSID, 1, DMDA); 498 PetscValidHeaderSpecific(vec, VEC_CLASSID, 2); 499 PetscAssertPointer(array, 3); 500 PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 501 PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 502 PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 503 504 /* Handle case where user passes in global vector as opposed to local */ 505 PetscCall(VecGetLocalSize(vec, &N)); 506 if (N == xm * ym * zm * dof) { 507 gxm = xm; 508 gym = ym; 509 gzm = zm; 510 gxs = xs; 511 gys = ys; 512 gzs = zs; 513 } 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); 514 515 if (dim == 1) { 516 PetscCall(VecRestoreArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array)); 517 } else if (dim == 2) { 518 PetscCall(VecRestoreArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array)); 519 } else if (dim == 3) { 520 PetscCall(VecRestoreArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array)); 521 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 522 PetscFunctionReturn(PETSC_SUCCESS); 523 } 524 525 /*@C 526 DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with 527 the underlying vector and is indexed using the global or local dimensions of a `DMDA` 528 529 Not Collective 530 531 Input Parameters: 532 + da - the `DMDA` 533 - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 534 535 Output Parameter: 536 . array - the array 537 538 Level: intermediate 539 540 Notes: 541 Call `DMDAVecRestoreArrayDOFRead()` once you have finished accessing the vector entries. 542 543 In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]! 544 545 The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from 546 `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 547 548 Fortran Notes: 549 Use `DMDAVecGetArrayRead()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate 550 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 551 dimension of the `DMDA`. 552 553 The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise 554 `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 555 `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 556 557 .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, 558 `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()` 559 @*/ 560 PetscErrorCode DMDAVecGetArrayDOFRead(DM da, Vec vec, void *array) 561 { 562 PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 563 564 PetscFunctionBegin; 565 PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 566 PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 567 PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 568 569 /* Handle case where user passes in global vector as opposed to local */ 570 PetscCall(VecGetLocalSize(vec, &N)); 571 if (N == xm * ym * zm * dof) { 572 gxm = xm; 573 gym = ym; 574 gzm = zm; 575 gxs = xs; 576 gys = ys; 577 gzs = zs; 578 } 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); 579 580 if (dim == 1) { 581 PetscCall(VecGetArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array)); 582 } else if (dim == 2) { 583 PetscCall(VecGetArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array)); 584 } else if (dim == 3) { 585 PetscCall(VecGetArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array)); 586 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 587 PetscFunctionReturn(PETSC_SUCCESS); 588 } 589 590 /*@ 591 DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFRead()` 592 593 Not Collective 594 595 Input Parameters: 596 + da - the `DMDA` 597 . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 598 - array - the `array` pointer 599 600 Level: intermediate 601 602 .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`, 603 `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()` 604 @*/ 605 PetscErrorCode DMDAVecRestoreArrayDOFRead(DM da, Vec vec, void *array) 606 { 607 PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 608 609 PetscFunctionBegin; 610 PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 611 PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 612 PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 613 614 /* Handle case where user passes in global vector as opposed to local */ 615 PetscCall(VecGetLocalSize(vec, &N)); 616 if (N == xm * ym * zm * dof) { 617 gxm = xm; 618 gym = ym; 619 gzm = zm; 620 gxs = xs; 621 gys = ys; 622 gzs = zs; 623 } 624 625 if (dim == 1) { 626 PetscCall(VecRestoreArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array)); 627 } else if (dim == 2) { 628 PetscCall(VecRestoreArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array)); 629 } else if (dim == 3) { 630 PetscCall(VecRestoreArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array)); 631 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 632 PetscFunctionReturn(PETSC_SUCCESS); 633 } 634 635 /*@C 636 DMDAVecGetArrayDOFWrite - Returns a multiple dimension array that shares data with 637 the underlying vector and is indexed using the global or local dimensions of a `DMDA` 638 639 Not Collective 640 641 Input Parameters: 642 + da - the `DMDA` 643 - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 644 645 Output Parameter: 646 . array - the array 647 648 Level: intermediate 649 650 Notes: 651 Call `DMDAVecRestoreArrayDOFWrite()` once you have finished accessing the vector entries. 652 653 In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]! 654 655 The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1][0:dof-1]` where the values are obtained from 656 `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 657 658 Fortran Notes: 659 Use `DMDAVecGetArrayWrite()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate 660 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 661 dimension of the `DMDA`. 662 663 The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise 664 `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from 665 `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector. 666 667 .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, 668 `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()` 669 @*/ 670 PetscErrorCode DMDAVecGetArrayDOFWrite(DM da, Vec vec, void *array) 671 { 672 PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 673 674 PetscFunctionBegin; 675 PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 676 PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 677 PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 678 679 /* Handle case where user passes in global vector as opposed to local */ 680 PetscCall(VecGetLocalSize(vec, &N)); 681 if (N == xm * ym * zm * dof) { 682 gxm = xm; 683 gym = ym; 684 gzm = zm; 685 gxs = xs; 686 gys = ys; 687 gzs = zs; 688 } 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); 689 690 if (dim == 1) { 691 PetscCall(VecGetArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array)); 692 } else if (dim == 2) { 693 PetscCall(VecGetArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array)); 694 } else if (dim == 3) { 695 PetscCall(VecGetArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array)); 696 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 697 PetscFunctionReturn(PETSC_SUCCESS); 698 } 699 700 /*@ 701 DMDAVecRestoreArrayDOFWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFWrite()` 702 703 Not Collective 704 705 Input Parameters: 706 + da - the `DMDA` 707 . vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()` 708 - array - the `array` pointer 709 710 Level: intermediate 711 712 .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`, 713 `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()` 714 @*/ 715 PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM da, Vec vec, void *array) 716 { 717 PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof; 718 719 PetscFunctionBegin; 720 PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm)); 721 PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm)); 722 PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL)); 723 724 /* Handle case where user passes in global vector as opposed to local */ 725 PetscCall(VecGetLocalSize(vec, &N)); 726 if (N == xm * ym * zm * dof) { 727 gxm = xm; 728 gym = ym; 729 gzm = zm; 730 gxs = xs; 731 gys = ys; 732 gzs = zs; 733 } 734 735 if (dim == 1) { 736 PetscCall(VecRestoreArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array)); 737 } else if (dim == 2) { 738 PetscCall(VecRestoreArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array)); 739 } else if (dim == 3) { 740 PetscCall(VecRestoreArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array)); 741 } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim); 742 PetscFunctionReturn(PETSC_SUCCESS); 743 } 744