Actual source code: ex12.c

petsc-3.5.4 2015-05-23
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;
 19:   DM             da;
 20:   Vec            local,global,copy;
 21:   PetscScalar    *localptr,*copyptr;
 22:   PetscReal      h,k;

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

 26:   PetscOptionsGetInt(NULL,"-M",&M,NULL);
 27:   PetscOptionsGetInt(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:   DMCreateLocalVector(da,&local);
 33:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
 34:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

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


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

 44:   /* Initialize the Array */
 45:   VecGetLocalSize (local,&localsize);
 46:   VecGetArray (local,&localptr);

 48:   localptr[0] = copyptr[0] = 0.0;

 50:   localptr[localsize-1] = copyptr[localsize-1] = 1.0;
 51:   for (i=1; i<localsize-1; i++) {
 52:     j = (i-1) + mybase;

 54:     localptr[i] = PetscSinReal((PETSC_PI*j*6)/((PetscReal)M) + 1.2 * PetscSinReal((PETSC_PI*j*2)/((PetscReal)M))) * 4+4;
 55:   }

 57:   VecRestoreArray (copy,&copyptr);
 58:   VecRestoreArray(local,&localptr);
 59:   DMLocalToGlobalBegin(da,local,INSERT_VALUES,global);
 60:   DMLocalToGlobalEnd(da,local,INSERT_VALUES,global);

 62:   /* Assign Parameters */
 63:   h= 1.0/M;
 64:   k= h*h/2.2;

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

 68:     /* Global to Local */
 69:     DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);
 70:     DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);

 72:     /*Extract local array */
 73:     VecGetArray(local,&localptr);
 74:     VecGetArray (copy,&copyptr);

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

 82:     VecRestoreArray (copy,&copyptr);
 83:     VecRestoreArray(local,&localptr);

 85:     /* Local to Global */
 86:     DMLocalToGlobalBegin(da,copy,INSERT_VALUES,global);
 87:     DMLocalToGlobalEnd(da,copy,INSERT_VALUES,global);

 89:     /* View Wave */
 90:     /* Set Up Display to Show Heat Graph */
 91: #if defined(PETSC_USE_SOCKET_VIEWER)
 92:     VecView(global,PETSC_VIEWER_SOCKET_WORLD);
 93: #endif
 94:   }

 96:   VecDestroy(&copy);
 97:   VecDestroy(&local);
 98:   VecDestroy(&global);
 99:   DMDestroy(&da);
100:   PetscFinalize();
101:   return 0;
102: }