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>

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

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

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

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

 31:   value = rank+1;
 32:   VecSet(x,value);
 33:   VecSet(y,zero);

 35:   VecScatterCreate(x,is1,y,is2,&ctx);
 36:   VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);
 37:   VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_FORWARD);
 38:   VecScatterDestroy(&ctx);

 40:   VecView(y,PETSC_VIEWER_STDOUT_WORLD);

 42:   VecDestroy(&x);
 43:   VecDestroy(&y);
 44:   ISDestroy(&is1);
 45:   ISDestroy(&is2);

 47:   PetscFinalize();
 48:   return 0;
 49: }

 51: /*TEST

 53:    test:
 54:       nsize: 2

 56: TEST*/