Actual source code: ex22.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:   VecScatterBegin(ctx, x, y, ADD_VALUES, SCATTER_FORWARD);
 44:   VecScatterEnd(ctx, x, y, ADD_VALUES, SCATTER_FORWARD);
 45:   VecScatterDestroy(&ctx);

 47:   VecView(y, PETSC_VIEWER_STDOUT_WORLD);

 49:   VecDestroy(&x);
 50:   VecDestroy(&y);
 51:   ISDestroy(&is1);
 52:   ISDestroy(&is2);

 54:   PetscFinalize();
 55:   return 0;
 56: }

 58: /*TEST

 60:    testset:
 61:       nsize: 4
 62:       output_file: output/ex22_1.out
 63:       filter: grep -v "  type:"
 64:       diff_args: -j
 65:       test:
 66:         suffix: standard
 67:         args: -vec_type standard
 68:       test:
 69:         requires: cuda
 70:         suffix: cuda
 71:         args: -vec_type cuda
 72:       test:
 73:         requires: viennacl
 74:         suffix:  viennacl
 75:         args: -vec_type viennacl
 76:       test:
 77:         requires: kokkos_kernels
 78:         suffix: kokkos
 79:         args: -vec_type kokkos
 80:       test:
 81:         requires: hip
 82:         suffix: hip
 83:         args: -vec_type hip

 85: TEST*/