Actual source code: ex1.c
2: static char help[] = "Tests 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: PetscMPIInt rank;
12: PetscInt M = -10,N = -8;
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: #if defined(PETSC_HAVE_MATLAB_ENGINE)
22: PetscViewer mviewer;
23: #endif
25: PetscInitialize(&argc,&argv,(char*)0,help);
26: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",300,0,300,300,&viewer);
27: #if defined(PETSC_HAVE_MATLAB_ENGINE)
28: PetscViewerMatlabOpen(PETSC_COMM_WORLD,"tmp.mat",FILE_MODE_WRITE,&mviewer);
29: #endif
31: PetscOptionsHasName(PETSC_NULL,"-star_stencil",&flg);
32: if (flg) stype = DA_STENCIL_STAR;
34: /* Use the options
35: -da_grid_x <nx> - number of grid points in x direction, if M < 0
36: -da_grid_y <ny> - number of grid points in y direction, if N < 0
37: -da_processors_x <MX> number of processors in x directio
38: -da_processors_y <MY> number of processors in x direction
39: */
40: /* Create distributed array and get vectors */
41: DACreate2d(PETSC_COMM_WORLD,ptype,stype,M,N,PETSC_DECIDE,PETSC_DECIDE,1,1,PETSC_NULL,PETSC_NULL,&da);
42: DACreateGlobalVector(da,&global);
43: DACreateLocalVector(da,&local);
45: value = -3.0;
46: VecSet(global,value);
47: DAGlobalToLocalBegin(da,global,INSERT_VALUES,local);
48: DAGlobalToLocalEnd(da,global,INSERT_VALUES,local);
50: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
51: value = rank+1;
52: VecScale(local,value);
53: DALocalToGlobal(da,local,ADD_VALUES,global);
55: DAView(da,viewer);
56: VecView(global,viewer);
57: #if defined(PETSC_HAVE_MATLAB_ENGINE)
58: DAView(da,mviewer);
59: VecView(global,mviewer);
60: #endif
62: /* Free memory */
63: #if defined(PETSC_HAVE_MATLAB_ENGINE)
64: PetscViewerDestroy(mviewer);
65: #endif
66: PetscViewerDestroy(viewer);
67: VecDestroy(local);
68: VecDestroy(global);
69: DADestroy(da);
70: PetscFinalize();
71: return 0;
72: }
73: