Actual source code: ex5.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[] = "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 tol;
44: PetscInt N,m=15,nev,maxit,its;
46:
47: SlepcInitialize(&argc,&argv,(char*)0,help);
49: PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
50: N = m*(m+1)/2;
51: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n\n",N,m);
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Compute the operator matrix that defines the eigensystem, Ax=kx
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57: MatCreate(PETSC_COMM_WORLD,&A);
58: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
59: MatSetFromOptions(A);
60: MatMarkovModel(m,A);
62: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
63: Create the eigensolver and set various options
64: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
66: /*
67: Create eigensolver context
68: */
69: EPSCreate(PETSC_COMM_WORLD,&eps);
71: /*
72: Set operators. In this case, it is a standard eigenvalue problem
73: */
74: EPSSetOperators(eps,A,PETSC_NULL);
75: EPSSetProblemType(eps,EPS_NHEP);
77: /*
78: Set solver parameters at runtime
79: */
80: EPSSetFromOptions(eps);
82: /*
83: Set the initial vector. This is optional, if not done the initial
84: vector is set to random values
85: */
86: MatGetVecs(A,&v0,PETSC_NULL);
87: VecSet(v0,1.0);
88: EPSSetInitialSpace(eps,1,&v0);
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Solve the eigensystem
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: EPSSolve(eps);
95: EPSGetIterationNumber(eps,&its);
96: PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);
98: /*
99: Optional: Get some information from the solver and display it
100: */
101: EPSGetType(eps,&type);
102: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
103: EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);
104: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
105: EPSGetTolerances(eps,&tol,&maxit);
106: PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4G, maxit=%D\n",tol,maxit);
108: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109: Display solution and clean up
110: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
112: EPSPrintSolution(eps,PETSC_NULL);
113: EPSDestroy(&eps);
114: MatDestroy(&A);
115: VecDestroy(&v0);
116: SlepcFinalize();
117: return 0;
118: }
122: /*
123: Matrix generator for a Markov model of a random walk on a triangular grid.
125: This subroutine generates a test matrix that models a random walk on a
126: triangular grid. This test example was used by G. W. Stewart ["{SRRIT} - a
127: FORTRAN subroutine to calculate the dominant invariant subspaces of a real
128: matrix", Tech. report. TR-514, University of Maryland (1978).] and in a few
129: papers on eigenvalue problems by Y. Saad [see e.g. LAA, vol. 34, pp. 269-295
130: (1980) ]. These matrices provide reasonably easy test problems for eigenvalue
131: algorithms. The transpose of the matrix is stochastic and so it is known
132: that one is an exact eigenvalue. One seeks the eigenvector of the transpose
133: associated with the eigenvalue unity. The problem is to calculate the steady
134: state probability distribution of the system, which is the eigevector
135: associated with the eigenvalue one and scaled in such a way that the sum all
136: the components is equal to one.
138: Note: the code will actually compute the transpose of the stochastic matrix
139: that contains the transition probabilities.
140: */
141: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
142: {
143: const PetscReal cst = 0.5/(PetscReal)(m-1);
144: PetscReal pd,pu;
145: PetscInt Istart,Iend,i,j,jmax,ix=0;
146: PetscErrorCode ierr;
149: MatGetOwnershipRange(A,&Istart,&Iend);
150: for (i=1;i<=m;i++) {
151: jmax = m-i+1;
152: for (j=1;j<=jmax;j++) {
153: ix = ix + 1;
154: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
155: if (j!=jmax) {
156: pd = cst*(PetscReal)(i+j-1);
157: /* north */
158: if (i==1) {
159: MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
160: } else {
161: MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
162: }
163: /* east */
164: if (j==1) {
165: MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
166: } else {
167: MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
168: }
169: }
170: /* south */
171: pu = 0.5 - cst*(PetscReal)(i+j-3);
172: if (j>1) {
173: MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
174: }
175: /* west */
176: if (i>1) {
177: MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
178: }
179: }
180: }
181: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
182: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
183: return(0);
184: }