Actual source code: ex42.c

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

  3: #include <petscvec.h>

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

 14:   PetscFunctionBeginUser;
 15:   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
 16:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
 17:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));

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

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

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

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

 39:   PetscCall(VecSet(y, zero));

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

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

 62:   PetscCall(PetscFinalize());
 63:   return 0;
 64: }