Actual source code: ex25.c


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

  5: #include <petscvec.h>

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

 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:   if (rank == 0) {
 26:     VecCreateSeq(PETSC_COMM_SELF,N,&x);
 27:   } else {
 28:     VecCreateSeq(PETSC_COMM_SELF,0,&x);
 29:   }

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

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

 50:   VecScatterCreate(y,is2,x,is1,&ctx);
 51:   VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_FORWARD);
 52:   VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_FORWARD);
 53:   VecScatterDestroy(&ctx);

 55:   if (rank == 0) {
 56:     PetscPrintf(PETSC_COMM_SELF,"----\n");
 57:     VecView(x,PETSC_VIEWER_STDOUT_SELF);
 58:   }

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

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

 69: /*TEST

 71:    test:
 72:       nsize: 3

 74: TEST*/