Actual source code: ex12.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: {
 10:   PetscInt       bs=1,n=5,i,low;
 11:   PetscInt       ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,4},iy1[3] = {0,1,3};
 12:   PetscMPIInt    size,rank;
 13:   PetscScalar    *array;
 14:   Vec            x,y;
 15:   IS             isx,isy;
 16:   VecScatter     ctx;

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

 22:   if (size <2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor");

 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:   VecSet(y,0.0);

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

 52:   if (rank == 10) {
 53:     PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank);
 54:     ISView(isx,PETSC_VIEWER_STDOUT_SELF);
 55:   }

 57:   VecScatterCreate(x,isx,y,isy,&ctx);
 58:   VecScatterSetFromOptions(ctx);

 60:   /* Test forward vecscatter */
 61:   VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);
 62:   VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);
 63:   if (rank == 0) {
 64:     PetscPrintf(PETSC_COMM_SELF,"[%d] y:\n",rank);
 65:     VecView(y,PETSC_VIEWER_STDOUT_SELF);
 66:   }

 68:   /* Test reverse vecscatter */
 69:   VecScale(y,-1.0);
 70:   if (rank) {
 71:     VecScale(y,1.0/(size - 1));
 72:   }

 74:   VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_REVERSE);
 75:   VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_REVERSE);
 76:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

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

 88: /*TEST

 90:    test:
 91:       nsize: 3

 93: TEST*/