Actual source code: ex5.c

  2: static char help[] = "Scatters from a parallel vector to a sequential vector.\n\
  3: This does case when we are merely selecting the local part of the\n\
  4: parallel vector.\n\n";

 6:  #include petscvec.h

 10: int main(int argc,char **argv)
 11: {
 13:   PetscMPIInt    size,rank;
 14:   PetscInt       n = 5,i;
 15:   PetscScalar    value;
 16:   Vec            x,y;
 17:   IS             is1,is2;
 18:   VecScatter     ctx = 0;

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

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

 30:   /* create two index sets */
 31:   ISCreateStride(PETSC_COMM_SELF,n,n*rank,1,&is1);
 32:   ISCreateStride(PETSC_COMM_SELF,n,0,1,&is2);

 34:   /* each processor inserts the entire vector */
 35:   /* this is redundant but tests assembly */
 36:   for (i=0; i<n*size; i++) {
 37:     value = (PetscScalar) i;
 38:     VecSetValues(x,1,&i,&value,INSERT_VALUES);
 39:   }
 40:   VecAssemblyBegin(x);
 41:   VecAssemblyEnd(x);
 42:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

 44:   VecScatterCreate(x,is1,y,is2,&ctx);
 45:   VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
 46:   VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
 47:   VecScatterDestroy(ctx);
 48: 
 49:   if (!rank)
 50:    {printf("----\n"); VecView(y,PETSC_VIEWER_STDOUT_SELF);}

 52:   VecDestroy(x);
 53:   VecDestroy(y);
 54:   ISDestroy(is1);
 55:   ISDestroy(is2);

 57:   PetscFinalize();
 58:   return 0;
 59: }
 60: