Actual source code: ex24.c
2: static char help[] = "Scatters from a parallel vector to a sequential vector.\n\
3: Tests where the local part of the scatter is a copy.\n\n";
5: #include <petscvec.h>
7: int main(int argc,char **argv)
8: {
9: PetscMPIInt size,rank;
10: PetscInt n = 5,i,*blks,bs = 1,m = 2;
11: PetscScalar value;
12: Vec x,y;
13: IS is1,is2;
14: VecScatter ctx = 0;
15: PetscViewer sviewer;
17: PetscInitialize(&argc,&argv,(char*)0,help);
19: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
20: PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);
22: MPI_Comm_size(PETSC_COMM_WORLD,&size);
23: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
25: /* create two vectors */
26: VecCreate(PETSC_COMM_WORLD,&x);
27: VecSetSizes(x,PETSC_DECIDE,size*bs*n);
28: VecSetFromOptions(x);
30: /* create two index sets */
31: if (rank < size-1) m = n + 2;
32: else m = n;
34: PetscMalloc1(m,&blks);
35: blks[0] = n*rank;
36: for (i=1; i<m; i++) blks[i] = blks[i-1] + 1;
37: ISCreateBlock(PETSC_COMM_SELF,bs,m,blks,PETSC_COPY_VALUES,&is1);
38: PetscFree(blks);
40: VecCreateSeq(PETSC_COMM_SELF,bs*m,&y);
41: ISCreateStride(PETSC_COMM_SELF,bs*m,0,1,&is2);
43: /* each processor inserts the entire vector */
44: /* this is redundant but tests assembly */
45: for (i=0; i<bs*n*size; i++) {
46: value = (PetscScalar) i;
47: VecSetValues(x,1,&i,&value,INSERT_VALUES);
48: }
49: VecAssemblyBegin(x);
50: VecAssemblyEnd(x);
51: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
53: VecScatterCreate(x,is1,y,is2,&ctx);
54: VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
55: VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
57: PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);
58: PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"----\n");
59: PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);
60: VecView(y,sviewer)); fflush(stdout;
61: PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);
62: PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
63: PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);
65: VecScatterDestroy(&ctx);
67: VecDestroy(&x);
68: VecDestroy(&y);
69: ISDestroy(&is1);
70: ISDestroy(&is2);
72: PetscFinalize();
73: return 0;
74: }
76: /*TEST
78: testset:
79: nsize: 3
80: output_file: output/ex24_1.out
81: filter: grep -v " type:"
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
94: TEST*/