Actual source code: ex12.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, i, low;
 10:   PetscInt     ix0[3] = {5, 7, 9}, iy0[3] = {1, 2, 4}, ix1[3] = {2, 3, 4}, iy1[3] = {0, 1, 3};
 11:   PetscMPIInt  size, rank;
 12:   PetscScalar *array;
 13:   Vec          x, y;
 14:   IS           isx, isy;
 15:   VecScatter   ctx;

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

 22:   PetscCheck(size >= 2, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Must run more than one processor");

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

 27:   /* Create vector x over shared memory */
 28:   PetscCall(VecCreate(PETSC_COMM_WORLD, &x));
 29:   PetscCall(VecSetSizes(x, n, PETSC_DECIDE));
 30:   PetscCall(VecSetFromOptions(x));

 32:   PetscCall(VecGetOwnershipRange(x, &low, NULL));
 33:   PetscCall(VecGetArray(x, &array));
 34:   for (i = 0; i < n; i++) array[i] = (PetscScalar)(i + low);
 35:   PetscCall(VecRestoreArray(x, &array));

 37:   /* Create a sequential vector y */
 38:   PetscCall(VecCreateSeq(PETSC_COMM_SELF, n, &y));
 39:   PetscCall(VecSet(y, 0.0));

 41:   /* Create two index sets */
 42:   if (rank == 0) {
 43:     PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, ix0, PETSC_COPY_VALUES, &isx));
 44:     PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, iy0, PETSC_COPY_VALUES, &isy));
 45:   } else {
 46:     PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, ix1, PETSC_COPY_VALUES, &isx));
 47:     PetscCall(ISCreateBlock(PETSC_COMM_SELF, bs, 3, iy1, PETSC_COPY_VALUES, &isy));
 48:   }

 50:   if (rank == 10) {
 51:     PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n[%d] isx:\n", rank));
 52:     PetscCall(ISView(isx, PETSC_VIEWER_STDOUT_SELF));
 53:   }

 55:   PetscCall(VecScatterCreate(x, isx, y, isy, &ctx));
 56:   PetscCall(VecScatterSetFromOptions(ctx));

 58:   /* Test forward vecscatter */
 59:   PetscCall(VecScatterBegin(ctx, x, y, ADD_VALUES, SCATTER_FORWARD));
 60:   PetscCall(VecScatterEnd(ctx, x, y, ADD_VALUES, SCATTER_FORWARD));
 61:   if (rank == 0) {
 62:     PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] y:\n", rank));
 63:     PetscCall(VecView(y, PETSC_VIEWER_STDOUT_SELF));
 64:   }

 66:   /* Test reverse vecscatter */
 67:   PetscCall(VecScale(y, -1.0));
 68:   if (rank) PetscCall(VecScale(y, 1.0 / (size - 1)));

 70:   PetscCall(VecScatterBegin(ctx, y, x, ADD_VALUES, SCATTER_REVERSE));
 71:   PetscCall(VecScatterEnd(ctx, y, x, ADD_VALUES, SCATTER_REVERSE));
 72:   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));

 74:   /* Free spaces */
 75:   PetscCall(VecScatterDestroy(&ctx));
 76:   PetscCall(ISDestroy(&isx));
 77:   PetscCall(ISDestroy(&isy));
 78:   PetscCall(VecDestroy(&x));
 79:   PetscCall(VecDestroy(&y));
 80:   PetscCall(PetscFinalize());
 81:   return 0;
 82: }

 84: /*TEST

 86:    test:
 87:       nsize: 3

 89: TEST*/