Actual source code: ex14.c


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

  5: #include <petscvec.h>

  7: int main(int argc, char **argv)
  8: {
  9:   PetscInt    n = 5, N;
 10:   PetscMPIInt size, rank;
 11:   PetscScalar value, zero = 0.0;
 12:   Vec         x, y;
 13:   IS          is1, is2;
 14:   VecScatter  ctx = 0;

 17:   PetscInitialize(&argc, &argv, (char *)0, help);
 18:   MPI_Comm_size(PETSC_COMM_WORLD, &size);
 19:   MPI_Comm_rank(PETSC_COMM_WORLD, &rank);

 21:   /* create two vectors */
 22:   N = size * n;
 23:   VecCreate(PETSC_COMM_WORLD, &y);
 24:   VecSetSizes(y, PETSC_DECIDE, N);
 25:   VecSetFromOptions(y);
 26:   VecCreateSeq(PETSC_COMM_SELF, N, &x);

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

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

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

 41:   VecView(y, PETSC_VIEWER_STDOUT_WORLD);

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

 48:   PetscFinalize();
 49:   return 0;
 50: }

 52: /*TEST

 54:    test:
 55:       nsize: 2

 57: TEST*/