Actual source code: ex5.c
petsc-3.4.2 2013-07-02
2: /* This file created by Peter Mell 6/30/95 */
4: static char help[] = "Solves the one dimensional heat equation.\n\n";
6: #include <petscdmda.h>
7: #include <petscdraw.h>
11: int main(int argc,char **argv)
12: {
13: PetscMPIInt rank,size;
15: PetscInt M = 14,time_steps = 1000,w=1,s=1,localsize,j,i,mybase,myend;
16: DM da;
17: PetscViewer viewer;
18: PetscDraw draw;
19: Vec local,global,copy;
20: PetscScalar *localptr,*copyptr;
21: PetscReal h,k;
23: PetscInitialize(&argc,&argv,(char*)0,help);
25: PetscOptionsGetInt(NULL,"-M",&M,NULL);
26: PetscOptionsGetInt(NULL,"-time",&time_steps,NULL);
28: /* Set up the array */
29: DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,M,w,s,NULL,&da);
30: DMCreateGlobalVector(da,&global);
31: DMCreateLocalVector(da,&local);
32: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
33: MPI_Comm_size(PETSC_COMM_WORLD,&size);
35: /* Make copy of local array for doing updates */
36: VecDuplicate(local,©);
38: /* Set Up Display to Show Heat Graph */
39: PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",80,480,500,160,&viewer);
40: PetscViewerDrawGetDraw(viewer,0,&draw);
41: PetscDrawSetDoubleBuffer(draw);
43: /* determine starting point of each processor */
44: VecGetOwnershipRange(global,&mybase,&myend);
46: /* Initialize the Array */
47: VecGetLocalSize (local,&localsize);
48: VecGetArray (local,&localptr);
49: VecGetArray (copy,©ptr);
51: localptr[0] = copyptr[0] = 0.0;
53: localptr[localsize-1] = copyptr[localsize-1] = 1.0;
54: for (i=1; i<localsize-1; i++) {
55: j = (i-1) + mybase;
57: localptr[i] = sin((PETSC_PI*j*6)/((PetscReal)M) + 1.2 * sin((PETSC_PI*j*2)/((PetscReal)M))) * 4+4;
58: }
60: VecRestoreArray(local,&localptr);
61: VecRestoreArray(copy,©ptr);
62: DMLocalToGlobalBegin(da,local,INSERT_VALUES,global);
63: DMLocalToGlobalEnd(da,local,INSERT_VALUES,global);
65: /* Assign Parameters */
66: h= 1.0/M;
67: k= h*h/2.2;
69: for (j=0; j<time_steps; j++) {
71: /* Global to Local */
72: DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);
73: DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);
75: /*Extract local array */
76: VecGetArray(local,&localptr);
77: VecGetArray (copy,©ptr);
79: /* Update Locally - Make array of new values */
80: /* Note: I don't do anything for the first and last entry */
81: for (i=1; i< localsize-1; i++) {
82: copyptr[i] = localptr[i] + (k/(h*h)) * (localptr[i+1]-2.0*localptr[i]+localptr[i-1]);
83: }
85: VecRestoreArray(copy,©ptr);
86: VecRestoreArray(local,&localptr);
88: /* Local to Global */
89: DMLocalToGlobalBegin(da,copy,INSERT_VALUES,global);
90: DMLocalToGlobalEnd(da,copy,INSERT_VALUES,global);
92: /* View Wave */
93: VecView(global,viewer);
95: }
97: PetscViewerDestroy(&viewer);
98: VecDestroy(©);
99: VecDestroy(&local);
100: VecDestroy(&global);
101: DMDestroy(&da);
102: PetscFinalize();
103: return 0;
104: }