Actual source code: ex12.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  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 <petscdm.h>
 10: #include <petscdmda.h>

 14: int main(int argc,char **argv)
 15: {
 16:   PetscMPIInt    rank,size;
 17:   PetscInt       M = 14,time_steps = 20,w=1,s=1,localsize,j,i,mybase,myend,globalsize;
 19:   DM             da;
 20:   Vec            global,local;
 21:   PetscScalar    *globalptr,*localptr;
 22:   PetscReal      h,k;

 24:   PetscInitialize(&argc,&argv,(char*)0,help);

 26:   PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);
 27:   PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL);

 29:   /* Set up the array */
 30:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,M,w,s,NULL,&da);
 31:   DMCreateGlobalVector(da,&global);
 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:   DMCreateLocalVector(da,&local);


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

 42:   /* Initialize the Array */
 43:   VecGetLocalSize (global,&globalsize);
 44:   VecGetArray (global,&globalptr);


 47:   for (i=0; i<globalsize; i++) {
 48:     j = i + mybase;

 50:     globalptr[i] = PetscSinReal((PETSC_PI*j*6)/((PetscReal)M) + 1.2 * PetscSinReal((PETSC_PI*j*2)/((PetscReal)M))) * 4+4;
 51:   }

 53:   VecRestoreArray(global,&localptr);

 55:   /* Assign Parameters */
 56:   h= 1.0/M;
 57:   k= h*h/2.2;
 58:   VecGetLocalSize(local,&localsize);

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

 62:     /* Global to Local */
 63:     DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);
 64:     DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);

 66:     /*Extract local array */
 67:     VecGetArray(local,&localptr);
 68:     VecGetArray (global,&globalptr);

 70:     /* Update Locally - Make array of new values */
 71:     /* Note: I don't do anything for the first and last entry */
 72:     for (i=1; i< localsize-1; i++) {
 73:       globalptr[i-1] = localptr[i] + (k/(h*h)) * (localptr[i+1]-2.0*localptr[i]+localptr[i-1]);
 74:     }

 76:     VecRestoreArray (global,&globalptr);
 77:     VecRestoreArray(local,&localptr);

 79:     /* View Wave */
 80:     /* Set Up Display to Show Heat Graph */
 81: #if defined(PETSC_USE_SOCKET_VIEWER)
 82:     VecView(global,PETSC_VIEWER_SOCKET_WORLD);
 83: #endif
 84:   }

 86:   VecDestroy(&local);
 87:   VecDestroy(&global);
 88:   DMDestroy(&da);
 89:   PetscFinalize();
 90:   return 0;
 91: }