Actual source code: ex12.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2011, Universitat 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[] = "Solves the same eigenproblem as in example ex5, but computing also left eigenvectors. "
23: "It is a Markov model of a random walk on a triangular grid. "
24: "A standard nonsymmetric eigenproblem with real eigenvalues. The rightmost eigenvalue is known to be 1.\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,w0; /* initial vector */
40: Vec *X,*Y; /* right and left eigenvectors */
41: Mat A; /* operator matrix */
42: EPS eps; /* eigenproblem solver context */
43: const EPSType type;
44: PetscReal error1,error2,tol,re,im;
45: PetscScalar kr,ki;
46: PetscInt nev,maxit,i,its,nconv,N,m=15;
49: SlepcInitialize(&argc,&argv,(char*)0,help);
51: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
52: N = m*(m+1)/2;
53: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n\n",N,m);
55: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56: Compute the operator matrix that defines the eigensystem, Ax=kx
57: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59: MatCreate(PETSC_COMM_WORLD,&A);
60: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
61: MatSetFromOptions(A);
62: MatMarkovModel(m,A);
64: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65: Create the eigensolver and set various options
66: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68: /*
69: Create eigensolver context
70: */
71: EPSCreate(PETSC_COMM_WORLD,&eps);
73: /*
74: Set operators. In this case, it is a standard eigenvalue problem.
75: Request also left eigenvectors
76: */
77: EPSSetOperators(eps,A,PETSC_NULL);
78: EPSSetProblemType(eps,EPS_NHEP);
79: EPSSetLeftVectorsWanted(eps,PETSC_TRUE);
81: /*
82: Set solver parameters at runtime
83: */
84: EPSSetFromOptions(eps);
86: /*
87: Set the initial vector. This is optional, if not done the initial
88: vector is set to random values
89: */
90: MatGetVecs(A,&v0,&w0);
91: VecSet(v0,1.0);
92: MatMult(A,v0,w0);
93: EPSSetInitialSpace(eps,1,&v0);
94: EPSSetInitialSpaceLeft(eps,1,&w0);
96: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: Solve the eigensystem
98: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100: EPSSolve(eps);
101: EPSGetIterationNumber(eps,&its);
102: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);
104: /*
105: Optional: Get some information from the solver and display it
106: */
107: EPSGetType(eps,&type);
108: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
109: EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);
110: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
111: EPSGetTolerances(eps,&tol,&maxit);
112: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4G, maxit=%D\n",tol,maxit);
114: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115: Display solution and clean up
116: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118: /*
119: Get number of converged approximate eigenpairs
120: */
121: EPSGetConverged(eps,&nconv);
122: PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %D\n\n",nconv);
124: if (nconv>0) {
125: /*
126: Display eigenvalues and relative errors
127: */
128: PetscPrintf(PETSC_COMM_WORLD,
129: " k ||Ax-kx||/||kx|| ||y'A-ky'||/||ky||\n"
130: " ----------------- ------------------ --------------------\n");
132: for (i=0;i<nconv;i++) {
133: /*
134: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
135: ki (imaginary part)
136: */
137: EPSGetEigenvalue(eps,i,&kr,&ki);
138: /*
139: Compute the relative errors associated to both right and left eigenvectors
140: */
141: EPSComputeRelativeError(eps,i,&error1);
142: EPSComputeRelativeErrorLeft(eps,i,&error2);
144: #if defined(PETSC_USE_COMPLEX)
145: re = PetscRealPart(kr);
146: im = PetscImaginaryPart(kr);
147: #else
148: re = kr;
149: im = ki;
150: #endif
151: if (im!=0.0) {
152: PetscPrintf(PETSC_COMM_WORLD," %9F%+9F j %12G%12G\n",re,im,error1,error2);
153: } else {
154: PetscPrintf(PETSC_COMM_WORLD," %12F %12G %12G\n",re,error1,error2);
155: }
156: }
157: PetscPrintf(PETSC_COMM_WORLD,"\n");
159: VecDuplicateVecs(v0,nconv,&X);
160: VecDuplicateVecs(w0,nconv,&Y);
161: for (i=0;i<nconv;i++) {
162: EPSGetEigenvector(eps,i,X[i],PETSC_NULL);
163: EPSGetEigenvectorLeft(eps,i,Y[i],PETSC_NULL);
164: }
165: PetscPrintf(PETSC_COMM_WORLD,
166: " Bi-orthogonality <x,y> \n"
167: " ---------------------------------------------------------\n");
169: SlepcCheckOrthogonality(X,nconv,Y,nconv,PETSC_NULL,PETSC_NULL);
170: PetscPrintf(PETSC_COMM_WORLD,"\n");
171: VecDestroyVecs(nconv,&X);
172: VecDestroyVecs(nconv,&Y);
174: }
175:
176: /*
177: Free work space
178: */
179: VecDestroy(&v0);
180: VecDestroy(&w0);
181: EPSDestroy(&eps);
182: MatDestroy(&A);
183: SlepcFinalize();
184: return 0;
185: }
189: /*
190: Matrix generator for a Markov model of a random walk on a triangular grid.
192: This subroutine generates a test matrix that models a random walk on a
193: triangular grid. This test example was used by G. W. Stewart ["{SRRIT} - a
194: FORTRAN subroutine to calculate the dominant invariant subspaces of a real
195: matrix", Tech. report. TR-514, University of Maryland (1978).] and in a few
196: papers on eigenvalue problems by Y. Saad [see e.g. LAA, vol. 34, pp. 269-295
197: (1980) ]. These matrices provide reasonably easy test problems for eigenvalue
198: algorithms. The transpose of the matrix is stochastic and so it is known
199: that one is an exact eigenvalue. One seeks the eigenvector of the transpose
200: associated with the eigenvalue unity. The problem is to calculate the steady
201: state probability distribution of the system, which is the eigevector
202: associated with the eigenvalue one and scaled in such a way that the sum all
203: the components is equal to one.
204: Note: the code will actually compute the transpose of the stochastic matrix
205: that contains the transition probabilities.
206: */
207: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
208: {
209: const PetscReal cst = 0.5/(PetscReal)(m-1);
210: PetscReal pd,pu;
211: PetscInt i,j,jmax,ix=0,Istart,Iend;
212: PetscErrorCode ierr;
215: MatGetOwnershipRange(A,&Istart,&Iend);
216: for (i=1;i<=m;i++) {
217: jmax = m-i+1;
218: for (j=1;j<=jmax;j++) {
219: ix = ix + 1;
220: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
221: if (j!=jmax) {
222: pd = cst*(PetscReal)(i+j-1);
223: /* north */
224: if (i==1) {
225: MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
226: } else {
227: MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
228: }
229: /* east */
230: if (j==1) {
231: MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
232: } else {
233: MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
234: }
235: }
236: /* south */
237: pu = 0.5 - cst*(PetscReal)(i+j-3);
238: if (j>1) {
239: MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
240: }
241: /* west */
242: if (i>1) {
243: MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
244: }
245: }
246: }
247: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
248: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
249: return(0);
250: }