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