Actual source code: ex10.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: static char help[]= "Scatters from a parallel vector to a sequential vector.\n\
  3: uses block index sets\n\n";

  5: #include <petscvec.h>

  9: int main(int argc,char **argv)
 10: {
 12:   PetscInt       bs = 1,n = 5,ix0[3] = {5,7,9},ix1[3] = {2,3,4},i,iy0[3] = {1,2,4},iy1[3] = {0,1,3};
 13:   PetscMPIInt    size,rank;
 14:   PetscScalar    value;
 15:   Vec            x,y;
 16:   IS             isx,isy;
 17:   VecScatter     ctx = 0,newctx;

 19:   PetscInitialize(&argc,&argv,(char*)0,help);
 20:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 21:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

 23:   if (size != 2) SETERRQ(PETSC_COMM_SELF,1,"Must run with 2 processors");

 25:   PetscOptionsGetInt(NULL,"-bs",&bs,NULL);
 26:   n    = bs*n;

 28:   /* create two vectors */
 29:   VecCreate(PETSC_COMM_WORLD,&x);
 30:   VecSetSizes(x,PETSC_DECIDE,size*n);
 31:   VecSetFromOptions(x);
 32:   VecCreateSeq(PETSC_COMM_SELF,n,&y);

 34:   /* create two index sets */
 35:   if (!rank) {
 36:     ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx);
 37:     ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy);
 38:   } else {
 39:     ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx);
 40:     ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy);
 41:   }

 43:   /* fill local part of parallel vector */
 44:   for (i=n*rank; i<n*(rank+1); i++) {
 45:     value = (PetscScalar) i;
 46:     VecSetValues(x,1,&i,&value,INSERT_VALUES);
 47:   }
 48:   VecAssemblyBegin(x);
 49:   VecAssemblyEnd(x);

 51:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

 53:   /* fill local part of parallel vector */
 54:   for (i=0; i<n; i++) {
 55:     value = -(PetscScalar) (i + 100*rank);
 56:     VecSetValues(y,1,&i,&value,INSERT_VALUES);
 57:   }
 58:   VecAssemblyBegin(y);
 59:   VecAssemblyEnd(y);


 62:   VecScatterCreate(x,isx,y,isy,&ctx);
 63:   VecScatterCopy(ctx,&newctx);
 64:   VecScatterDestroy(&ctx);

 66:   VecScatterBegin(newctx,y,x,INSERT_VALUES,SCATTER_REVERSE);
 67:   VecScatterEnd(newctx,y,x,INSERT_VALUES,SCATTER_REVERSE);
 68:   VecScatterDestroy(&newctx);

 70:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

 72:   ISDestroy(&isx);
 73:   ISDestroy(&isy);
 74:   VecDestroy(&x);
 75:   VecDestroy(&y);

 77:   PetscFinalize();
 78:   return 0;
 79: }