Actual source code: daltol.c
1: /*
2: Code for manipulating distributed regular arrays in parallel.
3: */
5: #include <petsc/private/dmdaimpl.h>
7: /*
8: DMLocalToLocalCreate_DA - Creates the local to local scatter
10: Collective
12: Input Parameter:
13: . da - the `DMDA`
15: */
16: PetscErrorCode DMLocalToLocalCreate_DA(DM da)
17: {
18: PetscInt *idx, left, j, count, up, down, i, bottom, top, k, dim = da->dim;
19: DM_DA *dd = (DM_DA *)da->data;
21: PetscFunctionBegin;
24: if (dd->ltol) PetscFunctionReturn(PETSC_SUCCESS);
25: /*
26: We simply remap the values in the from part of
27: global to local to read from an array with the ghost values
28: rather than from the plain array.
29: */
30: PetscCall(VecScatterCopy(dd->gtol, &dd->ltol));
31: if (dim == 1) {
32: left = dd->xs - dd->Xs;
33: PetscCall(PetscMalloc1(dd->xe - dd->xs, &idx));
34: for (j = 0; j < dd->xe - dd->xs; j++) idx[j] = left + j;
35: } else if (dim == 2) {
36: left = dd->xs - dd->Xs;
37: down = dd->ys - dd->Ys;
38: up = down + dd->ye - dd->ys;
39: PetscCall(PetscMalloc1((dd->xe - dd->xs) * (up - down), &idx));
40: count = 0;
41: for (i = down; i < up; i++) {
42: for (j = 0; j < dd->xe - dd->xs; j++) idx[count++] = left + i * (dd->Xe - dd->Xs) + j;
43: }
44: } else if (dim == 3) {
45: left = dd->xs - dd->Xs;
46: bottom = dd->ys - dd->Ys;
47: top = bottom + dd->ye - dd->ys;
48: down = dd->zs - dd->Zs;
49: up = down + dd->ze - dd->zs;
50: count = (dd->xe - dd->xs) * (top - bottom) * (up - down);
51: PetscCall(PetscMalloc1(count, &idx));
52: count = 0;
53: for (i = down; i < up; i++) {
54: for (j = bottom; j < top; j++) {
55: for (k = 0; k < dd->xe - dd->xs; k++) idx[count++] = (left + j * (dd->Xe - dd->Xs)) + i * (dd->Xe - dd->Xs) * (dd->Ye - dd->Ys) + k;
56: }
57: }
58: } else SETERRQ(PetscObjectComm((PetscObject)da), PETSC_ERR_ARG_CORRUPT, "DMDA has invalid dimension %" PetscInt_FMT, dim);
60: PetscCall(VecScatterRemap(dd->ltol, idx, NULL));
61: PetscCall(PetscFree(idx));
62: PetscFunctionReturn(PETSC_SUCCESS);
63: }
65: PetscErrorCode DMLocalToLocalBegin_DA(DM da, Vec g, InsertMode mode, Vec l)
66: {
67: DM_DA *dd = (DM_DA *)da->data;
69: PetscFunctionBegin;
71: if (!dd->ltol) PetscCall(DMLocalToLocalCreate_DA(da));
72: PetscCall(VecScatterBegin(dd->ltol, g, l, mode, SCATTER_FORWARD));
73: PetscFunctionReturn(PETSC_SUCCESS);
74: }
76: PetscErrorCode DMLocalToLocalEnd_DA(DM da, Vec g, InsertMode mode, Vec l)
77: {
78: DM_DA *dd = (DM_DA *)da->data;
80: PetscFunctionBegin;
83: PetscCall(VecScatterEnd(dd->ltol, g, l, mode, SCATTER_FORWARD));
84: PetscFunctionReturn(PETSC_SUCCESS);
85: }