Actual source code: ex5.c
2: static char help[] = "Scatters from a parallel vector to a sequential vector.\n\
3: This does case when we are merely selecting the local part of the\n\
4: parallel vector.\n\n";
6: #include <petscvec.h>
8: int main(int argc, char **argv)
9: {
10: PetscMPIInt size, rank;
11: PetscInt n = 5, i;
12: PetscScalar value;
13: Vec x, y;
14: IS is1, is2;
15: VecScatter ctx = 0;
18: PetscInitialize(&argc, &argv, (char *)0, help);
19: MPI_Comm_size(PETSC_COMM_WORLD, &size);
20: MPI_Comm_rank(PETSC_COMM_WORLD, &rank);
22: /* create two vectors */
23: VecCreate(PETSC_COMM_WORLD, &x);
24: VecSetSizes(x, PETSC_DECIDE, size * n);
25: VecSetFromOptions(x);
26: VecCreateSeq(PETSC_COMM_SELF, n, &y);
28: /* create two index sets */
29: ISCreateStride(PETSC_COMM_SELF, n, n * rank, 1, &is1);
30: ISCreateStride(PETSC_COMM_SELF, n, 0, 1, &is2);
32: /* each processor inserts the entire vector */
33: /* this is redundant but tests assembly */
34: for (i = 0; i < n * size; i++) {
35: value = (PetscScalar)i;
36: VecSetValues(x, 1, &i, &value, INSERT_VALUES);
37: }
38: VecAssemblyBegin(x);
39: VecAssemblyEnd(x);
40: VecView(x, PETSC_VIEWER_STDOUT_WORLD);
42: VecScatterCreate(x, is1, y, is2, &ctx);
43: VecScatterBegin(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD);
44: VecScatterEnd(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD);
45: VecScatterDestroy(&ctx);
47: if (rank == 0) {
48: PetscPrintf(PETSC_COMM_SELF, "----\n");
49: VecView(y, PETSC_VIEWER_STDOUT_SELF);
50: }
52: VecDestroy(&x);
53: VecDestroy(&y);
54: ISDestroy(&is1);
55: ISDestroy(&is2);
57: PetscFinalize();
58: return 0;
59: }
61: /*TEST
63: test:
64: nsize: 2
66: TEST*/