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: {
  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;

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


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

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

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

 38:   /* Create a sequential vector y */
 39:   VecCreateSeq(PETSC_COMM_SELF,n,&y);
 40:   VecSet(y,0.0);

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

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

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

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

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

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

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

 87: /*TEST

 89:    test:
 90:       nsize: 3

 92: TEST*/