Actual source code: ex9.c
petsc-3.12.0 2019-09-29
2: static char help[]= "Scatters from a parallel vector to a sequential vector.\n\n";
4: #include <petscvec.h>
6: int main(int argc,char **argv)
7: {
9: PetscInt n = 5,i,idx2[3] = {0,2,3},idx1[3] = {0,1,2};
10: PetscMPIInt size,rank;
11: PetscScalar value;
12: Vec x,y;
13: IS is1,is2;
14: VecScatter ctx = 0;
16: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
17: MPI_Comm_size(PETSC_COMM_WORLD,&size);
18: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
20: /* create two vectors */
21: VecCreate(PETSC_COMM_WORLD,&x);
22: VecSetSizes(x,PETSC_DECIDE,size*n);
23: VecSetFromOptions(x);
24: VecCreateSeq(PETSC_COMM_SELF,n,&y);
26: /* create two index sets */
27: ISCreateGeneral(PETSC_COMM_SELF,3,idx1,PETSC_COPY_VALUES,&is1);
28: ISCreateGeneral(PETSC_COMM_SELF,3,idx2,PETSC_COPY_VALUES,&is2);
30: /* fill local part of parallel vector */
31: for (i=n*rank; i<n*(rank+1); i++) {
32: value = (PetscScalar) i;
33: VecSetValues(x,1,&i,&value,INSERT_VALUES);
34: }
35: VecAssemblyBegin(x);
36: VecAssemblyEnd(x);
38: VecView(x,PETSC_VIEWER_STDOUT_WORLD);
40: VecSet(y,-1.0);
42: VecScatterCreate(x,is1,y,is2,&ctx);
43: VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
44: VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
45: VecScatterDestroy(&ctx);
47: if (!rank) {
48: PetscPrintf(PETSC_COMM_SELF,"scattered vector\n");
49: VecView(y,PETSC_VIEWER_STDOUT_SELF);
50: }
51: ISDestroy(&is1);
52: ISDestroy(&is2);
53: VecDestroy(&x);
54: VecDestroy(&y);
56: PetscFinalize();
57: return ierr;
58: }
62: /*TEST
64: test:
65: nsize: 2
67: TEST*/