xref: /petsc/src/vec/vec/tests/ex23.c (revision 882062120559b0524b4d5919e1e1b8157c6a8559)
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   PetscFunctionBeginUser;
25   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
26   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
27   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
28 
29   PetscCheck(size == 2,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run with 2 processors");
30 
31   /* create two vectors */
32   if (rank == 0) nlocal = 8;
33   else nlocal = 4;
34   PetscCall(VecCreate(PETSC_COMM_WORLD,&x));
35   PetscCall(VecSetSizes(x,nlocal,12));
36   PetscCall(VecSetFromOptions(x));
37   PetscCall(VecCreate(PETSC_COMM_SELF,&y));
38   PetscCall(VecSetSizes(y,8,PETSC_DECIDE));
39   PetscCall(VecSetFromOptions(y));
40 
41   /* create two index sets */
42   if (rank == 0) {
43     blocks[0] = 0; blocks[1] = 2;
44   } else {
45     blocks[0] = 1; blocks[1] = 2;
46   }
47   PetscCall(ISCreateBlock(PETSC_COMM_SELF,4,2,blocks,PETSC_COPY_VALUES,&is1));
48   PetscCall(ISCreateStride(PETSC_COMM_SELF,8,0,1,&is2));
49 
50   for (i=0; i<12; i++) {
51     value = i;
52     PetscCall(VecSetValues(x,1,&i,&value,INSERT_VALUES));
53   }
54   PetscCall(VecAssemblyBegin(x));
55   PetscCall(VecAssemblyEnd(x));
56 
57   PetscCall(VecScatterCreate(x,is1,y,is2,&ctx));
58   PetscCall(VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD));
59   PetscCall(VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD));
60   PetscCall(VecScatterDestroy(&ctx));
61 
62   PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&subviewer));
63   PetscCall(VecView(y,subviewer));
64   PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&subviewer));
65 
66   PetscCall(VecDestroy(&x));
67   PetscCall(VecDestroy(&y));
68   PetscCall(ISDestroy(&is1));
69   PetscCall(ISDestroy(&is2));
70 
71   PetscCall(PetscFinalize());
72   return 0;
73 }
74 
75 /*TEST
76 
77    testset:
78       nsize: 2
79       output_file: output/ex23_1.out
80       filter: grep -v "  type:"
81       diff_args: -j
82       test:
83         suffix: standard
84         args: -vec_type standard
85       test:
86         requires: cuda
87         suffix: cuda
88         args: -vec_type cuda
89       test:
90         requires: viennacl
91         suffix:  viennacl
92         args: -vec_type viennacl
93       test:
94         requires: !sycl kokkos_kernels
95         suffix: kokkos
96         args: -vec_type kokkos
97       test:
98         requires: hip
99         suffix: hip
100         args: -vec_type hip
101 
102 TEST*/
103