Actual source code: ex3.c

  2: static char help[] = "Solves the 1-dimensional wave equation.\n\n";

 4:  #include petscda.h
 5:  #include petscsys.h

  9: int main(int argc,char **argv)
 10: {
 11:   PetscMPIInt    rank,size;
 13:   PetscInt       M = 60,time_steps = 100, localsize,j,i,mybase,myend,width,xbase,*localnodes = PETSC_NULL;
 14:   DA             da;
 15:   PetscViewer    viewer,viewer_private;
 16:   PetscDraw      draw;
 17:   Vec            local,global,copy;
 18:   PetscScalar    *localptr,*copyptr;
 19:   PetscReal      a,h,k;
 20:   PetscTruth     flg;
 21: 
 22:   PetscInitialize(&argc,&argv,(char*)0,help);
 23:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 24:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 26:   PetscOptionsGetInt(PETSC_NULL,"-M",&M,PETSC_NULL);
 27:   PetscOptionsGetInt(PETSC_NULL,"-time",&time_steps,PETSC_NULL);
 28:   /*
 29:       Test putting two nodes on each processor, exact last processor gets the rest
 30:   */
 31:   PetscOptionsHasName(PETSC_NULL,"-distribute",&flg);
 32:   if (flg) {
 33:     PetscMalloc(size*sizeof(PetscInt),&localnodes);
 34:     for (i=0; i<size-1; i++) { localnodes[i] = 2;}
 35:     localnodes[size-1] = M - 2*(size-1);
 36:   }
 37: 
 38:   /* Set up the array */
 39:   DACreate1d(PETSC_COMM_WORLD,DA_XPERIODIC,M,1,1,localnodes,&da);
 40:   PetscFree(localnodes);
 41:   DACreateGlobalVector(da,&global);
 42:   DACreateLocalVector(da,&local);

 44:   /* Set up display to show combined wave graph */
 45:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"Entire Solution",20,480,800,200,&viewer);
 46:   PetscViewerDrawGetDraw(viewer,0,&draw);
 47:   PetscDrawSetDoubleBuffer(draw);

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

 52:   /* set up display to show my portion of the wave */
 53:   xbase = (int)((mybase)*((800.0 - 4.0*size)/M) + 4.0*rank);
 54:   width = (int)((myend-mybase)*800./M);
 55:   PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Local Portion of Solution",xbase,200,
 56:                          width,200,&viewer_private);
 57:   PetscViewerDrawGetDraw(viewer_private,0,&draw);
 58:   PetscDrawSetDoubleBuffer(draw);



 62:   /* Initialize the array */
 63:   VecGetLocalSize(local,&localsize);
 64:   VecGetArray(local,&localptr);
 65:   localptr[0] = 0.0;
 66:   localptr[localsize-1] = 0.0;
 67:   for (i=1; i<localsize-1; i++) {
 68:     j=(i-1)+mybase;
 69:     localptr[i] = sin((PETSC_PI*j*6)/((PetscReal)M)
 70:                         + 1.2 * sin((PETSC_PI*j*2)/((PetscReal)M))) * 2;
 71:   }

 73:   VecRestoreArray(local,&localptr);
 74:   DALocalToGlobal(da,local,INSERT_VALUES,global);

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

 79:   /* Assign Parameters */
 80:   a= 1.0;
 81:   h= 1.0/M;
 82:   k= h;

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

 86:     /* Global to Local */
 87:     DAGlobalToLocalBegin(da,global,INSERT_VALUES,local);
 88:     DAGlobalToLocalEnd(da,global,INSERT_VALUES,local);

 90:     /*Extract local array */
 91:     VecGetArray(local,&localptr);
 92:     VecGetArray(copy,&copyptr);

 94:     /* Update Locally - Make array of new values */
 95:     /* Note: I don't do anything for the first and last entry */
 96:     for (i=1; i< localsize-1; i++) {
 97:       copyptr[i] = .5*(localptr[i+1]+localptr[i-1]) -
 98:                     (k / (2.0*a*h)) * (localptr[i+1] - localptr[i-1]);
 99:     }
100:     VecRestoreArray(copy,&copyptr);
101:     VecRestoreArray(local,&localptr);

103:     /* Local to Global */
104:     DALocalToGlobal(da,copy,INSERT_VALUES,global);
105: 
106:     /* View my part of Wave */
107:     VecView(copy,viewer_private);

109:     /* View global Wave */
110:     VecView(global,viewer);
111:   }

113:   DADestroy(da);
114:   PetscViewerDestroy(viewer);
115:   PetscViewerDestroy(viewer_private);
116:   VecDestroy(copy);
117:   VecDestroy(local);
118:   VecDestroy(global);

120:   PetscFinalize();
121:   return 0;
122: }
123: