Actual source code: ex9.c

  1: static char help[] = "Scatters from a parallel vector to a sequential vector.\n\n";

  3: #include <petscvec.h>

  5: int main(int argc, char **argv)
  6: {
  7:   PetscInt    n = 5, i, idx2[3] = {0, 2, 3}, idx1[3] = {0, 1, 2};
  8:   PetscMPIInt size, rank;
  9:   PetscScalar value;
 10:   Vec         x, y;
 11:   IS          is1, is2;
 12:   VecScatter  ctx = 0;

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

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

 25:   /* create two index sets */
 26:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 3, idx1, PETSC_COPY_VALUES, &is1));
 27:   PetscCall(ISCreateGeneral(PETSC_COMM_SELF, 3, idx2, PETSC_COPY_VALUES, &is2));

 29:   /* fill local part of parallel vector */
 30:   for (i = n * rank; i < n * (rank + 1); i++) {
 31:     value = (PetscScalar)i;
 32:     PetscCall(VecSetValues(x, 1, &i, &value, INSERT_VALUES));
 33:   }
 34:   PetscCall(VecAssemblyBegin(x));
 35:   PetscCall(VecAssemblyEnd(x));

 37:   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));

 39:   PetscCall(VecSet(y, -1.0));

 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, "scattered vector\n"));
 48:     PetscCall(VecView(y, PETSC_VIEWER_STDOUT_SELF));
 49:   }
 50:   PetscCall(ISDestroy(&is1));
 51:   PetscCall(ISDestroy(&is2));
 52:   PetscCall(VecDestroy(&x));
 53:   PetscCall(VecDestroy(&y));

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

 59: /*TEST

 61:    test:
 62:       nsize: 2

 64: TEST*/