Actual source code: ex10.c

petsc-3.9.0 2018-04-07
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: /*T
  6:    requires: x
  7: T*/

  9:  #include <petscvec.h>

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

 21:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 22:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 23:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

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

 27:   PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);
 28:   n    = bs*n;

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

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

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

 53:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

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


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

 68:   VecScatterBegin(newctx,y,x,INSERT_VALUES,SCATTER_REVERSE);
 69:   VecScatterEnd(newctx,y,x,INSERT_VALUES,SCATTER_REVERSE);
 70:   VecScatterDestroy(&newctx);

 72:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

 74:   ISDestroy(&isx);
 75:   ISDestroy(&isy);
 76:   VecDestroy(&x);
 77:   VecDestroy(&y);

 79:   PetscFinalize();
 80:   return ierr;
 81: }



 85: /*TEST

 87:    test:
 88:       nsize: 2

 90: TEST*/