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.  Inn
  4: this case processor zero is as long as the entire parallel vector; rest are zero length.nn";

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

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

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

 22:   /* create two vectors */
 23:   N = size*n;
 24:   VecCreate(PETSC_COMM_WORLD,&y);
 25:   VecSetSizes(y,PETSC_DECIDE,N);
 26:   VecSetFromOptions(y);
 27:   if (!rank) {
 28:     VecCreateSeq(PETSC_COMM_SELF,N,&x);
 29:   } else {
 30:     VecCreateSeq(PETSC_COMM_SELF,0,&x);
 31:   }

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

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

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

 60:   VecDestroy(x);
 61:   VecDestroy(y);
 62:   ISDestroy(is1);
 63:   ISDestroy(is2);

 65:   PetscFinalize();
 66:   return 0;
 67: }
 68: