Actual source code: ex25.c

petsc-3.12.1 2019-10-22
Report Typos and Errors

  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: {
 10:   PetscInt       n = 5,N,low,high,iglobal,i;
 11:   PetscMPIInt    size,rank;
 12:   PetscScalar    value,zero = 0.0;
 13:   Vec            x,y;
 14:   IS             is1,is2;
 15:   VecScatter     ctx;

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

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

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

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

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

 56:   if (!rank) {
 57:     printf("----\n");
 58:     VecView(x,PETSC_VIEWER_STDOUT_SELF);
 59:   }

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

 66:   PetscFinalize();
 67:   return ierr;
 68: }



 72: /*TEST

 74:    test:
 75:       nsize: 3

 77: TEST*/