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