Actual source code: ex23.c

petsc-3.8.3 2017-12-09
Report Typos and Errors

  2: static char help[] = "Scatters from a parallel vector to a sequential vector.\n\
  3:   Using a blocked send and a strided receive.\n\n";

  5: /*
  6:         0 1 2 3 | 4 5 6 7 ||  8 9 10 11

  8:      Scatter first and third block to first processor and
  9:      second and third block to second processor
 10: */
 11:  #include <petscvec.h>

 13: int main(int argc,char **argv)
 14: {
 16:   PetscInt       i,blocks[2],nlocal;
 17:   PetscMPIInt    size,rank;
 18:   PetscScalar    value;
 19:   Vec            x,y;
 20:   IS             is1,is2;
 21:   VecScatter     ctx = 0;

 23:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 24:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 25:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

 27:   if (size != 2) SETERRQ(PETSC_COMM_SELF,1,"Must run with 2 processors");

 29:   /* create two vectors */
 30:   if (!rank) nlocal = 8;
 31:   else nlocal = 4;
 32:   VecCreate(PETSC_COMM_WORLD,&x);
 33:   VecSetSizes(x,nlocal,12);
 34:   VecSetFromOptions(x);
 35:   VecCreateSeq(PETSC_COMM_SELF,8,&y);

 37:   /* create two index sets */
 38:   if (!rank) {
 39:     blocks[0] = 0; blocks[1] = 2;
 40:   } else {
 41:     blocks[0] = 1; blocks[1] = 2;
 42:   }
 43:   ISCreateBlock(PETSC_COMM_SELF,4,2,blocks,PETSC_COPY_VALUES,&is1);
 44:   ISCreateStride(PETSC_COMM_SELF,8,0,1,&is2);

 46:   for (i=0; i<12; i++) {
 47:     value = i;
 48:     VecSetValues(x,1,&i,&value,INSERT_VALUES);
 49:   }
 50:   VecAssemblyBegin(x);
 51:   VecAssemblyEnd(x);

 53:   VecScatterCreate(x,is1,y,is2,&ctx);
 54:   VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
 55:   VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
 56:   VecScatterDestroy(&ctx);

 58:   PetscSleep(2.0*rank);
 59:   VecView(y,PETSC_VIEWER_STDOUT_SELF);

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

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