Actual source code: ex5.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7:
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
22: static char help[] = "Eigenvalue problem associated with a Markov model of a random walk on a triangular grid. "
23: "It is a standard nonsymmetric eigenproblem with real eigenvalues and the rightmost eigenvalue is known to be 1.\n"
24: "This example illustrates how the user can set the initial vector.\n\n"
25: "The command line options are:\n"
26: " -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";
28: #include slepceps.h
30: /*
31: User-defined routines
32: */
33: PetscErrorCode MatMarkovModel( PetscInt m, Mat A );
37: int main( int argc, char **argv )
38: {
39: Vec v0; /* initial vector */
40: Mat A; /* operator matrix */
41: EPS eps; /* eigenproblem solver context */
42: const EPSType type;
43: PetscReal error, tol, re, im;
44: PetscScalar kr, ki;
45: PetscInt N, m=15, nev, maxit, i, its, nconv;
47:
48: SlepcInitialize(&argc,&argv,(char*)0,help);
50: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
51: N = m*(m+1)/2;
52: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%d (m=%d)\n\n",N,m);
54: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55: Compute the operator matrix that defines the eigensystem, Ax=kx
56: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58: MatCreate(PETSC_COMM_WORLD,&A);
59: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
60: MatSetFromOptions(A);
61: MatMarkovModel( m, A );
63: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64: Create the eigensolver and set various options
65: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
67: /*
68: Create eigensolver context
69: */
70: EPSCreate(PETSC_COMM_WORLD,&eps);
72: /*
73: Set operators. In this case, it is a standard eigenvalue problem
74: */
75: EPSSetOperators(eps,A,PETSC_NULL);
76: EPSSetProblemType(eps,EPS_NHEP);
78: /*
79: Set solver parameters at runtime
80: */
81: EPSSetFromOptions(eps);
83: /*
84: Set the initial vector. This is optional, if not done the initial
85: vector is set to random values
86: */
87: MatGetVecs(A,&v0,PETSC_NULL);
88: VecSet(v0,1.0);
89: EPSSetInitialVector(eps,v0);
91: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92: Solve the eigensystem
93: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
95: EPSSolve(eps);
96: EPSGetIterationNumber(eps, &its);
97: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);
99: /*
100: Optional: Get some information from the solver and display it
101: */
102: EPSGetType(eps,&type);
103: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
104: EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);
105: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
106: EPSGetTolerances(eps,&tol,&maxit);
107: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);
109: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110: Display solution and clean up
111: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113: /*
114: Get number of converged approximate eigenpairs
115: */
116: EPSGetConverged(eps,&nconv);
117: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
119: if (nconv>0) {
120: /*
121: Display eigenvalues and relative errors
122: */
123: PetscPrintf(PETSC_COMM_WORLD,
124: " k ||Ax-kx||/||kx||\n"
125: " ----------------- ------------------\n" );
127: for( i=0; i<nconv; i++ ) {
128: /*
129: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
130: ki (imaginary part)
131: */
132: EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
133: /*
134: Compute the relative error associated to each eigenpair
135: */
136: EPSComputeRelativeError(eps,i,&error);
138: #ifdef PETSC_USE_COMPLEX
139: re = PetscRealPart(kr);
140: im = PetscImaginaryPart(kr);
141: #else
142: re = kr;
143: im = ki;
144: #endif
145: if (im!=0.0) {
146: PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12g\n",re,im,error);
147: } else {
148: PetscPrintf(PETSC_COMM_WORLD," %12f %12g\n",re,error);
149: }
150: }
151: PetscPrintf(PETSC_COMM_WORLD,"\n" );
152: }
153:
154: /*
155: Free work space
156: */
157: EPSDestroy(eps);
158: MatDestroy(A);
159: VecDestroy(v0);
160: SlepcFinalize();
161: return 0;
162: }
166: /*
167: Matrix generator for a Markov model of a random walk on a triangular grid.
169: This subroutine generates a test matrix that models a random walk on a
170: triangular grid. This test example was used by G. W. Stewart ["{SRRIT} - a
171: FORTRAN subroutine to calculate the dominant invariant subspaces of a real
172: matrix", Tech. report. TR-514, University of Maryland (1978).] and in a few
173: papers on eigenvalue problems by Y. Saad [see e.g. LAA, vol. 34, pp. 269-295
174: (1980) ]. These matrices provide reasonably easy test problems for eigenvalue
175: algorithms. The transpose of the matrix is stochastic and so it is known
176: that one is an exact eigenvalue. One seeks the eigenvector of the transpose
177: associated with the eigenvalue unity. The problem is to calculate the steady
178: state probability distribution of the system, which is the eigevector
179: associated with the eigenvalue one and scaled in such a way that the sum all
180: the components is equal to one.
182: Note: the code will actually compute the transpose of the stochastic matrix
183: that contains the transition probabilities.
184: */
185: PetscErrorCode MatMarkovModel( PetscInt m, Mat A )
186: {
187: const PetscReal cst = 0.5/(PetscReal)(m-1);
188: PetscReal pd, pu;
189: PetscErrorCode ierr;
190: PetscInt Istart, Iend, i, j, jmax, ix=0;
193: MatGetOwnershipRange(A,&Istart,&Iend);
194: for( i=1; i<=m; i++ ) {
195: jmax = m-i+1;
196: for( j=1; j<=jmax; j++ ) {
197: ix = ix + 1;
198: if( ix-1<Istart || ix>Iend ) continue; /* compute only owned rows */
199: if( j!=jmax ) {
200: pd = cst*(PetscReal)(i+j-1);
201: /* north */
202: if( i==1 ) {
203: MatSetValue( A, ix-1, ix, 2*pd, INSERT_VALUES );
204: }
205: else {
206: MatSetValue( A, ix-1, ix, pd, INSERT_VALUES );
207: }
208: /* east */
209: if( j==1 ) {
210: MatSetValue( A, ix-1, ix+jmax-1, 2*pd, INSERT_VALUES );
211: }
212: else {
213: MatSetValue( A, ix-1, ix+jmax-1, pd, INSERT_VALUES );
214: }
215: }
216: /* south */
217: pu = 0.5 - cst*(PetscReal)(i+j-3);
218: if( j>1 ) {
219: MatSetValue( A, ix-1, ix-2, pu, INSERT_VALUES );
220: }
221: /* west */
222: if( i>1 ) {
223: MatSetValue( A, ix-1, ix-jmax-2, pu, INSERT_VALUES );
224: }
225: }
226: }
227: MatAssemblyBegin( A, MAT_FINAL_ASSEMBLY );
228: MatAssemblyEnd( A, MAT_FINAL_ASSEMBLY );
229: return(0);
230: }