Actual source code: ex2.c

  1: /*$Id: ex2.c,v 1.15 2001/08/07 03:04:45 balay Exp $*/

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

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

  8: #undef __FUNCT__
 10: int main(int argc,char **argv)
 11: {
 12:   int            i,j,rank,M = 10,N = 8,m = PETSC_DECIDE,n = PETSC_DECIDE,ierr;
 13:   PetscTruth     flg;
 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;

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

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

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

 39:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 40:   value = 5.0*rank;
 41:   VecSet(&value,global);

 43:   VecView(global,viewer);

 45:   /*
 46:      Create Scatter from global DA parallel vector to local vector that
 47:    contains all entries
 48:   */
 49:   DAGlobalToNaturalAllCreate(da,&tolocalall);
 50:   DANaturalAllToGlobalCreate(da,&fromlocalall);

 52:   VecScatterBegin(global,localall,INSERT_VALUES,SCATTER_FORWARD,tolocalall);
 53:   VecScatterEnd(global,localall,INSERT_VALUES,SCATTER_FORWARD,tolocalall);

 55:   VecGetArray(localall,&vlocal);
 56:   for (j=0; j<N; j++) {
 57:     for (i=0; i<M; i++) {
 58:       *vlocal++ += i + j*M;
 59:     }
 60:   }
 61:   VecRestoreArray(localall,&vlocal);

 63:   /* scatter back to global vector */
 64:   VecScatterBegin(localall,global,INSERT_VALUES,SCATTER_FORWARD,fromlocalall);
 65:   VecScatterEnd(localall,global,INSERT_VALUES,SCATTER_FORWARD,fromlocalall);

 67:   VecView(global,viewer);

 69:   /* Free memory */
 70:   VecScatterDestroy(tolocalall);
 71:   VecScatterDestroy(fromlocalall);
 72:   PetscViewerDestroy(viewer);
 73:   VecDestroy(localall);
 74:   VecDestroy(global);
 75:   DADestroy(da);
 76:   PetscFinalize();
 77:   return 0;
 78: }
 79: