Actual source code: ex3.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  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;
 19:   PetscScalar    *localptr,*globalptr;
 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,width,200,&viewer_private);
 57:   PetscViewerDrawGetDraw(viewer_private,0,&draw);
 58:   PetscDrawSetDoubleBuffer(draw);



 62:   /* Initialize the array */
 63:   VecGetLocalSize(local,&localsize);
 64:   VecGetArray(global,&globalptr);

 66:   for (i=1; i<localsize-1; i++) {
 67:     j           = (i-1)+mybase;
 68:     globalptr[i-1] = PetscSinReal((PETSC_PI*j*6)/((PetscReal)M) + 1.2 * PetscSinReal((PETSC_PI*j*2)/((PetscReal)M))) * 2;
 69:   }

 71:   VecRestoreArray(global,&globalptr);

 73:   /* Assign Parameters */
 74:   a= 1.0;
 75:   h= 1.0/M;
 76:   k= h;

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

 80:     /* Global to Local */
 81:     DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);
 82:     DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);

 84:     /*Extract local array */
 85:     VecGetArray(local,&localptr);
 86:     VecGetArray(global,&globalptr);

 88:     /* Update Locally - Make array of new values */
 89:     /* Note: I don't do anything for the first and last entry */
 90:     for (i=1; i< localsize-1; i++) {
 91:       globalptr[i-1] = .5*(localptr[i+1]+localptr[i-1]) - (k / (2.0*a*h)) * (localptr[i+1] - localptr[i-1]);
 92:     }
 93:     VecRestoreArray(global,&globalptr);
 94:     VecRestoreArray(local,&localptr);

 96:     /* View my part of Wave */
 97:     VecView(global,viewer_private);

 99:     /* View global Wave */
100:     VecView(global,viewer);
101:   }

103:   DMDestroy(&da);
104:   PetscViewerDestroy(&viewer);
105:   PetscViewerDestroy(&viewer_private);
106:   VecDestroy(&local);
107:   VecDestroy(&global);

109:   PetscFinalize();
110:   return 0;
111: }