Actual source code: ex12.c

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

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

 9:  #include petscda.h

 13: int main(int argc,char **argv)
 14: {
 15:   PetscMPIInt    rank,size;
 16:   PetscInt       M = 14,time_steps = 20,w=1,s=1,localsize,j,i,mybase,myend;
 18:   DA             da;
 19:   Vec            local,global,copy;
 20:   PetscScalar    *localptr,*copyptr;
 21:   PetscReal      h,k;
 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);


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

 43:   /* Initialize the Array */
 44:   VecGetLocalSize (local,&localsize);
 45:   VecGetArray (local,&localptr);
 46:   localptr[0] = copyptr[0] = 0.0;
 47:   localptr[localsize-1] = copyptr[localsize-1] = 1.0;
 48:   for (i=1; i<localsize-1; i++) {
 49:     j=(i-1)+mybase;
 50:     localptr[i] = sin((PETSC_PI*j*6)/((PetscReal)M)
 51:                         + 1.2 * sin((PETSC_PI*j*2)/((PetscReal)M))) * 4+4;
 52:   }

 54:   VecRestoreArray (copy,&copyptr);
 55:   VecRestoreArray(local,&localptr);
 56:   DALocalToGlobal(da,local,INSERT_VALUES,global);

 58:   /* Assign Parameters */
 59:   h= 1.0/M;
 60:   k= h*h/2.2;

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

 64:     /* Global to Local */
 65:     DAGlobalToLocalBegin(da,global,INSERT_VALUES,local);
 66:     DAGlobalToLocalEnd(da,global,INSERT_VALUES,local);

 68:     /*Extract local array */
 69:     VecGetArray(local,&localptr);
 70:     VecGetArray (copy,&copyptr);

 72:     /* Update Locally - Make array of new values */
 73:     /* Note: I don't do anything for the first and last entry */
 74:     for (i=1; i< localsize-1; i++) {
 75:       copyptr[i] = localptr[i] + (k/(h*h)) *
 76:                            (localptr[i+1]-2.0*localptr[i]+localptr[i-1]);
 77:     }
 78: 
 79:     VecRestoreArray (copy,&copyptr);
 80:     VecRestoreArray(local,&localptr);

 82:     /* Local to Global */
 83:     DALocalToGlobal(da,copy,INSERT_VALUES,global);
 84: 
 85:     /* View Wave */
 86:   /* Set Up Display to Show Heat Graph */
 87: #if defined(PETSC_USE_SOCKET_VIEWER)
 88:     VecView(global,PETSC_VIEWER_SOCKET_WORLD);
 89: #endif
 90:   }

 92:   VecDestroy(copy);
 93:   VecDestroy(local);
 94:   VecDestroy(global);
 95:   DADestroy(da);
 96:   PetscFinalize();
 97:   return 0;
 98: }
 99: