Actual source code: ex10.c


  2: static char help[] = "Scatters from a parallel vector to a sequential vector.\n\
  3: uses block index sets\n\n";

  5: #include <petscvec.h>

  7: int main(int argc, char **argv)
  8: {
  9:   PetscInt    bs = 1, n = 5, ix0[3] = {5, 7, 9}, ix1[3] = {2, 3, 4}, i, iy0[3] = {1, 2, 4}, iy1[3] = {0, 1, 3};
 10:   PetscMPIInt size, rank;
 11:   PetscScalar value;
 12:   Vec         x, y;
 13:   IS          isx, isy;
 14:   VecScatter  ctx = 0, newctx;

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

 21:   PetscCheck(size == 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Must run with 2 processors");

 23:   PetscCall(PetscOptionsGetInt(NULL, NULL, "-bs", &bs, NULL));
 24:   n = bs * n;

 26:   /* create two vectors */
 27:   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
 28:   PetscCall(VecSetSizes(x, PETSC_DECIDE, size * n));
 29:   PetscCall(VecSetFromOptions(x));
 30:   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &y));

 32:   /* create two index sets */
 33:   if (rank == 0) {
 34:     PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, ix0, PETSC_COPY_VALUES, &isx));
 35:     PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, iy0, PETSC_COPY_VALUES, &isy));
 36:   } else {
 37:     PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, ix1, PETSC_COPY_VALUES, &isx));
 38:     PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, iy1, PETSC_COPY_VALUES, &isy));
 39:   }

 41:   /* fill local part of parallel vector */
 42:   for (i = n * rank; i < n * (rank + 1); i++) {
 43:     value = (PetscScalar)i;
 44:     PetscCall(VecSetValues(x, 1, &i, &value, INSERT_VALUES));
 45:   }
 46:   PetscCall(VecAssemblyBegin(x));
 47:   PetscCall(VecAssemblyEnd(x));

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

 51:   /* fill local part of parallel vector */
 52:   for (i = 0; i < n; i++) {
 53:     value = -(PetscScalar)(i + 100 * rank);
 54:     PetscCall(VecSetValues(y, 1, &i, &value, INSERT_VALUES));
 55:   }
 56:   PetscCall(VecAssemblyBegin(y));
 57:   PetscCall(VecAssemblyEnd(y));

 59:   PetscCall(VecScatterCreate(x, isx, y, isy, &ctx));
 60:   PetscCall(VecScatterCopy(ctx, &newctx));
 61:   PetscCall(VecScatterDestroy(&ctx));

 63:   PetscCall(VecScatterBegin(newctx, y, x, INSERT_VALUES, SCATTER_REVERSE));
 64:   PetscCall(VecScatterEnd(newctx, y, x, INSERT_VALUES, SCATTER_REVERSE));
 65:   PetscCall(VecScatterDestroy(&newctx));

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

 69:   PetscCall(ISDestroy(&isx));
 70:   PetscCall(ISDestroy(&isy));
 71:   PetscCall(VecDestroy(&x));
 72:   PetscCall(VecDestroy(&y));

 74:   PetscCall(PetscFinalize());
 75:   return 0;
 76: }

 78: /*TEST

 80:    test:
 81:       nsize: 2

 83: TEST*/