Actual source code: ex17.c


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

  5: #include <petscvec.h>

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

 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:   /* create two vectors */
 22:   N = size * n;
 23:   PetscCall(VecCreate(PETSC_COMM_WORLD, &y));
 24:   PetscCall(VecSetSizes(y, PETSC_DECIDE, N));
 25:   PetscCall(VecSetFromOptions(y));
 26:   PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &x));

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

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

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

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

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

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

 62: /*TEST

 64:    test:
 65:       nsize: 2

 67: TEST*/