Actual source code: ex10.c


  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>

  7: int main(int argc,char **argv)
  8: {
  9:   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};
 10:   PetscMPIInt    size,rank;
 11:   PetscScalar    value;
 12:   Vec            x,y;
 13:   IS             isx,isy;
 14:   VecScatter     ctx = 0,newctx;

 16:   PetscInitialize(&argc,&argv,(char*)0,help);
 17:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 18:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);


 22:   PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);
 23:   n    = bs*n;

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

 31:   /* create two index sets */
 32:   if (rank == 0) {
 33:     ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx);
 34:     ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy);
 35:   } else {
 36:     ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx);
 37:     ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy);
 38:   }

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

 48:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

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

 58:   VecScatterCreate(x,isx,y,isy,&ctx);
 59:   VecScatterCopy(ctx,&newctx);
 60:   VecScatterDestroy(&ctx);

 62:   VecScatterBegin(newctx,y,x,INSERT_VALUES,SCATTER_REVERSE);
 63:   VecScatterEnd(newctx,y,x,INSERT_VALUES,SCATTER_REVERSE);
 64:   VecScatterDestroy(&newctx);

 66:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

 68:   ISDestroy(&isx);
 69:   ISDestroy(&isy);
 70:   VecDestroy(&x);
 71:   VecDestroy(&y);

 73:   PetscFinalize();
 74:   return 0;
 75: }

 77: /*TEST

 79:    test:
 80:       nsize: 2

 82: TEST*/