Actual source code: ex13.c
1: static char help[] = "Demonstrates scattering with the indices specified by a process that is not sender or receiver.\n\n";
3: #include <petscvec.h>
5: int main(int argc, char **argv)
6: {
7: PetscMPIInt rank, size;
8: Vec x, y;
9: IS is1, is2;
10: PetscInt n, N, ix[2], iy[2];
11: VecScatter ctx;
13: PetscFunctionBeginUser;
14: PetscCall(PetscInitialize(&argc, &argv, NULL, help));
15: PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
16: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
17: PetscCheck(size >= 3, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This example needs at least 3 processes");
19: /* create two vectors */
20: n = 2;
21: N = 2 * size;
23: PetscCall(VecCreateFromOptions(PETSC_COMM_WORLD, NULL, 1, n, N, &x));
24: PetscCall(VecDuplicate(x, &y));
26: /* Specify indices to send from the next process in the ring */
27: ix[0] = ((rank + 1) * n + 0) % N;
28: ix[1] = ((rank + 1) * n + 1) % N;
29: /* And put them on the process after that in the ring */
30: iy[0] = ((rank + 2) * n + 0) % N;
31: iy[1] = ((rank + 2) * n + 1) % N;
33: /* create two index sets */
34: PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n, ix, PETSC_USE_POINTER, &is1));
35: PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, n, iy, PETSC_USE_POINTER, &is2));
37: PetscCall(VecSetValue(x, rank * n, rank * n, INSERT_VALUES));
38: PetscCall(VecSetValue(x, rank * n + 1, rank * n + 1, INSERT_VALUES));
39: PetscCall(VecAssemblyBegin(x));
40: PetscCall(VecAssemblyEnd(x));
42: PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
43: PetscCall(PetscPrintf(PETSC_COMM_WORLD, "----\n"));
45: PetscCall(VecScatterCreate(x, is1, y, is2, &ctx));
46: PetscCall(VecScatterBegin(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD));
47: PetscCall(VecScatterEnd(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD));
48: PetscCall(VecScatterDestroy(&ctx));
50: PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
52: PetscCall(ISDestroy(&is1));
53: PetscCall(ISDestroy(&is2));
54: PetscCall(VecDestroy(&x));
55: PetscCall(VecDestroy(&y));
56: PetscCall(PetscFinalize());
57: return 0;
58: }
60: /*TEST
62: test:
63: nsize: 4
65: TEST*/