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 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL)); 24 CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL)); 25 26 /* Set up the array */ 27 CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,M,w,s,NULL,&da)); 28 CHKERRQ(DMSetFromOptions(da)); 29 CHKERRQ(DMSetUp(da)); 30 CHKERRQ(DMCreateGlobalVector(da,&global)); 31 CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 32 CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 33 34 /* Make copy of local array for doing updates */ 35 CHKERRQ(DMCreateLocalVector(da,&local)); 36 37 /* determine starting point of each processor */ 38 CHKERRQ(VecGetOwnershipRange(global,&mybase,&myend)); 39 40 /* Initialize the Array */ 41 CHKERRQ(VecGetLocalSize (global,&globalsize)); 42 CHKERRQ(VecGetArray (global,&globalptr)); 43 44 for (i=0; i<globalsize; i++) { 45 j = i + mybase; 46 47 globalptr[i] = PetscSinReal((PETSC_PI*j*6)/((PetscReal)M) + 1.2 * PetscSinReal((PETSC_PI*j*2)/((PetscReal)M))) * 4+4; 48 } 49 50 CHKERRQ(VecRestoreArray(global,&localptr)); 51 52 /* Assign Parameters */ 53 h= 1.0/M; 54 k= h*h/2.2; 55 CHKERRQ(VecGetLocalSize(local,&localsize)); 56 57 for (j=0; j<time_steps; j++) { 58 59 /* Global to Local */ 60 CHKERRQ(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local)); 61 CHKERRQ(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local)); 62 63 /*Extract local array */ 64 CHKERRQ(VecGetArray(local,&localptr)); 65 CHKERRQ(VecGetArray (global,&globalptr)); 66 67 /* Update Locally - Make array of new values */ 68 /* Note: I don't do anything for the first and last entry */ 69 for (i=1; i< localsize-1; i++) { 70 globalptr[i-1] = localptr[i] + (k/(h*h)) * (localptr[i+1]-2.0*localptr[i]+localptr[i-1]); 71 } 72 73 CHKERRQ(VecRestoreArray (global,&globalptr)); 74 CHKERRQ(VecRestoreArray(local,&localptr)); 75 76 /* View Wave */ 77 /* Set Up Display to Show Heat Graph */ 78 #if defined(PETSC_USE_SOCKET_VIEWER) 79 CHKERRQ(VecView(global,PETSC_VIEWER_SOCKET_WORLD)); 80 #endif 81 } 82 83 CHKERRQ(VecDestroy(&local)); 84 CHKERRQ(VecDestroy(&global)); 85 CHKERRQ(DMDestroy(&da)); 86 ierr = PetscFinalize(); 87 return ierr; 88 } 89