Actual source code: ex2.c

  2: static char help[] = "Tests DAGlobalToNaturalAllCreate() using contour plotting for 2d DAs.\n\n";

 4:  #include petscda.h

  8: int main(int argc,char **argv)
  9: {
 10:   PetscInt       i,j,M = 10,N = 8,m = PETSC_DECIDE,n = PETSC_DECIDE;
 11:   PetscMPIInt    rank;
 13:   PetscTruth     flg = PETSC_FALSE;
 14:   DA             da;
 15:   PetscViewer    viewer;
 16:   Vec            localall,global;
 17:   PetscScalar    value,*vlocal;
 18:   DAPeriodicType ptype = DA_NONPERIODIC;
 19:   DAStencilType  stype = DA_STENCIL_BOX;
 20:   VecScatter     tolocalall,fromlocalall;
 21:   PetscInt       start,end;
 22: 

 24:   PetscInitialize(&argc,&argv,(char*)0,help);
 25:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",300,0,300,300,&viewer);

 27:   /* Read options */
 28:   PetscOptionsGetInt(PETSC_NULL,"-M",&M,PETSC_NULL);
 29:   PetscOptionsGetInt(PETSC_NULL,"-N",&N,PETSC_NULL);
 30:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 31:   PetscOptionsGetInt(PETSC_NULL,"-n",&n,PETSC_NULL);
 32:   PetscOptionsGetTruth(PETSC_NULL,"-star_stencil",&flg,PETSC_NULL);
 33:   if (flg) stype = DA_STENCIL_STAR;

 35:   /* Create distributed array and get vectors */
 36:   DACreate2d(PETSC_COMM_WORLD,ptype,stype,
 37:                     M,N,m,n,1,1,PETSC_NULL,PETSC_NULL,&da);
 38:   DACreateGlobalVector(da,&global);
 39:   VecCreateSeq(PETSC_COMM_SELF,M*N,&localall);

 41:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 42:   VecGetOwnershipRange(global,&start,&end);
 43:   for (i=start; i<end; i++) {
 44:     value = 5.0*rank;
 45:     VecSetValues(global,1,&i,&value,INSERT_VALUES);
 46:   }
 47:   VecView(global,viewer);

 49:   /*
 50:      Create Scatter from global DA parallel vector to local vector that
 51:    contains all entries
 52:   */
 53:   DAGlobalToNaturalAllCreate(da,&tolocalall);
 54:   DANaturalAllToGlobalCreate(da,&fromlocalall);

 56:   VecScatterBegin(tolocalall,global,localall,INSERT_VALUES,SCATTER_FORWARD);
 57:   VecScatterEnd(tolocalall,global,localall,INSERT_VALUES,SCATTER_FORWARD);

 59:   VecGetArray(localall,&vlocal);
 60:   for (j=0; j<N; j++) {
 61:     for (i=0; i<M; i++) {
 62:       *vlocal++ += i + j*M;
 63:     }
 64:   }
 65:   VecRestoreArray(localall,&vlocal);

 67:   /* scatter back to global vector */
 68:   VecScatterBegin(fromlocalall,localall,global,INSERT_VALUES,SCATTER_FORWARD);
 69:   VecScatterEnd(fromlocalall,localall,global,INSERT_VALUES,SCATTER_FORWARD);

 71:   VecView(global,viewer);

 73:   /* Free memory */
 74:   VecScatterDestroy(tolocalall);
 75:   VecScatterDestroy(fromlocalall);
 76:   PetscViewerDestroy(viewer);
 77:   VecDestroy(localall);
 78:   VecDestroy(global);
 79:   DADestroy(da);
 80:   PetscFinalize();
 81:   return 0;
 82: }
 83: