Actual source code: ex5.c

  1: static char help[] = "Scatters from a parallel vector to a sequential vector.\n\
  2: This does case when we are merely selecting the local part of the\n\
  3: parallel vector.\n\n";

  5: #include <petscvec.h>

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

 16:   PetscFunctionBeginUser;
 17:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
 18:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
 19:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));

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

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

 31:   /* each processor inserts the entire vector */
 32:   /* this is redundant but tests assembly */
 33:   for (i = 0; i < n * size; i++) {
 34:     value = (PetscScalar)i;
 35:     PetscCall(VecSetValues(x, 1, &i, &value, INSERT_VALUES));
 36:   }
 37:   PetscCall(VecAssemblyBegin(x));
 38:   PetscCall(VecAssemblyEnd(x));
 39:   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));

 41:   PetscCall(VecScatterCreate(x, is1, y, is2, &ctx));
 42:   PetscCall(VecScatterBegin(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD));
 43:   PetscCall(VecScatterEnd(ctx, x, y, INSERT_VALUES, SCATTER_FORWARD));
 44:   PetscCall(VecScatterDestroy(&ctx));

 46:   if (rank == 0) {
 47:     PetscCall(PetscPrintf(PETSC_COMM_SELF, "----\n"));
 48:     PetscCall(VecView(y, PETSC_VIEWER_STDOUT_SELF));
 49:   }

 51:   PetscCall(VecDestroy(&x));
 52:   PetscCall(VecDestroy(&y));
 53:   PetscCall(ISDestroy(&is1));
 54:   PetscCall(ISDestroy(&is2));

 56:   PetscCall(PetscFinalize());
 57:   return 0;
 58: }

 60: /*TEST

 62:      test:
 63:        nsize: 2

 65: TEST*/