Actual source code: ex13.c


  2: static char help[]= "Scatters from a sequential vector to a parallel 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,i,low;
 10:   PetscInt       ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,4},iy1[3] = {0,1,3};
 11:   PetscMPIInt    size,rank;
 12:   PetscScalar    *array;
 13:   Vec            x,y;
 14:   IS             isx,isy;
 15:   VecScatter     ctx;
 16:   PetscViewer    sviewer;

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


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

 27:   /* Create vector x over shared memory */
 28:   VecCreate(PETSC_COMM_WORLD,&x);
 29:   VecSetSizes(x,n,PETSC_DECIDE);
 30:   VecSetFromOptions(x);

 32:   VecGetOwnershipRange(x,&low,NULL);
 33:   VecGetArray(x,&array);
 34:   for (i=0; i<n; i++) {
 35:     array[i] = (PetscScalar)(i + low);
 36:   }
 37:   VecRestoreArray(x,&array);

 39:   /* Create a sequential vector y */
 40:   VecCreateSeq(PETSC_COMM_SELF,n,&y);
 41:   VecGetArray(y,&array);
 42:   for (i=0; i<n; i++) {
 43:     array[i] = (PetscScalar)(i + 100*rank);
 44:   }
 45:   VecRestoreArray(y,&array);

 47:   /* Create two index sets */
 48:   if (rank == 0) {
 49:     ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx);
 50:     ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy);
 51:   } else {
 52:     ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx);
 53:     ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy);
 54:   }

 56:   if (rank == 10) {
 57:     PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank);
 58:     ISView(isx,PETSC_VIEWER_STDOUT_SELF);
 59:   }

 61:   VecScatterCreate(y,isy,x,isx,&ctx);
 62:   VecScatterSetFromOptions(ctx);

 64:   /* Test forward vecscatter */
 65:   VecSet(x,0.0);
 66:   VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_FORWARD);
 67:   VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_FORWARD);
 68:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

 70:   /* Test reverse vecscatter */
 71:   VecScale(x,-1.0);
 72:   VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_REVERSE);
 73:   VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_REVERSE);
 74:   PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);
 75:   if (rank == 1) {
 76:     VecView(y,sviewer);
 77:   }
 78:   PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);

 80:   /* Free spaces */
 81:   VecScatterDestroy(&ctx);
 82:   ISDestroy(&isx);
 83:   ISDestroy(&isy);
 84:   VecDestroy(&x);
 85:   VecDestroy(&y);
 86:   PetscFinalize();
 87:   return 0;
 88: }

 90: /*TEST

 92:    test:
 93:       nsize: 3

 95: TEST*/