Actual source code: ex23.c
petsc-3.10.5 2019-03-28
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";
5: /*
6: 0 1 2 3 | 4 5 6 7 || 8 9 10 11
8: Scatter first and third block to first processor and
9: second and third block to second processor
10: */
11: /*T
12: requires: x
13: T*/
15: #include <petscvec.h>
17: int main(int argc,char **argv)
18: {
20: PetscInt i,blocks[2],nlocal;
21: PetscMPIInt size,rank;
22: PetscScalar value;
23: Vec x,y;
24: IS is1,is2;
25: VecScatter ctx = 0;
26: PetscViewer subviewer;
28: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
29: MPI_Comm_size(PETSC_COMM_WORLD,&size);
30: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
32: if (size != 2) SETERRQ(PETSC_COMM_SELF,1,"Must run with 2 processors");
34: /* create two vectors */
35: if (!rank) nlocal = 8;
36: else nlocal = 4;
37: VecCreate(PETSC_COMM_WORLD,&x);
38: VecSetSizes(x,nlocal,12);
39: VecSetFromOptions(x);
40: VecCreate(PETSC_COMM_SELF,&y);
41: VecSetSizes(y,8,PETSC_DECIDE);
42: VecSetFromOptions(y);
44: /* create two index sets */
45: if (!rank) {
46: blocks[0] = 0; blocks[1] = 2;
47: } else {
48: blocks[0] = 1; blocks[1] = 2;
49: }
50: ISCreateBlock(PETSC_COMM_SELF,4,2,blocks,PETSC_COPY_VALUES,&is1);
51: ISCreateStride(PETSC_COMM_SELF,8,0,1,&is2);
53: for (i=0; i<12; i++) {
54: value = i;
55: VecSetValues(x,1,&i,&value,INSERT_VALUES);
56: }
57: VecAssemblyBegin(x);
58: VecAssemblyEnd(x);
60: VecScatterCreate(x,is1,y,is2,&ctx);
61: VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
62: VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
63: VecScatterDestroy(&ctx);
65: PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&subviewer);
66: VecView(y,subviewer);
67: PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&subviewer);
69: VecDestroy(&x);
70: VecDestroy(&y);
71: ISDestroy(&is1);
72: ISDestroy(&is2);
74: PetscFinalize();
75: return ierr;
76: }
80: /*TEST
82: testset:
83: nsize: 2
84: output_file: output/ex23_1.out
85: filter: grep -v " type:"
86: test:
87: suffix: standard
88: args: -vec_type standard
89: test:
90: requires: veccuda
91: suffix: cuda
92: args: -vec_type cuda
93: test:
94: requires: viennacl
95: suffix: viennacl
96: args: -vec_type viennacl
97: TEST*/