xref: /petsc/src/vec/vec/tests/ex45.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 
2 static char help[] = "Demonstrates VecStrideSubSetScatter() and VecStrideSubSetGather().\n\n";
3 
4 /*T
5    Concepts: vectors^sub-vectors;
6    Processors: n
7 
8    Allows one to easily pull out some components of a multi-component vector and put them in another vector.
9 
10    Note that these are special cases of VecScatter
11 T*/
12 
13 /*
14   Include "petscvec.h" so that we can use vectors.  Note that this file
15   automatically includes:
16      petscsys.h       - base PETSc routines   petscis.h     - index sets
17      petscviewer.h - viewers
18 */
19 
20 #include <petscvec.h>
21 
22 int main(int argc,char **argv)
23 {
24   Vec            v,s;
25   PetscInt       i,start,end,n = 8;
26   PetscScalar    value;
27   const PetscInt vidx[] = {1,2},sidx[] = {1,0};
28   PetscInt       miidx[2];
29   PetscReal      mvidx[2];
30 
31   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
32   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
33 
34   /*
35       Create multi-component vector with 4 components
36   */
37   PetscCall(VecCreate(PETSC_COMM_WORLD,&v));
38   PetscCall(VecSetSizes(v,PETSC_DECIDE,n));
39   PetscCall(VecSetBlockSize(v,4));
40   PetscCall(VecSetFromOptions(v));
41 
42   /*
43       Create double-component vectors
44   */
45   PetscCall(VecCreate(PETSC_COMM_WORLD,&s));
46   PetscCall(VecSetSizes(s,PETSC_DECIDE,n/2));
47   PetscCall(VecSetBlockSize(s,2));
48   PetscCall(VecSetFromOptions(s));
49 
50   /*
51      Set the vector values
52   */
53   PetscCall(VecGetOwnershipRange(v,&start,&end));
54   for (i=start; i<end; i++) {
55     value = i;
56     PetscCall(VecSetValues(v,1,&i,&value,INSERT_VALUES));
57   }
58 
59   /*
60      Get the components from the large multi-component vector to the small multi-component vector,
61      scale the smaller vector and then move values back to the large vector
62   */
63   PetscCall(VecStrideSubSetGather(v,PETSC_DETERMINE,vidx,NULL,s,INSERT_VALUES));
64   PetscCall(VecView(s,PETSC_VIEWER_STDOUT_WORLD));
65   PetscCall(VecScale(s,100.0));
66 
67   PetscCall(VecStrideSubSetScatter(s,PETSC_DETERMINE,NULL,vidx,v,ADD_VALUES));
68   PetscCall(VecView(v,PETSC_VIEWER_STDOUT_WORLD));
69 
70   /*
71      Get the components from the large multi-component vector to the small multi-component vector,
72      scale the smaller vector and then move values back to the large vector
73   */
74   PetscCall(VecStrideSubSetGather(v,2,vidx,sidx,s,INSERT_VALUES));
75   PetscCall(VecView(s,PETSC_VIEWER_STDOUT_WORLD));
76   PetscCall(VecScale(s,100.0));
77 
78   PetscCall(VecStrideSubSetScatter(s,2,sidx,vidx,v,ADD_VALUES));
79   PetscCall(VecView(v,PETSC_VIEWER_STDOUT_WORLD));
80 
81   PetscCall(VecStrideMax(v,1,&miidx[0],&mvidx[0]));
82   PetscCall(VecStrideMin(v,1,&miidx[1],&mvidx[1]));
83   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Min/Max: %" PetscInt_FMT " %g, %" PetscInt_FMT " %g\n",miidx[0],(double)mvidx[0],miidx[1],(double)mvidx[1]));
84   /*
85      Free work space.  All PETSc objects should be destroyed when they
86      are no longer needed.
87   */
88   PetscCall(VecDestroy(&v));
89   PetscCall(VecDestroy(&s));
90   PetscCall(PetscFinalize());
91   return 0;
92 }
93 
94 /*TEST
95 
96    test:
97       filter: grep -v type | grep -v "MPI processes" | grep -v Process
98       diff_args: -j
99       nsize: 2
100 
101    test:
102       filter: grep -v type | grep -v "MPI processes" | grep -v Process
103       output_file: output/ex45_1.out
104       diff_args: -j
105       suffix: 2
106       nsize: 1
107       args: -vec_type {{seq mpi}}
108 
109 TEST*/
110