xref: /petsc/src/vec/vec/tests/ex23.c (revision cac4c232dc4f93991e342196e27ef7b0655dac7b)
1 
2 static char help[] = "Scatters from a parallel vector to a sequential vector.\n\
3   Using a blocked send and a strided receive.\n\n";
4 
5 /*
6         0 1 2 3 | 4 5 6 7 ||  8 9 10 11
7 
8      Scatter first and third block to first processor and
9      second and third block to second processor
10 */
11 
12 #include <petscvec.h>
13 
14 int main(int argc,char **argv)
15 {
16   PetscInt       i,blocks[2],nlocal;
17   PetscMPIInt    size,rank;
18   PetscScalar    value;
19   Vec            x,y;
20   IS             is1,is2;
21   VecScatter     ctx = 0;
22   PetscViewer    subviewer;
23 
24   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
25   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
26   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
27 
28   PetscCheck(size == 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run with 2 processors");
29 
30   /* create two vectors */
31   if (rank == 0) nlocal = 8;
32   else nlocal = 4;
33   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
34   PetscCall(VecSetSizes(x,nlocal,12));
35   PetscCall(VecSetFromOptions(x));
36   PetscCall(VecCreate(PETSC_COMM_SELF,&y));
37   PetscCall(VecSetSizes(y,8,PETSC_DECIDE));
38   PetscCall(VecSetFromOptions(y));
39 
40   /* create two index sets */
41   if (rank == 0) {
42     blocks[0] = 0; blocks[1] = 2;
43   } else {
44     blocks[0] = 1; blocks[1] = 2;
45   }
46   PetscCall(ISCreateBlock(PETSC_COMM_SELF,4,2,blocks,PETSC_COPY_VALUES,&is1));
47   PetscCall(ISCreateStride(PETSC_COMM_SELF,8,0,1,&is2));
48 
49   for (i=0; i<12; i++) {
50     value = i;
51     PetscCall(VecSetValues(x,1,&i,&value,INSERT_VALUES));
52   }
53   PetscCall(VecAssemblyBegin(x));
54   PetscCall(VecAssemblyEnd(x));
55 
56   PetscCall(VecScatterCreate(x,is1,y,is2,&ctx));
57   PetscCall(VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD));
58   PetscCall(VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD));
59   PetscCall(VecScatterDestroy(&ctx));
60 
61   PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&subviewer));
62   PetscCall(VecView(y,subviewer));
63   PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&subviewer));
64 
65   PetscCall(VecDestroy(&x));
66   PetscCall(VecDestroy(&y));
67   PetscCall(ISDestroy(&is1));
68   PetscCall(ISDestroy(&is2));
69 
70   PetscCall(PetscFinalize());
71   return 0;
72 }
73 
74 /*TEST
75 
76    testset:
77       nsize: 2
78       output_file: output/ex23_1.out
79       filter: grep -v "  type:"
80       diff_args: -j
81       test:
82         suffix: standard
83         args: -vec_type standard
84       test:
85         requires: cuda
86         suffix: cuda
87         args: -vec_type cuda
88       test:
89         requires: viennacl
90         suffix:  viennacl
91         args: -vec_type viennacl
92       test:
93         requires: !sycl kokkos_kernels
94         suffix: kokkos
95         args: -vec_type kokkos
96       test:
97         requires: hip
98         suffix: hip
99         args: -vec_type hip
100 
101 TEST*/
102