Actual source code: dagetarray.c
  1: #include <petsc/private/dmdaimpl.h>
  3: /*@C
  4:   DMDAVecGetArray - Returns a multiple dimension array that shares data with
  5:   the underlying vector and is indexed using the global or local dimensions of a `DMDA`.
  7:   Logically Collective
  9:   Input Parameters:
 10: + da  - the `DMDA`
 11: - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
 13:   Output Parameter:
 14: . array - the array
 16:   Level: intermediate
 18:   Notes:
 19:   Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
 21:   In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
 23:   If `vec` is a local vector (obtained with DMCreateLocalVector() etc) then the ghost point locations are accessible. If it is
 24:   a global vector then the ghost points are not accessible. Of course, with a local vector you will have had to do the
 25:   appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
 27:   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
 28:   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
 30:   Fortran Notes:
 31:   Use `DMDAVecGetArray()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
 32:   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
 33:   dimension of the `DMDA`.
 35:   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise
 36:   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
 37:   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
 39: .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecRestoreArrayDOF()`
 40:           `DMDAVecGetArrayDOF()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
 41:           `DMStagVecGetArray()`
 42: @*/
 43: PetscErrorCode DMDAVecGetArray(DM da, Vec vec, void *array)
 44: {
 45:   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
 47:   PetscFunctionBegin;
 50:   PetscAssertPointer(array, 3);
 51:   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
 52:   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
 53:   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
 55:   /* Handle case where user passes in global vector as opposed to local */
 56:   PetscCall(VecGetLocalSize(vec, &N));
 57:   if (N == xm * ym * zm * dof) {
 58:     gxm = xm;
 59:     gym = ym;
 60:     gzm = zm;
 61:     gxs = xs;
 62:     gys = ys;
 63:     gzs = zs;
 64:   } 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);
 66:   if (dim == 1) {
 67:     PetscCall(VecGetArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
 68:   } else if (dim == 2) {
 69:     PetscCall(VecGetArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
 70:   } else if (dim == 3) {
 71:     PetscCall(VecGetArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
 72:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
 73:   PetscFunctionReturn(PETSC_SUCCESS);
 74: }
 76: /*@
 77:   DMDAVecRestoreArray - Restores a multiple dimension array obtained with `DMDAVecGetArray()`
 79:   Logically Collective
 81:   Input Parameters:
 82: + da    - the `DMDA`
 83: . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
 84: - array - the `array` pointer
 86:   Level: intermediate
 88: .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`,
 89:           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`,
 90:           `DMStagVecRestoreArray()`
 91: @*/
 92: PetscErrorCode DMDAVecRestoreArray(DM da, Vec vec, void *array)
 93: {
 94:   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
 96:   PetscFunctionBegin;
 99:   PetscAssertPointer(array, 3);
100:   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
101:   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
102:   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
104:   /* Handle case where user passes in global vector as opposed to local */
105:   PetscCall(VecGetLocalSize(vec, &N));
106:   if (N == xm * ym * zm * dof) {
107:     gxm = xm;
108:     gym = ym;
109:     gzm = zm;
110:     gxs = xs;
111:     gys = ys;
112:     gzs = zs;
113:   } 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);
115:   if (dim == 1) {
116:     PetscCall(VecRestoreArray1d(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
117:   } else if (dim == 2) {
118:     PetscCall(VecRestoreArray2d(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
119:   } else if (dim == 3) {
120:     PetscCall(VecRestoreArray3d(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
121:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
122:   PetscFunctionReturn(PETSC_SUCCESS);
123: }
125: /*@C
126:   DMDAVecGetArrayWrite - Returns a multiple dimension array that shares data with
127:   the underlying vector and is indexed using the global or local dimensions of a `DMDA`.
129:   Logically Collective
131:   Input Parameters:
132: + da  - the `DMDA`
133: - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
135:   Output Parameter:
136: . array - the array
138:   Level: intermediate
140:   Notes:
141:   Call `DMDAVecRestoreArray()` once you have finished accessing the vector entries.
143:   In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
145:   if `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
146:   a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
147:   appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
149:   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
150:   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
152:   Fortran Notes:
153:   Use `DMDAVecGetArrayWrite()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
154:   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
155:   dimension of the `DMDA`.
157:   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise
158:   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
159:   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
161:   Developer Note:
162:   This has code duplication with `DMDAVecGetArray()` and `DMDAVecGetArrayRead()`
164: .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecRestoreArrayDOF()`
165:           `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
166: @*/
167: PetscErrorCode DMDAVecGetArrayWrite(DM da, Vec vec, void *array)
168: {
169:   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
171:   PetscFunctionBegin;
174:   PetscAssertPointer(array, 3);
175:   if (da->localSection) {
176:     PetscCall(VecGetArrayWrite(vec, (PetscScalar **)array));
177:     PetscFunctionReturn(PETSC_SUCCESS);
178:   }
179:   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
180:   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
181:   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
183:   /* Handle case where user passes in global vector as opposed to local */
184:   PetscCall(VecGetLocalSize(vec, &N));
185:   if (N == xm * ym * zm * dof) {
186:     gxm = xm;
187:     gym = ym;
188:     gzm = zm;
189:     gxs = xs;
190:     gys = ys;
191:     gzs = zs;
192:   } 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);
194:   if (dim == 1) {
195:     PetscCall(VecGetArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
196:   } else if (dim == 2) {
197:     PetscCall(VecGetArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
198:   } else if (dim == 3) {
199:     PetscCall(VecGetArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
200:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
201:   PetscFunctionReturn(PETSC_SUCCESS);
202: }
204: /*@
205:   DMDAVecRestoreArrayWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayWrite()`
207:   Logically Collective
209:   Input Parameters:
210: + da    - the `DMDA`
211: . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
212: - array - the `array` pointer
214:   Level: intermediate
216: .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayWrite()`,
217:           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
218: @*/
219: PetscErrorCode DMDAVecRestoreArrayWrite(DM da, Vec vec, void *array)
220: {
221:   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
223:   PetscFunctionBegin;
226:   PetscAssertPointer(array, 3);
227:   if (da->localSection) {
228:     PetscCall(VecRestoreArray(vec, (PetscScalar **)array));
229:     PetscFunctionReturn(PETSC_SUCCESS);
230:   }
231:   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
232:   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
233:   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
235:   /* Handle case where user passes in global vector as opposed to local */
236:   PetscCall(VecGetLocalSize(vec, &N));
237:   if (N == xm * ym * zm * dof) {
238:     gxm = xm;
239:     gym = ym;
240:     gzm = zm;
241:     gxs = xs;
242:     gys = ys;
243:     gzs = zs;
244:   } 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);
246:   if (dim == 1) {
247:     PetscCall(VecRestoreArray1dWrite(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
248:   } else if (dim == 2) {
249:     PetscCall(VecRestoreArray2dWrite(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
250:   } else if (dim == 3) {
251:     PetscCall(VecRestoreArray3dWrite(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
252:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
253:   PetscFunctionReturn(PETSC_SUCCESS);
254: }
256: /*@C
257:   DMDAVecGetArrayDOF - Returns a multiple dimension array that shares data with
258:   the underlying vector and is indexed using the global or local dimensions of a `DMDA`
260:   Logically Collective
262:   Input Parameters:
263: + da  - the `DMDA`
264: - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
266:   Output Parameter:
267: . array - the `array` pointer
269:   Level: intermediate
271:   Notes:
272:   Call `DMDAVecRestoreArrayDOF()` once you have finished accessing the vector entries.
274:   In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]
276:   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
277:   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
279:   Fortran Notes:
280:   Use `DMDAVecGetArray()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
281:   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
282:   dimension of the `DMDA`.
284:   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when ndof is 1) otherwise
285:   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
286:   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
288: .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecRestoreArrayDOF()`,
289:           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`, `DMDAVecGetArrayDOFRead()`
290: @*/
291: PetscErrorCode DMDAVecGetArrayDOF(DM da, Vec vec, void *array)
292: {
293:   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
295:   PetscFunctionBegin;
296:   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
297:   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
298:   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
300:   /* Handle case where user passes in global vector as opposed to local */
301:   PetscCall(VecGetLocalSize(vec, &N));
302:   if (N == xm * ym * zm * dof) {
303:     gxm = xm;
304:     gym = ym;
305:     gzm = zm;
306:     gxs = xs;
307:     gys = ys;
308:     gzs = zs;
309:   } 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);
311:   if (dim == 1) {
312:     PetscCall(VecGetArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
313:   } else if (dim == 2) {
314:     PetscCall(VecGetArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
315:   } else if (dim == 3) {
316:     PetscCall(VecGetArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
317:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
318:   PetscFunctionReturn(PETSC_SUCCESS);
319: }
321: /*@
322:   DMDAVecRestoreArrayDOF - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOF()`
324:   Logically Collective
326:   Input Parameters:
327: + da    - the `DMDA`
328: . vec   - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
329: - array - the `array` point
331:   Level: intermediate
333: .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
334:           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
335: @*/
336: PetscErrorCode DMDAVecRestoreArrayDOF(DM da, Vec vec, void *array)
337: {
338:   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
340:   PetscFunctionBegin;
341:   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
342:   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
343:   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
345:   /* Handle case where user passes in global vector as opposed to local */
346:   PetscCall(VecGetLocalSize(vec, &N));
347:   if (N == xm * ym * zm * dof) {
348:     gxm = xm;
349:     gym = ym;
350:     gzm = zm;
351:     gxs = xs;
352:     gys = ys;
353:     gzs = zs;
354:   }
356:   if (dim == 1) {
357:     PetscCall(VecRestoreArray2d(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
358:   } else if (dim == 2) {
359:     PetscCall(VecRestoreArray3d(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
360:   } else if (dim == 3) {
361:     PetscCall(VecRestoreArray4d(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
362:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
363:   PetscFunctionReturn(PETSC_SUCCESS);
364: }
366: /*@C
367:   DMDAVecGetArrayRead - Returns a multiple dimension array that shares data with
368:   the underlying vector and is indexed using the global or local dimensions of a `DMDA`.
370:   Not Collective
372:   Input Parameters:
373: + da  - the `DMDA`
374: - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
376:   Output Parameter:
377: . array - the array
379:   Level: intermediate
381:   Notes:
382:   Call `DMDAVecRestoreArrayRead()` once you have finished accessing the vector entries.
384:   In C, the indexing is "backwards" from what expects: array[k][j][i] NOT array[i][j][k]!
386:   If `vec` is a local vector (obtained with `DMCreateLocalVector()` etc) then the ghost point locations are accessible. If it is
387:   a global vector then the ghost points are not accessible. Of course with the local vector you will have had to do the
388:   appropriate `DMGlobalToLocalBegin()` and `DMGlobalToLocalEnd()` to have correct values in the ghost locations.
390:   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
391:   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
393:   Fortran Notes:
394:   Use `DMDAVecGetArrayRead()` and pass for the array type `PetscScalar`,pointer :: array(:,...,:) of the appropriate
395:   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
396:   dimension of the `DMDA`.
398:   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise
399:   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
400:   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
402: .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`,
403: `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArrayRead()`,
404: `DMDAVecRestoreArrayDOF()`, `DMDAVecGetArrayDOF()`, `DMDAVecGetArray()`,
405: `DMDAVecRestoreArray()`, `DMStagVecGetArrayRead()`
406: @*/
407: PetscErrorCode DMDAVecGetArrayRead(DM da, Vec vec, void *array)
408: {
409:   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
411:   PetscFunctionBegin;
414:   PetscAssertPointer(array, 3);
415:   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
416:   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
417:   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
419:   /* Handle case where user passes in global vector as opposed to local */
420:   PetscCall(VecGetLocalSize(vec, &N));
421:   if (N == xm * ym * zm * dof) {
422:     gxm = xm;
423:     gym = ym;
424:     gzm = zm;
425:     gxs = xs;
426:     gys = ys;
427:     gzs = zs;
428:   } 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);
430:   if (dim == 1) {
431:     PetscCall(VecGetArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
432:   } else if (dim == 2) {
433:     PetscCall(VecGetArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
434:   } else if (dim == 3) {
435:     PetscCall(VecGetArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
436:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
437:   PetscFunctionReturn(PETSC_SUCCESS);
438: }
440: /*@
441:   DMDAVecRestoreArrayRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayRead()`
443:   Not Collective
445:   Input Parameters:
446: + da    - the `DMDA`
447: . vec   - vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
448: - array - the `array` pointer
450:   Level: intermediate
452: .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArrayRead()`,
453:           `DMDAVecGetArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`,
454:           `DMStagVecRestoreArrayRead()`
455: @*/
456: PetscErrorCode DMDAVecRestoreArrayRead(DM da, Vec vec, void *array)
457: {
458:   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
460:   PetscFunctionBegin;
463:   PetscAssertPointer(array, 3);
464:   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
465:   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
466:   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
468:   /* Handle case where user passes in global vector as opposed to local */
469:   PetscCall(VecGetLocalSize(vec, &N));
470:   if (N == xm * ym * zm * dof) {
471:     gxm = xm;
472:     gym = ym;
473:     gzm = zm;
474:     gxs = xs;
475:     gys = ys;
476:     gzs = zs;
477:   } 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);
479:   if (dim == 1) {
480:     PetscCall(VecRestoreArray1dRead(vec, gxm * dof, gxs * dof, (PetscScalar **)array));
481:   } else if (dim == 2) {
482:     PetscCall(VecRestoreArray2dRead(vec, gym, gxm * dof, gys, gxs * dof, (PetscScalar ***)array));
483:   } else if (dim == 3) {
484:     PetscCall(VecRestoreArray3dRead(vec, gzm, gym, gxm * dof, gzs, gys, gxs * dof, (PetscScalar ****)array));
485:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
486:   PetscFunctionReturn(PETSC_SUCCESS);
487: }
489: /*@C
490:   DMDAVecGetArrayDOFRead - Returns a multiple dimension array that shares data with
491:   the underlying vector and is indexed using the global or local dimensions of a `DMDA`
493:   Not Collective
495:   Input Parameters:
496: + da  - the `DMDA`
497: - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
499:   Output Parameter:
500: . array - the array
502:   Level: intermediate
504:   Notes:
505:   Call `DMDAVecRestoreArrayDOFRead()` once you have finished accessing the vector entries.
507:   In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
509:   The accessible indices are `array[zs:zs+zm-1][ys:ys+ym-1][xs:xs+xm-1]` where the values are obtained from
510:   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
512:   Fortran Notes:
513:   Use  `DMDAVecGetArrayRead()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
514:   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
515:   dimension of the `DMDA`.
517:   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise
518:   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
519:   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
521: .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
522:           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
523: @*/
524: PetscErrorCode DMDAVecGetArrayDOFRead(DM da, Vec vec, void *array)
525: {
526:   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
528:   PetscFunctionBegin;
529:   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
530:   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
531:   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
533:   /* Handle case where user passes in global vector as opposed to local */
534:   PetscCall(VecGetLocalSize(vec, &N));
535:   if (N == xm * ym * zm * dof) {
536:     gxm = xm;
537:     gym = ym;
538:     gzm = zm;
539:     gxs = xs;
540:     gys = ys;
541:     gzs = zs;
542:   } 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);
544:   if (dim == 1) {
545:     PetscCall(VecGetArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
546:   } else if (dim == 2) {
547:     PetscCall(VecGetArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
548:   } else if (dim == 3) {
549:     PetscCall(VecGetArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
550:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
551:   PetscFunctionReturn(PETSC_SUCCESS);
552: }
554: /*@
555:   DMDAVecRestoreArrayDOFRead - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFRead()`
557:   Not Collective
559:   Input Parameters:
560: + da    - the `DMDA`
561: . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
562: - array - the `array` pointer
564:   Level: intermediate
566: .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
567:           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`, `DMDAVecGetArrayRead()`, `DMDAVecRestoreArrayRead()`
568: @*/
569: PetscErrorCode DMDAVecRestoreArrayDOFRead(DM da, Vec vec, void *array)
570: {
571:   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
573:   PetscFunctionBegin;
574:   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
575:   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
576:   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
578:   /* Handle case where user passes in global vector as opposed to local */
579:   PetscCall(VecGetLocalSize(vec, &N));
580:   if (N == xm * ym * zm * dof) {
581:     gxm = xm;
582:     gym = ym;
583:     gzm = zm;
584:     gxs = xs;
585:     gys = ys;
586:     gzs = zs;
587:   }
589:   if (dim == 1) {
590:     PetscCall(VecRestoreArray2dRead(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
591:   } else if (dim == 2) {
592:     PetscCall(VecRestoreArray3dRead(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
593:   } else if (dim == 3) {
594:     PetscCall(VecRestoreArray4dRead(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
595:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
596:   PetscFunctionReturn(PETSC_SUCCESS);
597: }
599: /*@C
600:   DMDAVecGetArrayDOFWrite - Returns a multiple dimension array that shares data with
601:   the underlying vector and is indexed using the global or local dimensions of a `DMDA`
603:   Not Collective
605:   Input Parameters:
606: + da  - the `DMDA`
607: - vec - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
609:   Output Parameter:
610: . array - the array
612:   Level: intermediate
614:   Notes:
615:   Call `DMDAVecRestoreArrayDOFWrite()` once you have finished accessing the vector entries.
617:   In C, the indexing is "backwards" from what expects: array[k][j][i][DOF] NOT array[i][j][k][DOF]!
619:   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
620:   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
622:   Fortran Notes:
623:   Use  `DMDAVecGetArrayWrite()` and pass for the array type PetscScalar,pointer :: array(:,...,:) of the appropriate
624:   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
625:   dimension of the `DMDA`.
627:   The order of the indices is `array(xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` (when dof is 1) otherwise
628:   `array(0:dof-1,xs:xs+xm-1,ys:ys+ym-1,zs:zs+zm-1)` where the values are obtained from
629:   `DMDAGetCorners()` for a global vector or `DMDAGetGhostCorners()` for a local vector.
631: .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`,
632:           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
633: @*/
634: PetscErrorCode DMDAVecGetArrayDOFWrite(DM da, Vec vec, void *array)
635: {
636:   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
638:   PetscFunctionBegin;
639:   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
640:   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
641:   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
643:   /* Handle case where user passes in global vector as opposed to local */
644:   PetscCall(VecGetLocalSize(vec, &N));
645:   if (N == xm * ym * zm * dof) {
646:     gxm = xm;
647:     gym = ym;
648:     gzm = zm;
649:     gxs = xs;
650:     gys = ys;
651:     gzs = zs;
652:   } 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);
654:   if (dim == 1) {
655:     PetscCall(VecGetArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
656:   } else if (dim == 2) {
657:     PetscCall(VecGetArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
658:   } else if (dim == 3) {
659:     PetscCall(VecGetArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
660:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
661:   PetscFunctionReturn(PETSC_SUCCESS);
662: }
664: /*@
665:   DMDAVecRestoreArrayDOFWrite - Restores a multiple dimension array obtained with `DMDAVecGetArrayDOFWrite()`
667:   Not Collective
669:   Input Parameters:
670: + da    - the `DMDA`
671: . vec   - a vector the same size as one obtained with `DMCreateGlobalVector()` or `DMCreateLocalVector()`
672: - array - the `array` pointer
674:   Level: intermediate
676: .seealso: [](sec_struct), [](sec_struct_set), `DM`, `DMDA`, `DMDAGetGhostCorners()`, `DMDAGetCorners()`, `VecGetArray()`, `VecRestoreArray()`, `DMDAVecGetArray()`, `DMDAVecGetArrayDOF()`, `DMDAVecRestoreArrayDOF()`,
677:           `DMDAVecGetArrayWrite()`, `DMDAVecRestoreArrayWrite()`
678: @*/
679: PetscErrorCode DMDAVecRestoreArrayDOFWrite(DM da, Vec vec, void *array)
680: {
681:   PetscInt xs, ys, zs, xm, ym, zm, gxs, gys, gzs, gxm, gym, gzm, N, dim, dof;
683:   PetscFunctionBegin;
684:   PetscCall(DMDAGetCorners(da, &xs, &ys, &zs, &xm, &ym, &zm));
685:   PetscCall(DMDAGetGhostCorners(da, &gxs, &gys, &gzs, &gxm, &gym, &gzm));
686:   PetscCall(DMDAGetInfo(da, &dim, NULL, NULL, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
688:   /* Handle case where user passes in global vector as opposed to local */
689:   PetscCall(VecGetLocalSize(vec, &N));
690:   if (N == xm * ym * zm * dof) {
691:     gxm = xm;
692:     gym = ym;
693:     gzm = zm;
694:     gxs = xs;
695:     gys = ys;
696:     gzs = zs;
697:   }
699:   if (dim == 1) {
700:     PetscCall(VecRestoreArray2dWrite(vec, gxm, dof, gxs, 0, (PetscScalar ***)array));
701:   } else if (dim == 2) {
702:     PetscCall(VecRestoreArray3dWrite(vec, gym, gxm, dof, gys, gxs, 0, (PetscScalar ****)array));
703:   } else if (dim == 3) {
704:     PetscCall(VecRestoreArray4dWrite(vec, gzm, gym, gxm, dof, gzs, gys, gxs, 0, (PetscScalar *****)array));
705:   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_CORRUPT, "DMDA dimension not 1, 2, or 3, it is %" PetscInt_FMT, dim);
706:   PetscFunctionReturn(PETSC_SUCCESS);
707: }