Actual source code: ex24.c
petsc-3.9.1 2018-04-29
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: /*T
6: requires: x
7: T*/
9: #include <petscvec.h>
11: int main(int argc,char **argv)
12: {
14: PetscMPIInt size,rank;
15: PetscInt n = 5,i,*blks,bs = 1,m = 2;
16: PetscScalar value;
17: Vec x,y;
18: IS is1,is2;
19: VecScatter ctx = 0;
20: PetscViewer sviewer;
22: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
24: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
25: PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);
27: MPI_Comm_size(PETSC_COMM_WORLD,&size);
28: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
30: /* create two vectors */
31: VecCreate(PETSC_COMM_WORLD,&x);
32: VecSetSizes(x,PETSC_DECIDE,size*bs*n);
33: VecSetFromOptions(x);
35: /* create two index sets */
36: if (rank < size-1) m = n + 2;
37: else m = n;
39: PetscMalloc1(m,&blks);
40: blks[0] = n*rank;
41: for (i=1; i<m; i++) blks[i] = blks[i-1] + 1;
42: ISCreateBlock(PETSC_COMM_SELF,bs,m,blks,PETSC_COPY_VALUES,&is1);
43: PetscFree(blks);
45: VecCreateSeq(PETSC_COMM_SELF,bs*m,&y);
46: ISCreateStride(PETSC_COMM_SELF,bs*m,0,1,&is2);
48: /* each processor inserts the entire vector */
49: /* this is redundant but tests assembly */
50: for (i=0; i<bs*n*size; i++) {
51: value = (PetscScalar) i;
52: VecSetValues(x,1,&i,&value,INSERT_VALUES);
53: }
54: VecAssemblyBegin(x);
55: VecAssemblyEnd(x);
56: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
58: VecScatterCreate(x,is1,y,is2,&ctx);
59: VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
60: VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
62: PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);
63: PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"----\n");
64: PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);
65: VecView(y,sviewer); fflush(stdout);
66: PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);
67: PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);
68: PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);
70: VecScatterDestroy(&ctx);
72: VecDestroy(&x);
73: VecDestroy(&y);
74: ISDestroy(&is1);
75: ISDestroy(&is2);
77: PetscFinalize();
78: return ierr;
79: }
83: /*TEST
85: test:
86: nsize: 3
88: TEST*/