Actual source code: ex3.c

  1: /*$Id: ex3.c,v 1.49 2001/08/07 21:31:49 bsmith Exp $*/

  3: static char help[] = "Solves the 1-dimensional wave equation.nn";

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

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

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

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

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

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



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

 71:   VecRestoreArray(local,&localptr);
 72:   DALocalToGlobal(da,local,INSERT_VALUES,global);

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

 77:   /* Assign Parameters */
 78:   a= 1.0;
 79:   h= 1.0/M;
 80:   k= h;

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

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

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

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

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

107:     /* View global Wave */
108:     VecView(global,viewer);
109:   }

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

118:   PetscFinalize();
119:   return 0;
120: }
121: