Actual source code: ex5.c

  2: static char help[] = "Tests DAGetElements() and VecView() 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       M = 10,N = 8,m = PETSC_DECIDE,n = PETSC_DECIDE,ne,i;
 12:   const PetscInt *e;
 14:   PetscTruth     flg;
 15:   DA             da;
 16:   PetscViewer    viewer;
 17:   Vec            local,global;
 18:   PetscScalar    value;
 19:   DAPeriodicType ptype = DA_NONPERIODIC;
 20:   DAStencilType  stype = DA_STENCIL_BOX;
 21:   PetscScalar    *lv;

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

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

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

 39:   value = -3.0;
 40:   VecSet(global,value);
 41:   DAGlobalToLocalBegin(da,global,INSERT_VALUES,local);
 42:   DAGlobalToLocalEnd(da,global,INSERT_VALUES,local);

 44:   DAGetElements(da,&ne,&e);
 45:   VecGetArray(local,&lv);
 46:   for (i=0; i<ne; i++) {
 47:     PetscPrintf(PETSC_COMM_WORLD,"i %D e[3*i] %D %D %D\n",i,e[3*i],e[3*i+1],e[3*i+2]);
 48:     lv[e[3*i]] = i;
 49:   }
 50:   VecRestoreArray(local,&lv);
 51:   DARestoreElements(da,&ne,&e);

 53:   DALocalToGlobal(da,local,ADD_VALUES,global);

 55:   DAView(da,viewer);
 56:   VecView(global,viewer);

 58:   /* Free memory */
 59:   PetscViewerDestroy(viewer);
 60:   VecDestroy(local);
 61:   VecDestroy(global);
 62:   DADestroy(da);
 63:   PetscFinalize();
 64:   return 0;
 65: }
 66: