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 petscda.h
10: #include petscsys.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: DA da;
20: Vec local,global,copy;
21: PetscScalar *localptr,*copyptr;
22: PetscReal h,k;
23:
24: PetscInitialize(&argc,&argv,(char*)0,help);
26: PetscOptionsGetInt(PETSC_NULL,"-M",&M,PETSC_NULL);
27: PetscOptionsGetInt(PETSC_NULL,"-time",&time_steps,PETSC_NULL);
28:
29: /* Set up the array */
30: DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,M,w,s,PETSC_NULL,&da);
31: DACreateGlobalVector(da,&global);
32: DACreateLocalVector(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,©);
38: VecGetArray (copy,©ptr);
41: /* determine starting point of each processor */
42: VecGetOwnershipRange(global,&mybase,&myend);
44: /* Initialize the Array */
45: VecGetLocalSize (local,&localsize);
46: VecGetArray (local,&localptr);
47: localptr[0] = copyptr[0] = 0.0;
48: localptr[localsize-1] = copyptr[localsize-1] = 1.0;
49: for (i=1; i<localsize-1; i++) {
50: j=(i-1)+mybase;
51: localptr[i] = sin((PETSC_PI*j*6)/((PetscReal)M)
52: + 1.2 * sin((PETSC_PI*j*2)/((PetscReal)M))) * 4+4;
53: }
55: VecRestoreArray (copy,©ptr);
56: VecRestoreArray(local,&localptr);
57: DALocalToGlobal(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: DAGlobalToLocalBegin(da,global,INSERT_VALUES,local);
67: DAGlobalToLocalEnd(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: DALocalToGlobal(da,copy,INSERT_VALUES,global);
85:
86: /* View Wave */
87: /* Set Up Display to Show Heat Graph */
88: #if defined(PETSC_USE_SOCKET_VIEWER)
89: VecView(global,PETSC_VIEWER_SOCKET_WORLD);
90: #endif
91: }
93: VecDestroy(copy);
94: VecDestroy(local);
95: VecDestroy(global);
96: DADestroy(da);
97: PetscFinalize();
98: return 0;
99: }
100: