Actual source code: ex5.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  2: /* This file created by Peter Mell   6/30/95 */

  4: static char help[] = "Solves the one dimensional heat equation.\n\n";

  6: #include <petscdm.h>
  7: #include <petscdmda.h>
  8: #include <petscdraw.h>

 12: int main(int argc,char **argv)
 13: {
 14:   PetscMPIInt    rank,size;
 16:   PetscInt       M = 14,time_steps = 1000,w=1,s=1,localsize,j,i,mybase,myend;
 17:   DM             da;
 18:   PetscViewer    viewer;
 19:   PetscDraw      draw;
 20:   Vec            local,global,copy;
 21:   PetscScalar    *localptr,*copyptr;
 22:   PetscReal      h,k;

 24:   PetscInitialize(&argc,&argv,(char*)0,help);

 26:   PetscOptionsGetInt(NULL,"-M",&M,NULL);
 27:   PetscOptionsGetInt(NULL,"-time",&time_steps,NULL);

 29:   /* Set up the array */
 30:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,M,w,s,NULL,&da);
 31:   DMCreateGlobalVector(da,&global);
 32:   DMCreateLocalVector(da,&local);
 33:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 34:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

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

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

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

 47:   /* Initialize the Array */
 48:   VecGetLocalSize (local,&localsize);
 49:   VecGetArray (local,&localptr);
 50:   VecGetArray (copy,&copyptr);

 52:   localptr[0] = copyptr[0] = 0.0;

 54:   localptr[localsize-1] = copyptr[localsize-1] = 1.0;
 55:   for (i=1; i<localsize-1; i++) {
 56:     j = (i-1) + mybase;

 58:     localptr[i] = PetscSinReal((PETSC_PI*j*6)/((PetscReal)M) + 1.2 * PetscSinReal((PETSC_PI*j*2)/((PetscReal)M))) * 4+4;
 59:   }

 61:   VecRestoreArray(local,&localptr);
 62:   VecRestoreArray(copy,&copyptr);
 63:   DMLocalToGlobalBegin(da,local,INSERT_VALUES,global);
 64:   DMLocalToGlobalEnd(da,local,INSERT_VALUES,global);

 66:   /* Assign Parameters */
 67:   h= 1.0/M;
 68:   k= h*h/2.2;

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

 72:     /* Global to Local */
 73:     DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);
 74:     DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);

 76:     /*Extract local array */
 77:     VecGetArray(local,&localptr);
 78:     VecGetArray (copy,&copyptr);

 80:     /* Update Locally - Make array of new values */
 81:     /* Note: I don't do anything for the first and last entry */
 82:     for (i=1; i< localsize-1; i++) {
 83:       copyptr[i] = localptr[i] + (k/(h*h)) * (localptr[i+1]-2.0*localptr[i]+localptr[i-1]);
 84:     }

 86:     VecRestoreArray(copy,&copyptr);
 87:     VecRestoreArray(local,&localptr);

 89:     /* Local to Global */
 90:     DMLocalToGlobalBegin(da,copy,INSERT_VALUES,global);
 91:     DMLocalToGlobalEnd(da,copy,INSERT_VALUES,global);

 93:     /* View Wave */
 94:     VecView(global,viewer);

 96:   }

 98:   PetscViewerDestroy(&viewer);
 99:   VecDestroy(&copy);
100:   VecDestroy(&local);
101:   VecDestroy(&global);
102:   DMDestroy(&da);
103:   PetscFinalize();
104:   return 0;
105: }