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[] = "Solves the same problem as in ex5, but with a user-defined sorting criterion."
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 a custom spectrum selection.\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: */
34: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx);
35: PetscErrorCode MatMarkovModel(PetscInt m,Mat A);
39: int main(int argc,char **argv) 40: {
41: Vec v0; /* initial vector */
42: Mat A; /* operator matrix */
43: EPS eps; /* eigenproblem solver context */
44: ST st; /* spectral transformation associated */
45: EPSType type;
46: PetscReal tol=1000*PETSC_MACHINE_EPSILON;
47: PetscScalar target=0.5;
48: PetscInt N,m=15,nev;
50: char str[50];
52: SlepcInitialize(&argc,&argv,(char*)0,help);
54: PetscOptionsGetInt(NULL,"-m",&m,NULL);
55: N = m*(m+1)/2;
56: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n",N,m);
57: PetscOptionsGetScalar(NULL,"-target",&target,NULL);
58: SlepcSNPrintfScalar(str,50,target,PETSC_FALSE);
59: PetscPrintf(PETSC_COMM_WORLD,"Searching closest eigenvalues to the right of %s.\n\n",str);
61: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: Compute the operator matrix that defines the eigensystem, Ax=kx
63: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65: MatCreate(PETSC_COMM_WORLD,&A);
66: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
67: MatSetFromOptions(A);
68: MatSetUp(A);
69: MatMarkovModel(m,A);
71: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72: Create the eigensolver and set various options
73: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
75: /*
76: Create eigensolver context
77: */
78: EPSCreate(PETSC_COMM_WORLD,&eps);
80: /*
81: Set operators. In this case, it is a standard eigenvalue problem
82: */
83: EPSSetOperators(eps,A,NULL);
84: EPSSetProblemType(eps,EPS_NHEP);
85: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
87: /*
88: Set the custom comparing routine in order to obtain the eigenvalues
89: closest to the target on the right only
90: */
91: EPSSetEigenvalueComparison(eps,MyEigenSort,&target);
93: /*
94: Set solver parameters at runtime
95: */
96: EPSSetFromOptions(eps);
98: /*
99: Set the preconditioner based on A - target * I
100: */
101: EPSGetST(eps,&st);
102: STSetShift(st,target);
104: /*
105: Set the initial vector. This is optional, if not done the initial
106: vector is set to random values
107: */
108: MatCreateVecs(A,&v0,NULL);
109: VecSet(v0,1.0);
110: EPSSetInitialSpace(eps,1,&v0);
112: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: Solve the eigensystem
114: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116: EPSSolve(eps);
118: /*
119: Optional: Get some information from the solver and display it
120: */
121: EPSGetType(eps,&type);
122: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
123: EPSGetDimensions(eps,&nev,NULL,NULL);
124: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
126: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127: Display solution and clean up
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
131: EPSDestroy(&eps);
132: MatDestroy(&A);
133: VecDestroy(&v0);
134: SlepcFinalize();
135: return 0;
136: }
140: /*
141: Matrix generator for a Markov model of a random walk on a triangular grid.
143: This subroutine generates a test matrix that models a random walk on a
144: triangular grid. This test example was used by G. W. Stewart ["{SRRIT} - a
145: FORTRAN subroutine to calculate the dominant invariant subspaces of a real
146: matrix", Tech. report. TR-514, University of Maryland (1978).] and in a few
147: papers on eigenvalue problems by Y. Saad [see e.g. LAA, vol. 34, pp. 269-295
148: (1980) ]. These matrices provide reasonably easy test problems for eigenvalue
149: algorithms. The transpose of the matrix is stochastic and so it is known
150: that one is an exact eigenvalue. One seeks the eigenvector of the transpose
151: associated with the eigenvalue unity. The problem is to calculate the steady
152: state probability distribution of the system, which is the eigevector
153: associated with the eigenvalue one and scaled in such a way that the sum all
154: the components is equal to one.
156: Note: the code will actually compute the transpose of the stochastic matrix
157: that contains the transition probabilities.
158: */
159: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)160: {
161: const PetscReal cst = 0.5/(PetscReal)(m-1);
162: PetscReal pd,pu;
163: PetscInt Istart,Iend,i,j,jmax,ix=0;
164: PetscErrorCode ierr;
167: MatGetOwnershipRange(A,&Istart,&Iend);
168: for (i=1;i<=m;i++) {
169: jmax = m-i+1;
170: for (j=1;j<=jmax;j++) {
171: ix = ix + 1;
172: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
173: if (j!=jmax) {
174: pd = cst*(PetscReal)(i+j-1);
175: /* north */
176: if (i==1) {
177: MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
178: } else {
179: MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
180: }
181: /* east */
182: if (j==1) {
183: MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
184: } else {
185: MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
186: }
187: }
188: /* south */
189: pu = 0.5 - cst*(PetscReal)(i+j-3);
190: if (j>1) {
191: MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
192: }
193: /* west */
194: if (i>1) {
195: MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
196: }
197: }
198: }
199: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
200: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
201: return(0);
202: }
206: /*
207: Function for user-defined eigenvalue ordering criterion.
209: Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
210: one of them as the preferred one according to the criterion.
211: In this example, the preferred value is the one closest to the target,
212: but on the right side.
213: */
214: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)215: {
216: PetscScalar target = *(PetscScalar*)ctx;
217: PetscReal da,db;
218: PetscBool aisright,bisright;
221: if (PetscRealPart(target) < PetscRealPart(ar)) aisright = PETSC_TRUE;
222: else aisright = PETSC_FALSE;
223: if (PetscRealPart(target) < PetscRealPart(br)) bisright = PETSC_TRUE;
224: else bisright = PETSC_FALSE;
225: if (aisright == bisright) {
226: /* both are on the same side of the target */
227: da = SlepcAbsEigenvalue(ar-target,ai);
228: db = SlepcAbsEigenvalue(br-target,bi);
229: if (da < db) *r = -1;
230: else if (da > db) *r = 1;
231: else *r = 0;
232: } else if (aisright && !bisright) *r = -1; /* 'a' is on the right */
233: else *r = 1; /* 'b' is on the right */
234: return(0);
235: }