Actual source code: ex12.c

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

  3: /*
  4:    Simple example to show how PETSc programs can be run from Matlab. 
  5:   See ex12.m.
  6: */

  8: static char help[] = "Solves the one dimensional heat equation.nn";

 10:  #include petscda.h
 11:  #include petscsys.h

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

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

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

 39:   /* Set Up Display to Show Heat Graph */
 40:   viewer = PETSC_VIEWER_SOCKET_WORLD;

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

 45:   /* Initialize the Array */
 46:   VecGetLocalSize (local,&localsize);
 47:   VecGetArray (local,&localptr);
 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 (copy,&copyptr);
 57:   VecRestoreArray(local,&localptr);
 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:   VecDestroy(copy);
 93:   VecDestroy(local);
 94:   VecDestroy(global);
 95:   DADestroy(da);
 96:   PetscFinalize();
 97:   return 0;
 98: }
 99: