Actual source code: ex33.c


  2: static char help[] = "Tests the routines VecScatterCreateToAll(), VecScatterCreateToZero()\n\n";

  4: #include <petscvec.h>

  6: int main(int argc,char **argv)
  7: {
  8:   PetscInt       n = 3,i,len,start,end;
  9:   PetscMPIInt    size,rank;
 10:   PetscScalar    value,*yy;
 11:   Vec            x,y,z,y_t;
 12:   VecScatter     toall,tozero;

 14:   PetscInitialize(&argc,&argv,(char*)0,help);
 15:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 16:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);

 18:   /* create two vectors */
 19:   VecCreateMPI(PETSC_COMM_WORLD,PETSC_DECIDE,size*n,&x);

 21:   /* each processor inserts its values */

 23:   VecGetOwnershipRange(x,&start,&end);
 24:   for (i=start; i<end; i++) {
 25:     value = (PetscScalar) i;
 26:     VecSetValues(x,1,&i,&value,INSERT_VALUES);
 27:   }
 28:   VecAssemblyBegin(x);
 29:   VecAssemblyEnd(x);
 30:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

 32:   VecScatterCreateToAll(x,&toall,&y);
 33:   VecScatterBegin(toall,x,y,INSERT_VALUES,SCATTER_FORWARD);
 34:   VecScatterEnd(toall,x,y,INSERT_VALUES,SCATTER_FORWARD);
 35:   VecScatterDestroy(&toall);

 37:   /* Cannot view the above vector with VecView(), so place it in an MPI Vec
 38:      and do a VecView() */
 39:   VecGetArray(y,&yy);
 40:   VecGetLocalSize(y,&len);
 41:   VecCreateMPIWithArray(PETSC_COMM_WORLD,1,len,PETSC_DECIDE,yy,&y_t);
 42:   VecView(y_t,PETSC_VIEWER_STDOUT_WORLD);
 43:   VecDestroy(&y_t);
 44:   VecRestoreArray(y,&yy);

 46:   VecScatterCreateToAll(x,&tozero,&z);
 47:   VecScatterBegin(tozero,x,z,INSERT_VALUES,SCATTER_FORWARD);
 48:   VecScatterEnd(tozero,x,z,INSERT_VALUES,SCATTER_FORWARD);
 49:   VecScatterDestroy(&tozero);
 50:   if (rank == 0) {
 51:     VecView(z,PETSC_VIEWER_STDOUT_SELF);
 52:   }
 53:   VecDestroy(&z);

 55:   VecScatterCreateToZero(x,&tozero,&z);
 56:   VecScatterBegin(tozero,x,z,INSERT_VALUES,SCATTER_FORWARD);
 57:   VecScatterEnd(tozero,x,z,INSERT_VALUES,SCATTER_FORWARD);
 58:   VecScatterDestroy(&tozero);
 59:   VecDestroy(&z);

 61:   VecDestroy(&x);
 62:   VecDestroy(&y);

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

 68: /*TEST

 70:    test:
 71:       nsize: 4

 73: TEST*/