Actual source code: ex33.c

petsc-3.7.3 2016-07-24
Report Typos and Errors
  2: static char help[] = "Tests the routines VecScatterCreateToAll(), VecScatterCreateToZero()\n\n";

  4: #include <petscvec.h>

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

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

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

 24:   /* each processor inserts its values */

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

 35:   VecScatterCreateToAll(x,&toall,&y);
 36:   VecScatterBegin(toall,x,y,INSERT_VALUES,SCATTER_FORWARD);
 37:   VecScatterEnd(toall,x,y,INSERT_VALUES,SCATTER_FORWARD);
 38:   VecScatterDestroy(&toall);

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

 49:   VecScatterCreateToAll(x,&tozero,&z);
 50:   VecScatterBegin(tozero,x,z,INSERT_VALUES,SCATTER_FORWARD);
 51:   VecScatterEnd(tozero,x,z,INSERT_VALUES,SCATTER_FORWARD);
 52:   VecScatterDestroy(&tozero);
 53:   if (!rank) {
 54:     VecView(z,PETSC_VIEWER_STDOUT_SELF);
 55:   }
 56:   VecDestroy(&z);

 58:   VecScatterCreateToZero(x,&tozero,&z);
 59:   VecScatterBegin(tozero,x,z,INSERT_VALUES,SCATTER_FORWARD);
 60:   VecScatterEnd(tozero,x,z,INSERT_VALUES,SCATTER_FORWARD);
 61:   VecScatterDestroy(&tozero);
 62:   VecDestroy(&z);

 64:   VecDestroy(&x);
 65:   VecDestroy(&y);

 67:   PetscFinalize();
 68:   return 0;
 69: }