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