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