Actual source code: ex5.c

  1: /*$Id: ex5.c,v 1.50 2001/03/23 23:21:30 balay Exp $*/

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

  7: #include "petscvec.h"
  8: #include "petscsys.h"

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

 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:   VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,size*n,&x);
 24:   VecCreateSeq(PETSC_COMM_SELF,n,&y);

 26:   /* create two index sets */
 27:   ISCreateStride(PETSC_COMM_SELF,n,n*rank,1,&is1);
 28:   ISCreateStride(PETSC_COMM_SELF,n,0,1,&is2);

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

 40:   VecScatterCreate(x,is1,y,is2,&ctx);
 41:   VecScatterBegin(x,y,INSERT_VALUES,SCATTER_FORWARD,ctx);
 42:   VecScatterEnd(x,y,INSERT_VALUES,SCATTER_FORWARD,ctx);
 43:   VecScatterDestroy(ctx);
 44: 
 45:   if (!rank)
 46:    {printf("----n"); VecView(y,PETSC_VIEWER_STDOUT_SELF);}

 48:   VecDestroy(x);
 49:   VecDestroy(y);
 50:   ISDestroy(is1);
 51:   ISDestroy(is2);

 53:   PetscFinalize();
 54:   return 0;
 55: }
 56: