1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
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);
34: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx);
38: int main(int argc,char **argv) 39: {
40: Vec v0; /* initial vector */
41: Mat A; /* operator matrix */
42: EPS eps; /* eigenproblem solver context */
43: EPSType type;
44: PetscReal tol=1000*PETSC_MACHINE_EPSILON;
45: PetscInt N,m=15,nev;
46: PetscScalar origin=0.0;
49: SlepcInitialize(&argc,&argv,(char*)0,help);
51: PetscOptionsGetInt(NULL,"-m",&m,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: MatSetUp(A);
63: MatMarkovModel(m,A);
65: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66: Create the eigensolver and set various options
67: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69: /*
70: Create eigensolver context
71: */
72: EPSCreate(PETSC_COMM_WORLD,&eps);
74: /*
75: Set operators. In this case, it is a standard eigenvalue problem
76: */
77: EPSSetOperators(eps,A,NULL);
78: EPSSetProblemType(eps,EPS_NHEP);
79: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
81: /*
82: Set the custom comparing routine in order to obtain the eigenvalues
83: closest to the target on the right only
84: */
85: EPSSetEigenvalueComparison(eps,MyEigenSort,&origin);
88: /*
89: Set solver parameters at runtime
90: */
91: EPSSetFromOptions(eps);
93: /*
94: Set the initial vector. This is optional, if not done the initial
95: vector is set to random values
96: */
97: MatCreateVecs(A,&v0,NULL);
98: VecSetValue(v0,0,-1.5,INSERT_VALUES);
99: VecSetValue(v0,1,2.1,INSERT_VALUES);
100: VecAssemblyBegin(v0);
101: VecAssemblyEnd(v0);
102: EPSSetInitialSpace(eps,1,&v0);
104: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: Solve the eigensystem
106: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108: EPSSolve(eps);
110: /*
111: Optional: Get some information from the solver and display it
112: */
113: EPSGetType(eps,&type);
114: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
115: EPSGetDimensions(eps,&nev,NULL,NULL);
116: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
118: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119: Display solution and clean up
120: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
123: EPSDestroy(&eps);
124: MatDestroy(&A);
125: VecDestroy(&v0);
126: SlepcFinalize();
127: return 0;
128: }
132: /*
133: Matrix generator for a Markov model of a random walk on a triangular grid.
135: This subroutine generates a test matrix that models a random walk on a
136: triangular grid. This test example was used by G. W. Stewart ["{SRRIT} - a
137: FORTRAN subroutine to calculate the dominant invariant subspaces of a real
138: matrix", Tech. report. TR-514, University of Maryland (1978).] and in a few
139: papers on eigenvalue problems by Y. Saad [see e.g. LAA, vol. 34, pp. 269-295
140: (1980) ]. These matrices provide reasonably easy test problems for eigenvalue
141: algorithms. The transpose of the matrix is stochastic and so it is known
142: that one is an exact eigenvalue. One seeks the eigenvector of the transpose
143: associated with the eigenvalue unity. The problem is to calculate the steady
144: state probability distribution of the system, which is the eigevector
145: associated with the eigenvalue one and scaled in such a way that the sum all
146: the components is equal to one.
148: Note: the code will actually compute the transpose of the stochastic matrix
149: that contains the transition probabilities.
150: */
151: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)152: {
153: const PetscReal cst = 0.5/(PetscReal)(m-1);
154: PetscReal pd,pu;
155: PetscInt Istart,Iend,i,j,jmax,ix=0;
156: PetscErrorCode ierr;
159: MatGetOwnershipRange(A,&Istart,&Iend);
160: for (i=1;i<=m;i++) {
161: jmax = m-i+1;
162: for (j=1;j<=jmax;j++) {
163: ix = ix + 1;
164: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
165: if (j!=jmax) {
166: pd = cst*(PetscReal)(i+j-1);
167: /* north */
168: if (i==1) {
169: MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
170: } else {
171: MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
172: }
173: /* east */
174: if (j==1) {
175: MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
176: } else {
177: MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
178: }
179: }
180: /* south */
181: pu = 0.5 - cst*(PetscReal)(i+j-3);
182: if (j>1) {
183: MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
184: }
185: /* west */
186: if (i>1) {
187: MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
188: }
189: }
190: }
191: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
192: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
193: return(0);
194: }
198: /*
199: Function for user-defined eigenvalue ordering criterion.
201: Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
202: one of them as the preferred one according to the criterion.
203: In this example, the preferred value is the one furthest to the origin.
204: */
205: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)206: {
207: PetscScalar origin = *(PetscScalar*)ctx;
208: PetscReal d;
211: d = (SlepcAbsEigenvalue(br-origin,bi) - SlepcAbsEigenvalue(ar-origin,ai))/PetscMax(SlepcAbsEigenvalue(ar-origin,ai),SlepcAbsEigenvalue(br-origin,bi));
212: *r = d > PETSC_SQRT_MACHINE_EPSILON ? 1 : (d < -PETSC_SQRT_MACHINE_EPSILON ? -1 : PetscSign(PetscRealPart(br)));
213: return(0);
214: }