Actual source code: ex13.c

petsc-3.15.0 2021-03-30
Report Typos and Errors

  2: static char help[]= "Scatters from a sequential vector to a parallel vector. \n\
  3: uses block index sets\n\n";

  5: #include <petscvec.h>

  7: int main(int argc,char **argv)
  8: {
 10:   PetscInt       bs=1,n=5,i,low;
 11:   PetscInt       ix0[3] = {5,7,9},iy0[3] = {1,2,4},ix1[3] = {2,3,4},iy1[3] = {0,1,3};
 12:   PetscMPIInt    size,rank;
 13:   PetscScalar    *array;
 14:   Vec            x,y;
 15:   IS             isx,isy;
 16:   VecScatter     ctx;
 17:   PetscViewer    sviewer;

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

 23:   if (size <2) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Must run more than one processor");

 25:   PetscOptionsGetInt(NULL,NULL,"-bs",&bs,NULL);
 26:   n    = bs*n;

 28:   /* Create vector x over shared memory */
 29:   VecCreate(PETSC_COMM_WORLD,&x);
 30:   VecSetSizes(x,n,PETSC_DECIDE);
 31:   VecSetFromOptions(x);

 33:   VecGetOwnershipRange(x,&low,NULL);
 34:   VecGetArray(x,&array);
 35:   for (i=0; i<n; i++) {
 36:     array[i] = (PetscScalar)(i + low);
 37:   }
 38:   VecRestoreArray(x,&array);

 40:   /* Create a sequential vector y */
 41:   VecCreateSeq(PETSC_COMM_SELF,n,&y);
 42:   VecGetArray(y,&array);
 43:   for (i=0; i<n; i++) {
 44:     array[i] = (PetscScalar)(i + 100*rank);
 45:   }
 46:   VecRestoreArray(y,&array);

 48:   /* Create two index sets */
 49:   if (!rank) {
 50:     ISCreateBlock(PETSC_COMM_SELF,bs,3,ix0,PETSC_COPY_VALUES,&isx);
 51:     ISCreateBlock(PETSC_COMM_SELF,bs,3,iy0,PETSC_COPY_VALUES,&isy);
 52:   } else {
 53:     ISCreateBlock(PETSC_COMM_SELF,bs,3,ix1,PETSC_COPY_VALUES,&isx);
 54:     ISCreateBlock(PETSC_COMM_SELF,bs,3,iy1,PETSC_COPY_VALUES,&isy);
 55:   }

 57:   if (rank == 10) {
 58:     PetscPrintf(PETSC_COMM_SELF,"\n[%d] isx:\n",rank);
 59:     ISView(isx,PETSC_VIEWER_STDOUT_SELF);
 60:   }

 62:   VecScatterCreate(y,isy,x,isx,&ctx);
 63:   VecScatterSetFromOptions(ctx);

 65:   /* Test forward vecscatter */
 66:   VecSet(x,0.0);
 67:   VecScatterBegin(ctx,y,x,ADD_VALUES,SCATTER_FORWARD);
 68:   VecScatterEnd(ctx,y,x,ADD_VALUES,SCATTER_FORWARD);
 69:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

 71:   /* Test reverse vecscatter */
 72:   VecScale(x,-1.0);
 73:   VecScatterBegin(ctx,x,y,ADD_VALUES,SCATTER_REVERSE);
 74:   VecScatterEnd(ctx,x,y,ADD_VALUES,SCATTER_REVERSE);
 75:   PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);
 76:   if (rank == 1) {
 77:     VecView(y,sviewer);
 78:   }
 79:   PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);

 81:   /* Free spaces */
 82:   VecScatterDestroy(&ctx);
 83:   ISDestroy(&isx);
 84:   ISDestroy(&isy);
 85:   VecDestroy(&x);
 86:   VecDestroy(&y);
 87:   PetscFinalize();
 88:   return ierr;
 89: }

 91: /*TEST

 93:    test:
 94:       nsize: 3

 96: TEST*/