Actual source code: ex25.c
1: static char help[] = "Scatters from a parallel vector to a sequential vector. In\n\
2: this case processor zero is as long as the entire parallel vector; rest are zero length.\n\n";
4: #include <petscvec.h>
6: int main(int argc, char **argv)
7: {
8: PetscInt n = 5, N, low, high, iglobal, i;
9: PetscMPIInt size, rank;
10: PetscScalar value, zero = 0.0;
11: Vec x, y;
12: IS is1, is2;
13: VecScatter ctx;
15: PetscFunctionBeginUser;
16: PetscCall(PetscInitialize(&argc, &argv, NULL, help));
17: PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
18: PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
20: /* create two vectors */
21: N = size * n;
22: PetscCall(VecCreate(PETSC_COMM_WORLD, &y));
23: PetscCall(VecSetSizes(y, PETSC_DECIDE, N));
24: PetscCall(VecSetFromOptions(y));
25: if (rank == 0) {
26: PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &x));
27: } else {
28: PetscCall(VecCreateSeq(PETSC_COMM_SELF, 0, &x));
29: }
31: /* create two index sets */
32: if (rank == 0) {
33: PetscCall(ISCreateStride(PETSC_COMM_SELF, N, 0, 1, &is1));
34: PetscCall(ISCreateStride(PETSC_COMM_SELF, N, 0, 1, &is2));
35: } else {
36: PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &is1));
37: PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &is2));
38: }
40: PetscCall(VecSet(x, zero));
41: PetscCall(VecGetOwnershipRange(y, &low, &high));
42: for (i = 0; i < n; i++) {
43: iglobal = i + low;
44: value = (PetscScalar)(i + 10 * rank);
45: PetscCall(VecSetValues(y, 1, &iglobal, &value, INSERT_VALUES));
46: }
47: PetscCall(VecAssemblyBegin(y));
48: PetscCall(VecAssemblyEnd(y));
49: PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
51: PetscCall(VecScatterCreate(y, is2, x, is1, &ctx));
52: PetscCall(VecScatterBegin(ctx, y, x, ADD_VALUES, SCATTER_FORWARD));
53: PetscCall(VecScatterEnd(ctx, y, x, ADD_VALUES, SCATTER_FORWARD));
54: PetscCall(VecScatterDestroy(&ctx));
56: if (rank == 0) {
57: PetscCall(PetscPrintf(PETSC_COMM_SELF, "----\n"));
58: PetscCall(VecView(x, PETSC_VIEWER_STDOUT_SELF));
59: }
61: PetscCall(VecDestroy(&x));
62: PetscCall(VecDestroy(&y));
63: PetscCall(ISDestroy(&is1));
64: PetscCall(ISDestroy(&is2));
66: PetscCall(PetscFinalize());
67: return 0;
68: }
70: /*TEST
72: test:
73: nsize: 3
75: TEST*/