Actual source code: ex17.c

  1: static char help[] = "Scatters from a parallel vector to a sequential vector.  In\n\
  2: this case each local vector is as long as the entire parallel vector.\n\n";

  4: #include <petscvec.h>

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

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

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

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

 31:   PetscCall(VecSet(x, zero));
 32:   PetscCall(VecGetOwnershipRange(y, &low, &high));
 33:   for (i = 0; i < n; i++) {
 34:     iglobal = i + low;
 35:     value   = (PetscScalar)(i + 10 * rank);
 36:     PetscCall(VecSetValues(y, 1, &iglobal, &value, INSERT_VALUES));
 37:   }
 38:   PetscCall(VecAssemblyBegin(y));
 39:   PetscCall(VecAssemblyEnd(y));
 40:   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));

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

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

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

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

 61: /*TEST

 63:    test:
 64:       nsize: 2

 66: TEST*/