Actual source code: ex5.c

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

 4:  #include petscda.h

  8: int main(int argc,char **argv)
  9: {
 10:   PetscInt       M = 10,N = 8,m = PETSC_DECIDE,n = PETSC_DECIDE,ne,i;
 11:   const PetscInt *e;
 13:   PetscTruth     flg = PETSC_FALSE;
 14:   DA             da;
 15:   PetscViewer    viewer;
 16:   Vec            local,global;
 17:   PetscScalar    value;
 18:   DAPeriodicType ptype = DA_NONPERIODIC;
 19:   DAStencilType  stype = DA_STENCIL_BOX;
 20:   PetscScalar    *lv;

 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:   PetscOptionsGetTruth(PETSC_NULL,"-star_stencil",&flg,PETSC_NULL);
 31:   if (flg) stype = DA_STENCIL_STAR;

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

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

 43:   DAGetElements(da,&ne,&e);
 44:   VecGetArray(local,&lv);
 45:   for (i=0; i<ne; i++) {
 46:     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]);
 47:     lv[e[3*i]] = i;
 48:   }
 49:   VecRestoreArray(local,&lv);
 50:   DARestoreElements(da,&ne,&e);

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

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

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