Actual source code: ex3.c
petsc-3.5.4 2015-05-23
2: static char help[] = "Solves the 1-dimensional wave equation.\n\n";
4: #include <petscdm.h>
5: #include <petscdmda.h>
6: #include <petscdraw.h>
10: int main(int argc,char **argv)
11: {
12: PetscMPIInt rank,size;
14: PetscInt M = 60,time_steps = 100, localsize,j,i,mybase,myend,width,xbase,*localnodes = NULL;
15: DM da;
16: PetscViewer viewer,viewer_private;
17: PetscDraw draw;
18: Vec local,global,copy;
19: PetscScalar *localptr,*copyptr;
20: PetscReal a,h,k;
21: PetscBool flg = PETSC_FALSE;
23: PetscInitialize(&argc,&argv,(char*)0,help);
24: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
25: MPI_Comm_size(PETSC_COMM_WORLD,&size);
27: PetscOptionsGetInt(NULL,"-M",&M,NULL);
28: PetscOptionsGetInt(NULL,"-time",&time_steps,NULL);
29: /*
30: Test putting two nodes on each processor, exact last processor gets the rest
31: */
32: PetscOptionsGetBool(NULL,"-distribute",&flg,NULL);
33: if (flg) {
34: PetscMalloc1(size,&localnodes);
35: for (i=0; i<size-1; i++) localnodes[i] = 2;
36: localnodes[size-1] = M - 2*(size-1);
37: }
39: /* Set up the array */
40: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC,M,1,1,localnodes,&da);
41: PetscFree(localnodes);
42: DMCreateGlobalVector(da,&global);
43: DMCreateLocalVector(da,&local);
45: /* Set up display to show combined wave graph */
46: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"Entire Solution",20,480,800,200,&viewer);
47: PetscViewerDrawGetDraw(viewer,0,&draw);
48: PetscDrawSetDoubleBuffer(draw);
50: /* determine starting point of each processor */
51: VecGetOwnershipRange(global,&mybase,&myend);
53: /* set up display to show my portion of the wave */
54: xbase = (int)((mybase)*((800.0 - 4.0*size)/M) + 4.0*rank);
55: width = (int)((myend-mybase)*800./M);
56: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Local Portion of Solution",xbase,200,
57: width,200,&viewer_private);
58: PetscViewerDrawGetDraw(viewer_private,0,&draw);
59: PetscDrawSetDoubleBuffer(draw);
63: /* Initialize the array */
64: VecGetLocalSize(local,&localsize);
65: VecGetArray(local,&localptr);
67: localptr[0] = 0.0;
68: localptr[localsize-1] = 0.0;
69: for (i=1; i<localsize-1; i++) {
70: j = (i-1)+mybase;
71: localptr[i] = PetscSinReal((PETSC_PI*j*6)/((PetscReal)M) + 1.2 * PetscSinReal((PETSC_PI*j*2)/((PetscReal)M))) * 2;
72: }
74: VecRestoreArray(local,&localptr);
75: DMLocalToGlobalBegin(da,local,INSERT_VALUES,global);
76: DMLocalToGlobalEnd(da,local,INSERT_VALUES,global);
78: /* Make copy of local array for doing updates */
79: VecDuplicate(local,©);
81: /* Assign Parameters */
82: a= 1.0;
83: h= 1.0/M;
84: k= h;
86: for (j=0; j<time_steps; j++) {
88: /* Global to Local */
89: DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);
90: DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);
92: /*Extract local array */
93: VecGetArray(local,&localptr);
94: VecGetArray(copy,©ptr);
96: /* Update Locally - Make array of new values */
97: /* Note: I don't do anything for the first and last entry */
98: for (i=1; i< localsize-1; i++) {
99: copyptr[i] = .5*(localptr[i+1]+localptr[i-1]) - (k / (2.0*a*h)) * (localptr[i+1] - localptr[i-1]);
100: }
101: VecRestoreArray(copy,©ptr);
102: VecRestoreArray(local,&localptr);
104: /* Local to Global */
105: DMLocalToGlobalBegin(da,copy,INSERT_VALUES,global);
106: DMLocalToGlobalEnd(da,copy,INSERT_VALUES,global);
108: /* View my part of Wave */
109: VecView(copy,viewer_private);
111: /* View global Wave */
112: VecView(global,viewer);
113: }
115: DMDestroy(&da);
116: PetscViewerDestroy(&viewer);
117: PetscViewerDestroy(&viewer_private);
118: VecDestroy(©);
119: VecDestroy(&local);
120: VecDestroy(&global);
122: PetscFinalize();
123: return 0;
124: }