Actual source code: ex25.c


  2: static char help[] = "Scatters from a parallel vector to a sequential vector.  In\n\
  3: this case processor zero is as long as the entire parallel vector; rest are zero length.\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:   if (rank == 0) {
 27:     PetscCall(VecCreateSeq(PETSC_COMM_SELF, N, &x));
 28:   } else {
 29:     PetscCall(VecCreateSeq(PETSC_COMM_SELF, 0, &x));
 30:   }

 32:   /* create two index sets */
 33:   if (rank == 0) {
 34:     PetscCall(ISCreateStride(PETSC_COMM_SELF, N, 0, 1, &is1));
 35:     PetscCall(ISCreateStride(PETSC_COMM_SELF, N, 0, 1, &is2));
 36:   } else {
 37:     PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &is1));
 38:     PetscCall(ISCreateStride(PETSC_COMM_SELF, 0, 0, 1, &is2));
 39:   }

 41:   PetscCall(VecSet(x, zero));
 42:   PetscCall(VecGetOwnershipRange(y, &low, &high));
 43:   for (i = 0; i < n; i++) {
 44:     iglobal = i + low;
 45:     value   = (PetscScalar)(i + 10 * rank);
 46:     PetscCall(VecSetValues(y, 1, &iglobal, &value, INSERT_VALUES));
 47:   }
 48:   PetscCall(VecAssemblyBegin(y));
 49:   PetscCall(VecAssemblyEnd(y));
 50:   PetscCall(VecView(y, PETSC_VIEWER_STDOUT_WORLD));

 52:   PetscCall(VecScatterCreate(y, is2, x, is1, &ctx));
 53:   PetscCall(VecScatterBegin(ctx, y, x, ADD_VALUES, SCATTER_FORWARD));
 54:   PetscCall(VecScatterEnd(ctx, y, x, ADD_VALUES, SCATTER_FORWARD));
 55:   PetscCall(VecScatterDestroy(&ctx));

 57:   if (rank == 0) {
 58:     PetscCall(PetscPrintf(PETSC_COMM_SELF, "----\n"));
 59:     PetscCall(VecView(x, PETSC_VIEWER_STDOUT_SELF));
 60:   }

 62:   PetscCall(VecDestroy(&x));
 63:   PetscCall(VecDestroy(&y));
 64:   PetscCall(ISDestroy(&is1));
 65:   PetscCall(ISDestroy(&is2));

 67:   PetscCall(PetscFinalize());
 68:   return 0;
 69: }

 71: /*TEST

 73:    test:
 74:       nsize: 3

 76: TEST*/