xref: /petsc/src/vec/is/sf/tests/ex13.c (revision 97929ea760a70c780fc4b78d3065eb9e422113fe)
1*97929ea7SJunchao Zhang 
2*97929ea7SJunchao Zhang static char help[]= "Scatters from a sequential vector to a parallel vector. \n\
3*97929ea7SJunchao Zhang uses block index sets\n\n";
4*97929ea7SJunchao Zhang 
5*97929ea7SJunchao Zhang #include <petscvec.h>
6*97929ea7SJunchao Zhang 
7*97929ea7SJunchao Zhang int main(int argc,char **argv)
8*97929ea7SJunchao Zhang {
9*97929ea7SJunchao Zhang   PetscErrorCode ierr;
10*97929ea7SJunchao Zhang   PetscInt       bs=1,n=5,i,low;
11*97929ea7SJunchao Zhang   PetscInt       ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,4},iy1[3] = {0,1,3};
12*97929ea7SJunchao Zhang   PetscMPIInt    size,rank;
13*97929ea7SJunchao Zhang   PetscScalar    *array;
14*97929ea7SJunchao Zhang   Vec            x,y;
15*97929ea7SJunchao Zhang   IS             isx,isy;
16*97929ea7SJunchao Zhang   VecScatter     ctx;
17*97929ea7SJunchao Zhang   PetscViewer    sviewer;
18*97929ea7SJunchao Zhang 
19*97929ea7SJunchao Zhang   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
20*97929ea7SJunchao Zhang   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
21*97929ea7SJunchao Zhang   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
22*97929ea7SJunchao Zhang 
23*97929ea7SJunchao Zhang   if (size <2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor");
24*97929ea7SJunchao Zhang 
25*97929ea7SJunchao Zhang   ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr);
26*97929ea7SJunchao Zhang   n    = bs*n;
27*97929ea7SJunchao Zhang 
28*97929ea7SJunchao Zhang   /* Create vector x over shared memory */
29*97929ea7SJunchao Zhang   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
30*97929ea7SJunchao Zhang   ierr = VecSetSizes(x,n,PETSC_DECIDE);CHKERRQ(ierr);
31*97929ea7SJunchao Zhang   ierr = VecSetType(x,VECNODE);CHKERRQ(ierr);
32*97929ea7SJunchao Zhang   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
33*97929ea7SJunchao Zhang 
34*97929ea7SJunchao Zhang   ierr = VecGetOwnershipRange(x,&low,NULL);CHKERRQ(ierr);
35*97929ea7SJunchao Zhang   ierr = VecGetArray(x,&array);CHKERRQ(ierr);
36*97929ea7SJunchao Zhang   for (i=0; i<n; i++) {
37*97929ea7SJunchao Zhang     array[i] = (PetscScalar)(i + low);
38*97929ea7SJunchao Zhang   }
39*97929ea7SJunchao Zhang   ierr = VecRestoreArray(x,&array);CHKERRQ(ierr);
40*97929ea7SJunchao Zhang 
41*97929ea7SJunchao Zhang   /* Create a sequential vector y */
42*97929ea7SJunchao Zhang   ierr = VecCreateSeq(PETSC_COMM_SELF,n,&y);CHKERRQ(ierr);
43*97929ea7SJunchao Zhang   ierr = VecGetArray(y,&array);CHKERRQ(ierr);
44*97929ea7SJunchao Zhang   for (i=0; i<n; i++) {
45*97929ea7SJunchao Zhang     array[i] = (PetscScalar)(i + 100*rank);
46*97929ea7SJunchao Zhang   }
47*97929ea7SJunchao Zhang   ierr = VecRestoreArray(y,&array);CHKERRQ(ierr);
48*97929ea7SJunchao Zhang 
49*97929ea7SJunchao Zhang   /* Create two index sets */
50*97929ea7SJunchao Zhang   if (!rank) {
51*97929ea7SJunchao Zhang     ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr);
52*97929ea7SJunchao Zhang     ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr);
53*97929ea7SJunchao Zhang   } else {
54*97929ea7SJunchao Zhang     ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx);CHKERRQ(ierr);
55*97929ea7SJunchao Zhang     ierr = ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy);CHKERRQ(ierr);
56*97929ea7SJunchao Zhang   }
57*97929ea7SJunchao Zhang 
58*97929ea7SJunchao Zhang   if (rank == 10) {
59*97929ea7SJunchao Zhang     ierr = PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank);CHKERRQ(ierr);
60*97929ea7SJunchao Zhang     ierr = ISView(isx,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr);
61*97929ea7SJunchao Zhang   }
62*97929ea7SJunchao Zhang 
63*97929ea7SJunchao Zhang   ierr = VecScatterCreate(y,isy,x,isx,&ctx);CHKERRQ(ierr);
64*97929ea7SJunchao Zhang   ierr = VecScatterSetFromOptions(ctx);CHKERRQ(ierr);
65*97929ea7SJunchao Zhang 
66*97929ea7SJunchao Zhang   /* Test forward vecscatter */
67*97929ea7SJunchao Zhang   ierr = VecSet(x,0.0);CHKERRQ(ierr);
68*97929ea7SJunchao Zhang   ierr = VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
69*97929ea7SJunchao Zhang   ierr = VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
70*97929ea7SJunchao Zhang   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
71*97929ea7SJunchao Zhang 
72*97929ea7SJunchao Zhang   /* Test reverse vecscatter */
73*97929ea7SJunchao Zhang   ierr = VecScale(x,-1.0);CHKERRQ(ierr);
74*97929ea7SJunchao Zhang   ierr = VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
75*97929ea7SJunchao Zhang   ierr = VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
76*97929ea7SJunchao Zhang   ierr = PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
77*97929ea7SJunchao Zhang   if (rank == 1) {
78*97929ea7SJunchao Zhang     ierr = VecView(y,sviewer);CHKERRQ(ierr);
79*97929ea7SJunchao Zhang   }
80*97929ea7SJunchao Zhang   ierr = PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
81*97929ea7SJunchao Zhang 
82*97929ea7SJunchao Zhang   /* Free spaces */
83*97929ea7SJunchao Zhang   ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr);
84*97929ea7SJunchao Zhang   ierr = ISDestroy(&isx);CHKERRQ(ierr);
85*97929ea7SJunchao Zhang   ierr = ISDestroy(&isy);CHKERRQ(ierr);
86*97929ea7SJunchao Zhang   ierr = VecDestroy(&x);CHKERRQ(ierr);
87*97929ea7SJunchao Zhang   ierr = VecDestroy(&y);CHKERRQ(ierr);
88*97929ea7SJunchao Zhang   ierr = PetscFinalize();
89*97929ea7SJunchao Zhang   return ierr;
90*97929ea7SJunchao Zhang }
91*97929ea7SJunchao Zhang 
92*97929ea7SJunchao Zhang /*TEST
93*97929ea7SJunchao Zhang 
94*97929ea7SJunchao Zhang    test:
95*97929ea7SJunchao Zhang       nsize: 3
96*97929ea7SJunchao Zhang 
97*97929ea7SJunchao Zhang TEST*/
98