Actual source code: ex23.c

  2: static char help[] = "Scatters from a parallel vector to a sequential vector.\n\
  3:   Using a blocked send and a strided receive.\n\n";

  5: /*
  6:         0 1 2 3 | 4 5 6 7 ||  8 9 10 11 

  8:      Scatter first and third block to first processor and 
  9:      second and third block to second processor
 10: */
 11:  #include petscvec.h
 12:  #include petscsys.h

 16: int main(int argc,char **argv)
 17: {
 19:   PetscInt       i,blocks[2],nlocal;
 20:   PetscMPIInt    size,rank;
 21:   PetscScalar    value;
 22:   Vec            x,y;
 23:   IS             is1,is2;
 24:   VecScatter     ctx = 0;

 26:   PetscInitialize(&argc,&argv,(char*)0,help);
 27:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 28:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

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

 32:   /* create two vectors */
 33:   if (!rank) nlocal = 8;
 34:   else nlocal = 4;
 35:   VecCreate(PETSC_COMM_WORLD,&x);
 36:   VecSetSizes(x,nlocal,12);
 37:   VecSetFromOptions(x);
 38:   VecCreateSeq(PETSC_COMM_SELF,8,&y);

 40:   /* create two index sets */
 41:   if (!rank) {
 42:     blocks[0] = 0; blocks[1] = 8;
 43:   } else {
 44:     blocks[0] = 4; blocks[1] = 8;
 45:   }
 46:   ISCreateBlock(PETSC_COMM_SELF,4,2,blocks,&is1);
 47:   ISCreateStride(PETSC_COMM_SELF,8,0,1,&is2);

 49:   for (i=0; i<12; i++) {
 50:     value = i;
 51:     VecSetValues(x,1,&i,&value,INSERT_VALUES);
 52:   }
 53:   VecAssemblyBegin(x);
 54:   VecAssemblyEnd(x);

 56:   VecScatterCreate(x,is1,y,is2,&ctx);
 57:   VecScatterBegin(x,y,INSERT_VALUES,SCATTER_FORWARD,ctx);
 58:   VecScatterEnd(x,y,INSERT_VALUES,SCATTER_FORWARD,ctx);
 59:   VecScatterDestroy(ctx);
 60: 
 61:   PetscSleep(2*rank);
 62:   VecView(y,PETSC_VIEWER_STDOUT_SELF);

 64:   VecDestroy(x);
 65:   VecDestroy(y);
 66:   ISDestroy(is1);
 67:   ISDestroy(is2);

 69:   PetscFinalize();
 70:   return 0;
 71: }
 72: