Actual source code: ex5.c
2: static char help[] = "Eigenvalue problem associated with a Markov model of a random walk on a triangular grid. "
3: "It is a standard symmetric eigenproblem and the rightmost eigenvalue is known to be 1.\n\n"
4: "This example illustrates how the user can set the initial vector.\n\n"
5: "The command line options are:\n\n"
6: " -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";
8: #include slepceps.h
10: /*
11: User-defined routines
12: */
13: extern int MatMarkovModel( int m, Mat A );
17: int main( int argc, char **argv )
18: {
19: Vec v0; /* initial vector */
20: Mat A; /* operator matrix */
21: EPS eps; /* eigenproblem solver context */
22: EPSType type;
23: PetscReal error, tol, re, im;
24: PetscScalar kr, ki;
25: int N, M, m=15, nev, ierr, maxit, i, its, nconv;
26: PetscScalar alpha;
28: SlepcInitialize(&argc,&argv,(char*)0,help);
30: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
31: N = m*(m+1)/2;
32: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%d (m=%d)\n\n",N,m);
34: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35: Compute the operator matrix that defines the eigensystem, Ax=kx
36: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
38: MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&A);
39: MatSetFromOptions(A);
40: MatMarkovModel( m, A );
42: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43: Create the eigensolver and set various options
44: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
46: /*
47: Create eigensolver context
48: */
49: EPSCreate(PETSC_COMM_WORLD,&eps);
51: /*
52: Set operators. In this case, it is a standard eigenvalue problem
53: */
54: EPSSetOperators(eps,A,PETSC_NULL);
56: /*
57: Set solver parameters at runtime
58: */
59: EPSSetFromOptions(eps);
61: /*
62: Set the initial vector. This is optional, if not done the initial
63: vector is set to random values
64: */
65: MatGetLocalSize(A,PETSC_NULL,&M);
66: VecCreate(PETSC_COMM_WORLD,&v0);
67: VecSetSizes(v0,M,PETSC_DECIDE);
68: VecSetFromOptions(v0);
69: alpha = 1.0;
70: VecSet(&alpha,v0);
71: EPSSetInitialVector(eps,v0);
73: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
74: Solve the eigensystem
75: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77: EPSSolve(eps);
78: EPSGetIterationNumber(eps, &its);
79: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
81: /*
82: Optional: Get some information from the solver and display it
83: */
84: EPSGetType(eps,&type);
85: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
86: EPSGetDimensions(eps,&nev,PETSC_NULL);
87: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
88: EPSGetTolerances(eps,&tol,&maxit);
89: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
91: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92: Display solution and clean up
93: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95: /*
96: Get number of converged approximate eigenpairs
97: */
98: EPSGetConverged(eps,&nconv);
99: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
100:
102: if (nconv>0) {
103: /*
104: Display eigenvalues and relative errors
105: */
106: PetscPrintf(PETSC_COMM_WORLD,
107: " k ||Ax-kx||/||kx||\n"
108: " ----------------- ------------------\n" );
110: for( i=0; i<nconv; i++ ) {
111: /*
112: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
113: ki (imaginary part)
114: */
115: EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
116: /*
117: Compute the relative error associated to each eigenpair
118: */
119: EPSComputeRelativeError(eps,i,&error);
121: #ifdef PETSC_USE_COMPLEX
122: re = PetscRealPart(kr);
123: im = PetscImaginaryPart(kr);
124: #else
125: re = kr;
126: im = ki;
127: #endif
128: if (im!=0.0) {
129: PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12f\n",re,im,error);
130: } else {
131: PetscPrintf(PETSC_COMM_WORLD," %12f %12f\n",re,error);
132: }
133: }
134: PetscPrintf(PETSC_COMM_WORLD,"\n" );
135: }
136:
137: /*
138: Free work space
139: */
140: EPSDestroy(eps);
141: MatDestroy(A);
142: SlepcFinalize();
143: return 0;
144: }
148: /*
149: Matrix generator for a Markov model of a random walk on a triangular grid.
151: This subroutine generates a test matrix that models a random walk on a
152: triangular grid. This test example was used by G. W. Stewart ["{SRRIT} - a
153: FORTRAN subroutine to calculate the dominant invariant subspaces of a real
154: matrix", Tech. report. TR-514, University of Maryland (1978).] and in a few
155: papers on eigenvalue problems by Y. Saad [see e.g. LAA, vol. 34, pp. 269-295
156: (1980) ]. These matrices provide reasonably easy test problems for eigenvalue
157: algorithms. The transpose of the matrix is stochastic and so it is known
158: that one is an exact eigenvalue. One seeks the eigenvector of the transpose
159: associated with the eigenvalue unity. The problem is to calculate the steady
160: state probability distribution of the system, which is the eigevector
161: associated with the eigenvalue one and scaled in such a way that the sum all
162: the components is equal to one.
164: Note: the code will actually compute the transpose of the stochastic matrix
165: that contains the transition probabilities.
166: */
167: int MatMarkovModel( int m, Mat A )
168: {
169: const PetscReal cst = 0.5/(PetscReal)(m-1);
170: PetscReal pd, pu;
171: int ierr, i, j, jmax, ix=0, Istart, Iend;
173: MatGetOwnershipRange(A,&Istart,&Iend);
174: for( i=1; i<=m; i++ ) {
175: jmax = m-i+1;
176: for( j=1; j<=jmax; j++ ) {
177: ix = ix + 1;
178: if( ix-1<Istart || ix>Iend ) continue; /* compute only owned rows */
179: if( j!=jmax ) {
180: pd = cst*(PetscReal)(i+j-1);
181: /* north */
182: if( i==1 ) {
183: MatSetValue( A, ix-1, ix, 2*pd, INSERT_VALUES );
184:
185: }
186: else {
187: MatSetValue( A, ix-1, ix, pd, INSERT_VALUES );
188:
189: }
190: /* east */
191: if( j==1 ) {
192: MatSetValue( A, ix-1, ix+jmax-1, 2*pd, INSERT_VALUES );
193:
194: }
195: else {
196: MatSetValue( A, ix-1, ix+jmax-1, pd, INSERT_VALUES );
197:
198: }
199: }
200: /* south */
201: pu = 0.5 - cst*(PetscReal)(i+j-3);
202: if( j>1 ) {
203: MatSetValue( A, ix-1, ix-2, pu, INSERT_VALUES );
204:
205: }
206: /* west */
207: if( i>1 ) {
208: MatSetValue( A, ix-1, ix-jmax-2, pu, INSERT_VALUES );
209:
210: }
211: }
212: }
213: MatAssemblyBegin( A, MAT_FINAL_ASSEMBLY );
214: MatAssemblyEnd( A, MAT_FINAL_ASSEMBLY );
215: return 0;
216: }