Actual source code: ex11.c
2: static char help[] = "Scatters between parallel vectors. \n\
3: uses block index sets\n\n";
5: #include <petscvec.h>
7: int main(int argc, char **argv)
8: {
9: PetscInt bs = 1, n = 5, N, i, low;
10: PetscInt ix0[3] = {5, 7, 9}, iy0[3] = {1, 2, 4}, ix1[3] = {2, 3, 1}, iy1[3] = {0, 3, 9};
11: PetscMPIInt size, rank;
12: PetscScalar *array;
13: Vec x, x1, y;
14: IS isx, isy;
15: VecScatter ctx;
16: VecScatterType type;
17: PetscBool flg;
19: PetscFunctionBeginUser;
20: PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
21: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
22: PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
24: PetscCheck(size >= 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Must run more than one processor");
26: PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL));
27: n = bs * n;
29: /* Create vector x over shared memory */
30: PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
31: PetscCall(VecSetSizes(x, n, PETSC_DECIDE));
32: PetscCall(VecSetFromOptions(x));
34: PetscCall(VecGetOwnershipRange(x, &low, NULL));
35: PetscCall(VecGetArray(x, &array));
36: for (i = 0; i < n; i++) array[i] = (PetscScalar)(i + low);
37: PetscCall(VecRestoreArray(x, &array));
38: /* PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); */
40: /* Test some vector functions */
41: PetscCall(VecAssemblyBegin(x));
42: PetscCall(VecAssemblyEnd(x));
44: PetscCall(VecGetSize(x, &N));
45: PetscCall(VecGetLocalSize(x, &n));
47: PetscCall(VecDuplicate(x, &x1));
48: PetscCall(VecCopy(x, x1));
49: PetscCall(VecEqual(x, x1, &flg));
50: PetscCheck(flg, PetscObjectComm((PetscObject)x), PETSC_ERR_ARG_WRONG, "x1 != x");
52: PetscCall(VecScale(x1, 2.0));
53: PetscCall(VecSet(x1, 10.0));
54: /* PetscCall(VecView(x1,PETSC_VIEWER_STDOUT_WORLD)); */
56: /* Create vector y over shared memory */
57: PetscCall(VecCreate(PETSC_COMM_WORLD, &y));
58: PetscCall(VecSetSizes(y, n, PETSC_DECIDE));
59: PetscCall(VecSetFromOptions(y));
60: PetscCall(VecGetArray(y, &array));
61: for (i = 0; i < n; i++) array[i] = -(PetscScalar)(i + 100 * rank);
62: PetscCall(VecRestoreArray(y, &array));
63: PetscCall(VecAssemblyBegin(y));
64: PetscCall(VecAssemblyEnd(y));
65: /* PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD)); */
67: /* Create two index sets */
68: if (rank == 0) {
69: PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, ix0, PETSC_COPY_VALUES, &isx));
70: PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, iy0, PETSC_COPY_VALUES, &isy));
71: } else {
72: PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, ix1, PETSC_COPY_VALUES, &isx));
73: PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, iy1, PETSC_COPY_VALUES, &isy));
74: }
76: if (rank == 10) {
77: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n[%d] isx:\n", rank));
78: PetscCall(ISView(isx, PETSC_VIEWER_STDOUT_SELF));
79: PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n[%d] isy:\n", rank));
80: PetscCall(ISView(isy, PETSC_VIEWER_STDOUT_SELF));
81: }
83: /* Create Vector scatter */
84: PetscCall(VecScatterCreate(x, isx, y, isy, &ctx));
85: PetscCall(VecScatterSetFromOptions(ctx));
86: PetscCall(VecScatterGetType(ctx, &type));
87: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "scatter type %s\n", type));
89: /* Test forward vecscatter */
90: PetscCall(VecSet(y, 0.0));
91: PetscCall(VecScatterBegin(ctx, x, y, ADD_VALUES, SCATTER_FORWARD));
92: PetscCall(VecScatterEnd(ctx, x, y, ADD_VALUES, SCATTER_FORWARD));
93: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nSCATTER_FORWARD y:\n"));
94: PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
96: /* Test reverse vecscatter */
97: PetscCall(VecSet(x, 0.0));
98: PetscCall(VecSet(y, 0.0));
99: PetscCall(VecGetOwnershipRange(y, &low, NULL));
100: PetscCall(VecGetArray(y, &array));
101: for (i = 0; i < n; i++) array[i] = (PetscScalar)(i + low);
102: PetscCall(VecRestoreArray(y, &array));
103: PetscCall(VecScatterBegin(ctx, y, x, ADD_VALUES, SCATTER_REVERSE));
104: PetscCall(VecScatterEnd(ctx, y, x, ADD_VALUES, SCATTER_REVERSE));
105: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\nSCATTER_REVERSE x:\n"));
106: PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
108: /* Free objects */
109: PetscCall(VecScatterDestroy(&ctx));
110: PetscCall(ISDestroy(&isx));
111: PetscCall(ISDestroy(&isy));
112: PetscCall(VecDestroy(&x));
113: PetscCall(VecDestroy(&x1));
114: PetscCall(VecDestroy(&y));
115: PetscCall(PetscFinalize());
116: return 0;
117: }
119: /*TEST
121: test:
122: nsize: 2
124: TEST*/