xref: /petsc/src/vec/vec/tests/ex24.c (revision efa12513287cff49a2b9648ae83199dcbfaad71a)
1 
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";
4 
5 #include <petscvec.h>
6 
7 int main(int argc,char **argv)
8 {
9   PetscErrorCode ierr;
10   PetscMPIInt    size,rank;
11   PetscInt       n = 5,i,*blks,bs = 1,m = 2;
12   PetscScalar    value;
13   Vec            x,y;
14   IS             is1,is2;
15   VecScatter     ctx = 0;
16   PetscViewer    sviewer;
17 
18   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
19 
20   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
21   ierr = PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);CHKERRQ(ierr);
22 
23   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
24   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
25 
26   /* create two vectors */
27   ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr);
28   ierr = VecSetSizes(x,PETSC_DECIDE,size*bs*n);CHKERRQ(ierr);
29   ierr = VecSetFromOptions(x);CHKERRQ(ierr);
30 
31   /* create two index sets */
32   if (rank < size-1) m = n + 2;
33   else m = n;
34 
35   ierr = PetscMalloc1(m,&blks);CHKERRQ(ierr);
36   blks[0] = n*rank;
37   for (i=1; i<m; i++) blks[i] = blks[i-1] + 1;
38   ierr = ISCreateBlock(PETSC_COMM_SELF,bs,m,blks,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr);
39   ierr = PetscFree(blks);CHKERRQ(ierr);
40 
41   ierr = VecCreateSeq(PETSC_COMM_SELF,bs*m,&y);CHKERRQ(ierr);
42   ierr = ISCreateStride(PETSC_COMM_SELF,bs*m,0,1,&is2);CHKERRQ(ierr);
43 
44   /* each processor inserts the entire vector */
45   /* this is redundant but tests assembly */
46   for (i=0; i<bs*n*size; i++) {
47     value = (PetscScalar) i;
48     ierr  = VecSetValues(x,1,&i,&value,INSERT_VALUES);CHKERRQ(ierr);
49   }
50   ierr = VecAssemblyBegin(x);CHKERRQ(ierr);
51   ierr = VecAssemblyEnd(x);CHKERRQ(ierr);
52   ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
53 
54   ierr = VecScatterCreate(x,is1,y,is2,&ctx);CHKERRQ(ierr);
55   ierr = VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
56   ierr = VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
57 
58   ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
59   ierr = PetscViewerASCIISynchronizedPrintf(PETSC_VIEWER_STDOUT_WORLD,"----\n");CHKERRQ(ierr);
60   ierr = PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
61   ierr = VecView(y,sviewer);CHKERRQ(ierr); fflush(stdout);
62   ierr = PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
63   ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
64   ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
65 
66   ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr);
67 
68   ierr = VecDestroy(&x);CHKERRQ(ierr);
69   ierr = VecDestroy(&y);CHKERRQ(ierr);
70   ierr = ISDestroy(&is1);CHKERRQ(ierr);
71   ierr = ISDestroy(&is2);CHKERRQ(ierr);
72 
73   ierr = PetscFinalize();
74   return ierr;
75 }
76 
77 
78 
79 /*TEST
80 
81    testset:
82       nsize: 3
83       output_file: output/ex24_1.out
84       filter: grep -v "  type:"
85       test:
86         suffix: standard
87         args: -vec_type standard
88       test:
89         requires: cuda
90         suffix: cuda
91         args: -vec_type cuda
92       test:
93         requires: viennacl
94         suffix:  viennacl
95         args: -vec_type viennacl
96 
97 TEST*/
98