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.nn";

 7:  #include petscda.h
 8:  #include petscsys.h

 10: int main(int argc,char **argv)
 11: {
 12:   int         rank,size,M = 14,ierr,time_steps = 1000,w=1,s=1;
 13:   DA          da;
 14:   PetscViewer viewer;
 15:   PetscDraw   draw;
 16:   Vec         local,global,copy;
 17:   PetscScalar *localptr,*copyptr;
 18:   PetscReal   h,k;
 19:   int         localsize,j,i,mybase,myend;
 20: 
 21:   PetscInitialize(&argc,&argv,(char*)0,help);

 23:   PetscOptionsGetInt(PETSC_NULL,"-M",&M,PETSC_NULL);
 24:   PetscOptionsGetInt(PETSC_NULL,"-time",&time_steps,PETSC_NULL);
 25: 
 26:   /* Set up the array */
 27:   DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,M,w,s,PETSC_NULL,&da);
 28:   DACreateGlobalVector(da,&global);
 29:   DACreateLocalVector(da,&local);
 30:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 31:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 33:   /* Make copy of local array for doing updates */
 34:   VecDuplicate(local,&copy);

 36:   /* Set Up Display to Show Heat Graph */
 37:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",80,480,500,160,&viewer);
 38:   PetscViewerDrawGetDraw(viewer,0,&draw);
 39:   PetscDrawSetDoubleBuffer(draw);

 41:   /* determine starting point of each processor */
 42:   VecGetOwnershipRange(global,&mybase,&myend);

 44:   /* Initialize the Array */
 45:   VecGetLocalSize (local,&localsize);
 46:   VecGetArray (local,&localptr);
 47:   VecGetArray (copy,&copyptr);
 48:   localptr[0] = copyptr[0] = 0.0;
 49:   localptr[localsize-1] = copyptr[localsize-1] = 1.0;
 50:   for (i=1; i<localsize-1; i++) {
 51:     j=(i-1)+mybase;
 52:     localptr[i] = sin((PETSC_PI*j*6)/((PetscReal)M)
 53:                         + 1.2 * sin((PETSC_PI*j*2)/((PetscReal)M))) * 4+4;
 54:   }

 56:   VecRestoreArray(local,&localptr);
 57:   VecRestoreArray(copy,&copyptr);
 58:   DALocalToGlobal(da,local,INSERT_VALUES,global);

 60:   /* Assign Parameters */
 61:   h= 1.0/M;
 62:   k= h*h/2.2;

 64:   for (j=0; j<time_steps; j++) {

 66:     /* Global to Local */
 67:     DAGlobalToLocalBegin(da,global,INSERT_VALUES,local);
 68:     DAGlobalToLocalEnd(da,global,INSERT_VALUES,local);

 70:     /*Extract local array */
 71:     VecGetArray(local,&localptr);
 72:     VecGetArray (copy,&copyptr);

 74:     /* Update Locally - Make array of new values */
 75:     /* Note: I don't do anything for the first and last entry */
 76:     for (i=1; i< localsize-1; i++) {
 77:       copyptr[i] = localptr[i] + (k/(h*h)) *
 78:                            (localptr[i+1]-2.0*localptr[i]+localptr[i-1]);
 79:     }
 80: 
 81:     VecRestoreArray(copy,&copyptr);
 82:     VecRestoreArray(local,&localptr);

 84:     /* Local to Global */
 85:     DALocalToGlobal(da,copy,INSERT_VALUES,global);
 86: 
 87:     /* View Wave */
 88:     VecView(global,viewer);

 90:   }

 92:   PetscViewerDestroy(viewer);
 93:   VecDestroy(copy);
 94:   VecDestroy(local);
 95:   VecDestroy(global);
 96:   DADestroy(da);
 97:   PetscFinalize();
 98:   return 0;
 99: }
100: