Actual source code: ex14.c

  2: static char help[] = "Scatters from a sequential vector to a parallel vector.\n\
  3: This does the tricky case.\n\n";

 5:  #include petscvec.h

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

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

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

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

 34:   value = rank+1;
 35:   VecSet(x,value);
 36:   VecSet(y,zero);

 38:   VecScatterCreate(x,is1,y,is2,&ctx);
 39:   VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);
 40:   VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);
 41:   VecScatterDestroy(ctx);
 42: 
 43:   VecView(y,PETSC_VIEWER_STDOUT_WORLD);

 45:   VecDestroy(x);
 46:   VecDestroy(y);
 47:   ISDestroy(is1);
 48:   ISDestroy(is2);

 50:   PetscFinalize();
 51:   return 0;
 52: }
 53: