Actual source code: ex42.c


  2: static char help[] = "Scatters from a parallel vector to a parallel vector.\n\n";

  4: #include <petscvec.h>

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

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

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

 26:   VecCreate(PETSC_COMM_WORLD, &x);
 27:   VecSetSizes(x, n, PETSC_DECIDE);
 28:   VecSetFromOptions(x);

 30:   /* create two index sets */
 31:   ISCreateStride(PETSC_COMM_WORLD, n, n * rank, 1, &is1);
 32:   ISCreateStride(PETSC_COMM_WORLD, n, (n * (rank + 1)) % N, 1, &is2);

 34:   /* fill local part of parallel vector x */
 35:   value = (PetscScalar)(rank + 1);
 36:   for (i = n * rank; i < n * (rank + 1); i++) VecSetValues(x, 1, &i, &value, INSERT_VALUES);
 37:   VecAssemblyBegin(x);
 38:   VecAssemblyEnd(x);

 40:   VecSet(y, zero);

 42:   VecScatterCreate(x, is1, y, is2, &ctx);
 43:   for (i = 0; i < 100; i++) {
 44:     PetscReal ynorm;
 45:     PetscInt  j;
 46:     VecNormBegin(y, NORM_2, &ynorm);
 47:     PetscCommSplitReductionBegin(PetscObjectComm((PetscObject)y));
 48:     for (j = 0; j < 3; j++) {
 49:       VecScatterBegin(ctx, x, y, ADD_VALUES, SCATTER_FORWARD);
 50:       VecScatterEnd(ctx, x, y, ADD_VALUES, SCATTER_FORWARD);
 51:     }
 52:     VecNormEnd(y, NORM_2, &ynorm);
 53:     /* PetscPrintf(PETSC_COMM_WORLD,"ynorm = %8.2G\n",ynorm); */
 54:   }
 55:   VecScatterDestroy(&ctx);
 56:   VecView(y, PETSC_VIEWER_STDOUT_WORLD);

 58:   VecDestroy(&x);
 59:   VecDestroy(&y);
 60:   ISDestroy(&is1);
 61:   ISDestroy(&is2);

 63:   PetscFinalize();
 64:   return 0;
 65: }