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 <petscdmda.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: DM 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: DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,M,w,s,PETSC_NULL,&da);
30: DMCreateGlobalVector(da,&global);
31: DMCreateLocalVector(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,©);
37: VecGetArray (copy,©ptr);
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,©ptr);
55: VecRestoreArray(local,&localptr);
56: DMLocalToGlobalBegin(da,local,INSERT_VALUES,global);
57: DMLocalToGlobalEnd(da,local,INSERT_VALUES,global);
59: /* Assign Parameters */
60: h= 1.0/M;
61: k= h*h/2.2;
63: for (j=0; j<time_steps; j++) {
65: /* Global to Local */
66: DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);
67: DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);
69: /*Extract local array */
70: VecGetArray(local,&localptr);
71: VecGetArray (copy,©ptr);
73: /* Update Locally - Make array of new values */
74: /* Note: I don't do anything for the first and last entry */
75: for (i=1; i< localsize-1; i++) {
76: copyptr[i] = localptr[i] + (k/(h*h)) *
77: (localptr[i+1]-2.0*localptr[i]+localptr[i-1]);
78: }
79:
80: VecRestoreArray (copy,©ptr);
81: VecRestoreArray(local,&localptr);
83: /* Local to Global */
84: DMLocalToGlobalBegin(da,copy,INSERT_VALUES,global);
85: DMLocalToGlobalEnd(da,copy,INSERT_VALUES,global);
86:
87: /* View Wave */
88: /* Set Up Display to Show Heat Graph */
89: #if defined(PETSC_USE_SOCKET_VIEWER)
90: VecView(global,PETSC_VIEWER_SOCKET_WORLD);
91: #endif
92: }
94: VecDestroy(©);
95: VecDestroy(&local);
96: VecDestroy(&global);
97: DMDestroy(&da);
98: PetscFinalize();
99: return 0;
100: }
101: