Actual source code: ex10.c

petsc-3.12.2 2019-11-22
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>

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

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

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

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

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

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

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

 49:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

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


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

 64:   VecScatterBegin(newctx,y,x,INSERT_VALUES,SCATTER_REVERSE);
 65:   VecScatterEnd(newctx,y,x,INSERT_VALUES,SCATTER_REVERSE);
 66:   VecScatterDestroy(&newctx);

 68:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

 70:   ISDestroy(&isx);
 71:   ISDestroy(&isy);
 72:   VecDestroy(&x);
 73:   VecDestroy(&y);

 75:   PetscFinalize();
 76:   return ierr;
 77: }



 81: /*TEST

 83:    test:
 84:       nsize: 2

 86: TEST*/