Actual source code: ex11.c
petsc-3.15.0 2021-03-30
2: static char help[]= "Scatters between parallel vectors. \n\
3: uses block index sets\n\n";
5: #include <petscvec.h>
7: int main(int argc,char **argv)
8: {
10: PetscInt bs=1,n=5,N,i,low;
11: PetscInt ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,1},iy1[3] = {0,3,9};
12: PetscMPIInt size,rank;
13: PetscScalar *array;
14: Vec x,x1,y;
15: IS isx,isy;
16: VecScatter ctx;
17: VecScatterType type;
18: PetscBool flg;
20: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
21: MPI_Comm_size(PETSC_COMM_WORLD,&size);
22: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
24: if (size <2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor");
26: PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);
27: n = bs*n;
29: /* Create vector x over shared memory */
30: VecCreate(PETSC_COMM_WORLD,&x);
31: VecSetSizes(x,n,PETSC_DECIDE);
32: VecSetFromOptions(x);
34: VecGetOwnershipRange(x,&low,NULL);
35: VecGetArray(x,&array);
36: for (i=0; i<n; i++) {
37: array[i] = (PetscScalar)(i + low);
38: }
39: VecRestoreArray(x,&array);
40: /* VecView(x,PETSC_VIEWER_STDOUT_WORLD); */
42: /* Test some vector functions */
43: VecAssemblyBegin(x);
44: VecAssemblyEnd(x);
46: VecGetSize(x,&N);
47: VecGetLocalSize(x,&n);
49: VecDuplicate(x,&x1);
50: VecCopy(x,x1);
51: VecEqual(x,x1,&flg);
52: if (!flg) SETERRQ(PetscObjectComm((PetscObject)x),PETSC_ERR_ARG_WRONG,"x1 != x");
54: VecScale(x1,2.0);
55: VecSet(x1,10.0);
56: /* VecView(x1,PETSC_VIEWER_STDOUT_WORLD); */
58: /* Create vector y over shared memory */
59: VecCreate(PETSC_COMM_WORLD,&y);
60: VecSetSizes(y,n,PETSC_DECIDE);
61: VecSetFromOptions(y);
62: VecGetArray(y,&array);
63: for (i=0; i<n; i++) {
64: array[i] = -(PetscScalar) (i + 100*rank);
65: }
66: VecRestoreArray(y,&array);
67: VecAssemblyBegin(y);
68: VecAssemblyEnd(y);
69: /* VecView(y,PETSC_VIEWER_STDOUT_WORLD); */
71: /* Create two index sets */
72: if (!rank) {
73: ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx);
74: ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy);
75: } else {
76: ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx);
77: ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy);
78: }
80: if (rank == 10) {
81: PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank);
82: ISView(isx,PETSC_VIEWER_STDOUT_SELF);
83: PetscPrintf(PETSC_COMM_SELF,"\n[%d] isy:\n",rank);
84: ISView(isy,PETSC_VIEWER_STDOUT_SELF);
85: }
87: /* Create Vector scatter */
88: VecScatterCreate(x,isx,y,isy,&ctx);
89: VecScatterSetFromOptions(ctx);
90: VecScatterGetType(ctx,&type);
91: PetscPrintf(PETSC_COMM_WORLD,"scatter type %s\n",type);
93: /* Test forward vecscatter */
94: VecSet(y,0.0);
95: VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);
96: VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);
97: PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_FORWARD y:\n");
98: VecView(y,PETSC_VIEWER_STDOUT_WORLD);
100: /* Test reverse vecscatter */
101: VecSet(x,0.0);
102: VecSet(y,0.0);
103: VecGetOwnershipRange(y,&low,NULL);
104: VecGetArray(y,&array);
105: for (i=0; i<n; i++) {
106: array[i] = (PetscScalar)(i+ low);
107: }
108: VecRestoreArray(y,&array);
109: VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_REVERSE);
110: VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_REVERSE);
111: PetscPrintf(PETSC_COMM_WORLD,"\nSCATTER_REVERSE x:\n");
112: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
114: /* Free objects */
115: VecScatterDestroy(&ctx);
116: ISDestroy(&isx);
117: ISDestroy(&isy);
118: VecDestroy(&x);
119: VecDestroy(&x1);
120: VecDestroy(&y);
121: PetscFinalize();
122: return ierr;
123: }
125: /*TEST
127: test:
128: nsize: 2
130: TEST*/