Actual source code: ex9.c
2: static char help[] = "Scatters from a parallel vector to a sequential vector.\n\n";
4: #include <petscvec.h>
6: int main(int argc, char **argv)
7: {
8: PetscInt n = 5, i, idx2[3] = {0, 2, 3}, idx1[3] = {0, 1, 2};
9: PetscMPIInt size, rank;
10: PetscScalar value;
11: Vec x, y;
12: IS is1, is2;
13: VecScatter ctx = 0;
15: PetscFunctionBeginUser;
16: PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
17: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
18: PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
20: /* create two vectors */
21: PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
22: PetscCall(VecSetSizes(x, PETSC_DECIDE, size * n));
23: PetscCall(VecSetFromOptions(x));
24: PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &y));
26: /* create two index sets */
27: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 3, idx1, PETSC_COPY_VALUES, &is1));
28: PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 3, idx2, PETSC_COPY_VALUES, &is2));
30: /* fill local part of parallel vector */
31: for (i = n * rank; i < n * (rank + 1); i++) {
32: value = (PetscScalar)i;
33: PetscCall(VecSetValues(x, 1, &i, &value, INSERT_VALUES));
34: }
35: PetscCall(VecAssemblyBegin(x));
36: PetscCall(VecAssemblyEnd(x));
38: PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
40: PetscCall(VecSet(y, -1.0));
42: PetscCall(VecScatterCreate(x, is1, y, is2, &ctx));
43: PetscCall(VecScatterBegin(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD));
44: PetscCall(VecScatterEnd(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD));
45: PetscCall(VecScatterDestroy(&ctx));
47: if (rank == 0) {
48: PetscCall(PetscPrintf(PETSC_COMM_SELF, "scattered vector\n"));
49: PetscCall(VecView(y, PETSC_VIEWER_STDOUT_SELF));
50: }
51: PetscCall(ISDestroy(&is1));
52: PetscCall(ISDestroy(&is2));
53: PetscCall(VecDestroy(&x));
54: PetscCall(VecDestroy(&y));
56: PetscCall(PetscFinalize());
57: return 0;
58: }
60: /*TEST
62: test:
63: nsize: 2
65: TEST*/