Actual source code: ex14.c

  1: static char help[] = "Scatters from a sequential vector to a parallel vector.\n\
  2: This does the tricky case.\n\n";

  4: #include <petscvec.h>

  6: int main(int argc, char **argv)
  7: {
  8:   PetscInt    n = 5, N;
  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, 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:   PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &x));

 27:   /* create two index sets */
 28:   PetscCall(ISCreateStride(PETSC_COMM_SELF, n, 0, 1, &is1));
 29:   PetscCall(ISCreateStride(PETSC_COMM_SELF, n, rank, 1, &is2));

 31:   value = rank + 1;
 32:   PetscCall(VecSet(x, value));
 33:   PetscCall(VecSet(y, zero));

 35:   PetscCall(VecScatterCreate(x, is1, y, is2, &ctx));
 36:   PetscCall(VecScatterBegin(ctx, x, y, ADD_VALUES, SCATTER_FORWARD));
 37:   PetscCall(VecScatterEnd(ctx, x, y, ADD_VALUES, SCATTER_FORWARD));
 38:   PetscCall(VecScatterDestroy(&ctx));

 40:   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));

 42:   PetscCall(VecDestroy(&x));
 43:   PetscCall(VecDestroy(&y));
 44:   PetscCall(ISDestroy(&is1));
 45:   PetscCall(ISDestroy(&is2));

 47:   PetscCall(PetscFinalize());
 48:   return 0;
 49: }

 51: /*TEST

 53:    test:
 54:       nsize: 2

 56: TEST*/