xref: /petsc/src/vec/vec/tests/ex22.c (revision efbe7e8a80d07327753dbe0b33efee01e046af3f)
1 
2 static char help[] = "Scatters from a parallel vector to a parallel vector.\n\n";
3 
4 #include <petscvec.h>
5 
6 int main(int argc,char **argv)
7 {
8   PetscErrorCode ierr;
9   PetscInt       n = 5,N,i;
10   PetscMPIInt    size,rank;
11   PetscScalar    value,zero = 0.0;
12   Vec            x,y;
13   IS             is1,is2;
14   VecScatter     ctx = 0;
15 
16   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
17   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
18   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
19 
20   /* create two vectors */
21   N    = size*n;
22   ierr = VecCreate(PETSC_COMM_WORLD,&y);CHKERRQ(ierr);
23   ierr = VecSetSizes(y,n,PETSC_DECIDE);CHKERRQ(ierr);
24   ierr = VecSetFromOptions(y);CHKERRQ(ierr);
25 
26   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
27   ierr = VecSetSizes(x,n,PETSC_DECIDE);CHKERRQ(ierr);
28   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
29 
30   /* create two index sets */
31   ierr = ISCreateStride(PETSC_COMM_WORLD,n,n*rank,1,&is1);CHKERRQ(ierr);
32   ierr = ISCreateStride(PETSC_COMM_WORLD,n,(n*(rank+1))%N,1,&is2);CHKERRQ(ierr);
33 
34   /* fill local part of parallel vector x */
35   value = (PetscScalar)(rank+1);
36   for (i=n*rank; i<n*(rank+1); i++) {
37     ierr = VecSetValues(x,1,&i,&value,INSERT_VALUES);CHKERRQ(ierr);
38   }
39   ierr = VecAssemblyBegin(x);CHKERRQ(ierr);
40   ierr = VecAssemblyEnd(x);CHKERRQ(ierr);
41 
42   ierr = VecSet(y,zero);CHKERRQ(ierr);
43 
44   ierr = VecScatterCreate(x,is1,y,is2,&ctx);CHKERRQ(ierr);
45   ierr = VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
46   ierr = VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
47   ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr);
48 
49   ierr = VecView(y,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
50 
51   ierr = VecDestroy(&x);CHKERRQ(ierr);
52   ierr = VecDestroy(&y);CHKERRQ(ierr);
53   ierr = ISDestroy(&is1);CHKERRQ(ierr);
54   ierr = ISDestroy(&is2);CHKERRQ(ierr);
55 
56   ierr = PetscFinalize();
57   return ierr;
58 }
59 
60 
61 
62 /*TEST
63 
64    testset:
65       nsize: 4
66       output_file: output/ex22_1.out
67       filter: grep -v "  type:"
68       test:
69         suffix: standard
70         args: -vec_type standard
71       test:
72         requires: cuda
73         suffix: cuda
74         args: -vec_type cuda
75       test:
76         requires: viennacl
77         suffix:  viennacl
78         args: -vec_type viennacl
79       test:
80         requires: kokkos_kernels
81         suffix: kokkos
82         args: -vec_type kokkos
83 
84 TEST*/
85