Actual source code: ex5.c
1: /*$Id: ex5.c,v 1.44 2001/08/07 21:31:49 bsmith Exp $*/
3: /* This file created by Peter Mell 6/30/95 */
5: static char help[] = "Solves the one dimensional heat equation.\n\n";
7: #include petscda.h
8: #include petscsys.h
12: int main(int argc,char **argv)
13: {
14: int rank,size,M = 14,ierr,time_steps = 1000,w=1,s=1;
15: DA da;
16: PetscViewer viewer;
17: PetscDraw draw;
18: Vec local,global,copy;
19: PetscScalar *localptr,*copyptr;
20: PetscReal h,k;
21: int localsize,j,i,mybase,myend;
22:
23: PetscInitialize(&argc,&argv,(char*)0,help);
25: PetscOptionsGetInt(PETSC_NULL,"-M",&M,PETSC_NULL);
26: PetscOptionsGetInt(PETSC_NULL,"-time",&time_steps,PETSC_NULL);
27:
28: /* Set up the array */
29: DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,M,w,s,PETSC_NULL,&da);
30: DACreateGlobalVector(da,&global);
31: DACreateLocalVector(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);
50: localptr[0] = copyptr[0] = 0.0;
51: localptr[localsize-1] = copyptr[localsize-1] = 1.0;
52: for (i=1; i<localsize-1; i++) {
53: j=(i-1)+mybase;
54: localptr[i] = sin((PETSC_PI*j*6)/((PetscReal)M)
55: + 1.2 * sin((PETSC_PI*j*2)/((PetscReal)M))) * 4+4;
56: }
58: VecRestoreArray(local,&localptr);
59: VecRestoreArray(copy,©ptr);
60: DALocalToGlobal(da,local,INSERT_VALUES,global);
62: /* Assign Parameters */
63: h= 1.0/M;
64: k= h*h/2.2;
66: for (j=0; j<time_steps; j++) {
68: /* Global to Local */
69: DAGlobalToLocalBegin(da,global,INSERT_VALUES,local);
70: DAGlobalToLocalEnd(da,global,INSERT_VALUES,local);
72: /*Extract local array */
73: VecGetArray(local,&localptr);
74: VecGetArray (copy,©ptr);
76: /* Update Locally - Make array of new values */
77: /* Note: I don't do anything for the first and last entry */
78: for (i=1; i< localsize-1; i++) {
79: copyptr[i] = localptr[i] + (k/(h*h)) *
80: (localptr[i+1]-2.0*localptr[i]+localptr[i-1]);
81: }
82:
83: VecRestoreArray(copy,©ptr);
84: VecRestoreArray(local,&localptr);
86: /* Local to Global */
87: DALocalToGlobal(da,copy,INSERT_VALUES,global);
88:
89: /* View Wave */
90: VecView(global,viewer);
92: }
94: PetscViewerDestroy(viewer);
95: VecDestroy(copy);
96: VecDestroy(local);
97: VecDestroy(global);
98: DADestroy(da);
99: PetscFinalize();
100: return 0;
101: }
102: