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