Actual source code: ex25.c

  1: /*$Id: ex25.c,v 1.17 2001/09/11 16:32:10 bsmith Exp $*/

  3: static char help[] = "Scatters from a parallel vector to a sequential vector.  In\n\
  4: this case processor zero is as long as the entire parallel vector; rest are zero length.\n\n";

 6:  #include petscvec.h
 7:  #include petscsys.h

 11: int main(int argc,char **argv)
 12: {
 13:   int           n = 5,ierr;
 14:   int           size,rank,N,low,high,iglobal,i;
 15:   PetscScalar   value,zero = 0.0;
 16:   Vec           x,y;
 17:   IS            is1,is2;
 18:   VecScatter    ctx;

 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:   N = size*n;
 26:   VecCreate(PETSC_COMM_WORLD,&y);
 27:   VecSetSizes(y,PETSC_DECIDE,N);
 28:   VecSetFromOptions(y);
 29:   if (!rank) {
 30:     VecCreateSeq(PETSC_COMM_SELF,N,&x);
 31:   } else {
 32:     VecCreateSeq(PETSC_COMM_SELF,0,&x);
 33:   }

 35:   /* create two index sets */
 36:   if (!rank) {
 37:     ISCreateStride(PETSC_COMM_SELF,N,0,1,&is1);
 38:     ISCreateStride(PETSC_COMM_SELF,N,0,1,&is2);
 39:   } else {
 40:     ISCreateStride(PETSC_COMM_SELF,0,0,1,&is1);
 41:     ISCreateStride(PETSC_COMM_SELF,0,0,1,&is2);
 42:   }

 44:   VecSet(&zero,x);
 45:   VecGetOwnershipRange(y,&low,&high);
 46:   for (i=0; i<n; i++) {
 47:     iglobal = i + low; value = (PetscScalar) (i + 10*rank);
 48:     VecSetValues(y,1,&iglobal,&value,INSERT_VALUES);
 49:   }
 50:   VecAssemblyBegin(y);
 51:   VecAssemblyEnd(y);
 52:   VecView(y,PETSC_VIEWER_STDOUT_WORLD);

 54:   VecScatterCreate(y,is2,x,is1,&ctx);
 55:   VecScatterBegin(y,x,ADD_VALUES,SCATTER_FORWARD,ctx);
 56:   VecScatterEnd(y,x,ADD_VALUES,SCATTER_FORWARD,ctx);
 57:   VecScatterDestroy(ctx);
 58: 
 59:   if (!rank)
 60:     {printf("----\n"); VecView(x,PETSC_VIEWER_STDOUT_SELF);}

 62:   VecDestroy(x);
 63:   VecDestroy(y);
 64:   ISDestroy(is1);
 65:   ISDestroy(is2);

 67:   PetscFinalize();
 68:   return 0;
 69: }
 70: