Actual source code: ex24.c

  2: static char help[] = "Scatters from a parallel vector to a sequential vector.\n\
  3: Tests where the local part of the scatter is a copy.\n\n";

 5:  #include petscvec.h
 6:  #include petscsys.h

 10: int main(int argc,char **argv)
 11: {
 13:   PetscMPIInt    size,rank;
 14:   PetscInt       n = 5,i,*blks,bs = 1,m = 2;
 15:   PetscScalar    value;
 16:   Vec            x,y;
 17:   IS             is1,is2;
 18:   VecScatter     ctx = 0;
 19:   PetscViewer    sviewer;

 21:   PetscInitialize(&argc,&argv,(char*)0,help);

 23:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 24:   PetscOptionsGetInt(PETSC_NULL,"-bs",&bs,PETSC_NULL);

 26:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 27:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

 29:   /* create two vectors */
 30:   VecCreate(PETSC_COMM_WORLD,&x);
 31:   VecSetSizes(x,PETSC_DECIDE,size*bs*n);
 32:   VecSetFromOptions(x);

 34:   /* create two index sets */
 35:   if (rank < size-1) {
 36:     m = n + 2;
 37:   } else {
 38:     m = n;
 39:   }
 40:   PetscMalloc((m)*sizeof(PetscInt),&blks);
 41:   blks[0] = n*rank*bs;
 42:   for (i=1; i<m; i++) {
 43:     blks[i] = blks[i-1] + bs;
 44:   }
 45:   ISCreateBlock(PETSC_COMM_SELF,bs,m,blks,&is1);
 46:   PetscFree(blks);

 48:   VecCreateSeq(PETSC_COMM_SELF,bs*m,&y);
 49:   ISCreateStride(PETSC_COMM_SELF,bs*m,0,1,&is2);

 51:   /* each processor inserts the entire vector */
 52:   /* this is redundant but tests assembly */
 53:   for (i=0; i<bs*n*size; i++) {
 54:     value = (PetscScalar) i;
 55:     VecSetValues(x,1,&i,&value,INSERT_VALUES);
 56:   }
 57:   VecAssemblyBegin(x);
 58:   VecAssemblyEnd(x);
 59:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

 61:   VecScatterCreate(x,is1,y,is2,&ctx);
 62:   VecScatterBegin(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);
 63:   VecScatterEnd(ctx,x,y,INSERT_VALUES,SCATTER_FORWARD);

 65:   PetscSynchronizedPrintf(PETSC_COMM_WORLD,"----\n");
 66:   PetscViewerGetSingleton(PETSC_VIEWER_STDOUT_WORLD,&sviewer);
 67:   VecView(y,sviewer); fflush(stdout);
 68:   PetscViewerRestoreSingleton(PETSC_VIEWER_STDOUT_WORLD,&sviewer);
 69:   PetscSynchronizedFlush(PETSC_COMM_WORLD);

 71:   VecScatterDestroy(ctx);

 73:   VecDestroy(x);
 74:   VecDestroy(y);
 75:   ISDestroy(is1);
 76:   ISDestroy(is2);

 78:   PetscFinalize();
 79:   return 0;
 80: }
 81: