1 2 /* 3 Simple example to show how PETSc programs can be run from MATLAB. 4 See ex12.m. 5 */ 6 7 static char help[] = "Solves the one dimensional heat equation.\n\n"; 8 9 #include <petscdm.h> 10 #include <petscdmda.h> 11 12 int main(int argc,char **argv) 13 { 14 PetscMPIInt rank,size; 15 PetscInt M = 14,time_steps = 20,w=1,s=1,localsize,j,i,mybase,myend,globalsize; 16 PetscErrorCode ierr; 17 DM da; 18 Vec global,local; 19 PetscScalar *globalptr,*localptr; 20 PetscReal h,k; 21 22 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 23 ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr); 24 ierr = PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL);CHKERRQ(ierr); 25 26 /* Set up the array */ 27 ierr = DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,M,w,s,NULL,&da);CHKERRQ(ierr); 28 ierr = DMSetFromOptions(da);CHKERRQ(ierr); 29 ierr = DMSetUp(da);CHKERRQ(ierr); 30 ierr = DMCreateGlobalVector(da,&global);CHKERRQ(ierr); 31 ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 32 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 33 34 /* Make copy of local array for doing updates */ 35 ierr = DMCreateLocalVector(da,&local);CHKERRQ(ierr); 36 37 38 /* determine starting point of each processor */ 39 ierr = VecGetOwnershipRange(global,&mybase,&myend);CHKERRQ(ierr); 40 41 /* Initialize the Array */ 42 ierr = VecGetLocalSize (global,&globalsize);CHKERRQ(ierr); 43 ierr = VecGetArray (global,&globalptr);CHKERRQ(ierr); 44 45 46 for (i=0; i<globalsize; i++) { 47 j = i + mybase; 48 49 globalptr[i] = PetscSinReal((PETSC_PI*j*6)/((PetscReal)M) + 1.2 * PetscSinReal((PETSC_PI*j*2)/((PetscReal)M))) * 4+4; 50 } 51 52 ierr = VecRestoreArray(global,&localptr);CHKERRQ(ierr); 53 54 /* Assign Parameters */ 55 h= 1.0/M; 56 k= h*h/2.2; 57 ierr = VecGetLocalSize(local,&localsize);CHKERRQ(ierr); 58 59 for (j=0; j<time_steps; j++) { 60 61 /* Global to Local */ 62 ierr = DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);CHKERRQ(ierr); 63 ierr = DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);CHKERRQ(ierr); 64 65 /*Extract local array */ 66 ierr = VecGetArray(local,&localptr);CHKERRQ(ierr); 67 ierr = VecGetArray (global,&globalptr);CHKERRQ(ierr); 68 69 /* Update Locally - Make array of new values */ 70 /* Note: I don't do anything for the first and last entry */ 71 for (i=1; i< localsize-1; i++) { 72 globalptr[i-1] = localptr[i] + (k/(h*h)) * (localptr[i+1]-2.0*localptr[i]+localptr[i-1]); 73 } 74 75 ierr = VecRestoreArray (global,&globalptr);CHKERRQ(ierr); 76 ierr = VecRestoreArray(local,&localptr);CHKERRQ(ierr); 77 78 /* View Wave */ 79 /* Set Up Display to Show Heat Graph */ 80 #if defined(PETSC_USE_SOCKET_VIEWER) 81 ierr = VecView(global,PETSC_VIEWER_SOCKET_WORLD);CHKERRQ(ierr); 82 #endif 83 } 84 85 ierr = VecDestroy(&local);CHKERRQ(ierr); 86 ierr = VecDestroy(&global);CHKERRQ(ierr); 87 ierr = DMDestroy(&da);CHKERRQ(ierr); 88 ierr = PetscFinalize(); 89 return ierr; 90 } 91 92