Actual source code: ex42.c
2: static char help[] = "Scatters from a parallel vector to a parallel vector.\n\n";
4: #include <petscvec.h>
6: int main(int argc, char **argv)
7: {
8: PetscInt n = 5, N, i;
9: PetscMPIInt size, rank;
10: PetscScalar value, zero = 0.0;
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: N = size * n;
22: PetscCall(VecCreate(PETSC_COMM_WORLD, &y));
23: PetscCall(VecSetSizes(y, n, PETSC_DECIDE));
24: PetscCall(VecSetFromOptions(y));
26: PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
27: PetscCall(VecSetSizes(x, n, PETSC_DECIDE));
28: PetscCall(VecSetFromOptions(x));
30: /* create two index sets */
31: PetscCall(ISCreateStride(PETSC_COMM_WORLD, n, n * rank, 1, &is1));
32: PetscCall(ISCreateStride(PETSC_COMM_WORLD, n, (n * (rank + 1)) % N, 1, &is2));
34: /* fill local part of parallel vector x */
35: value = (PetscScalar)(rank + 1);
36: for (i = n * rank; i < n * (rank + 1); i++) PetscCall(VecSetValues(x, 1, &i, &value, INSERT_VALUES));
37: PetscCall(VecAssemblyBegin(x));
38: PetscCall(VecAssemblyEnd(x));
40: PetscCall(VecSet(y, zero));
42: PetscCall(VecScatterCreate(x, is1, y, is2, &ctx));
43: for (i = 0; i < 100; i++) {
44: PetscReal ynorm;
45: PetscInt j;
46: PetscCall(VecNormBegin(y, NORM_2, &ynorm));
47: PetscCall(PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)y)));
48: for (j = 0; j < 3; j++) {
49: PetscCall(VecScatterBegin(ctx, x, y, ADD_VALUES, SCATTER_FORWARD));
50: PetscCall(VecScatterEnd(ctx, x, y, ADD_VALUES, SCATTER_FORWARD));
51: }
52: PetscCall(VecNormEnd(y, NORM_2, &ynorm));
53: /* PetscCall(PetscPrintf(PETSC_COMM_WORLD,"ynorm = %8.2G\n",ynorm)); */
54: }
55: PetscCall(VecScatterDestroy(&ctx));
56: PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));
58: PetscCall(VecDestroy(&x));
59: PetscCall(VecDestroy(&y));
60: PetscCall(ISDestroy(&is1));
61: PetscCall(ISDestroy(&is2));
63: PetscCall(PetscFinalize());
64: return 0;
65: }