xref: /petsc/src/vec/vec/tests/ex22.c (revision 11486bccf1aeea1ca5536228f99d437b39bdaca6)
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   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;
14 
15   PetscFunctionBeginUser;
16   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
17   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
18   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
19 
20   /* create two vectors */
21   N    = size*n;
22   PetscCall(VecCreate(PETSC_COMM_WORLD,&y));
23   PetscCall(VecSetSizes(y,n,PETSC_DECIDE));
24   PetscCall(VecSetFromOptions(y));
25 
26   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
27   PetscCall(VecSetSizes(x,n,PETSC_DECIDE));
28   PetscCall(VecSetFromOptions(x));
29 
30   /* create two index sets */
31   PetscCall(ISCreateStride(PETSC_COMM_WORLD,n,n*rank,1,&is1));
32   PetscCall(ISCreateStride(PETSC_COMM_WORLD,n,(n*(rank+1))%N,1,&is2));
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     PetscCall(VecSetValues(x,1,&i,&value,INSERT_VALUES));
38   }
39   PetscCall(VecAssemblyBegin(x));
40   PetscCall(VecAssemblyEnd(x));
41 
42   PetscCall(VecSet(y,zero));
43 
44   PetscCall(VecScatterCreate(x,is1,y,is2,&ctx));
45   PetscCall(VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD));
46   PetscCall(VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD));
47   PetscCall(VecScatterDestroy(&ctx));
48 
49   PetscCall(VecView(y,PETSC_VIEWER_STDOUT_WORLD));
50 
51   PetscCall(VecDestroy(&x));
52   PetscCall(VecDestroy(&y));
53   PetscCall(ISDestroy(&is1));
54   PetscCall(ISDestroy(&is2));
55 
56   PetscCall(PetscFinalize());
57   return 0;
58 }
59 
60 /*TEST
61 
62    testset:
63       nsize: 4
64       output_file: output/ex22_1.out
65       filter: grep -v "  type:"
66       diff_args: -j
67       test:
68         suffix: standard
69         args: -vec_type standard
70       test:
71         requires: cuda
72         suffix: cuda
73         args: -vec_type cuda
74       test:
75         requires: viennacl
76         suffix:  viennacl
77         args: -vec_type viennacl
78       test:
79         requires: !sycl kokkos_kernels
80         suffix: kokkos
81         args: -vec_type kokkos
82       test:
83         requires: hip
84         suffix: hip
85         args: -vec_type hip
86 
87 TEST*/
88