Actual source code: ex2.c

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

 4:  #include petscda.h
 5:  #include petscsys.h

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

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

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

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

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

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

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

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

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

 72:   VecView(global,viewer);

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