1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
11: static char help[] = "Solves the same problem as in ex5, but with a user-defined sorting criterion."
12: "It is a standard nonsymmetric eigenproblem with real eigenvalues and the rightmost eigenvalue is known to be 1.\n"
13: "This example illustrates how the user can set a custom spectrum selection.\n\n"
14: "The command line options are:\n"
15: " -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";
17: #include <slepceps.h>
19: /*
20: User-defined routines
21: */
23: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx);
24: PetscErrorCode MatMarkovModel(PetscInt m,Mat A);
26: int main(int argc,char **argv) 27: {
28: Mat A; /* operator matrix */
29: EPS eps; /* eigenproblem solver context */
30: EPSType type;
31: PetscScalar target=0.5;
32: PetscInt N,m=15,nev;
33: PetscBool terse;
34: PetscViewer viewer;
36: char str[50];
38: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
40: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
41: N = m*(m+1)/2;
42: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n",N,m);
43: PetscOptionsGetScalar(NULL,NULL,"-target",&target,NULL);
44: SlepcSNPrintfScalar(str,50,target,PETSC_FALSE);
45: PetscPrintf(PETSC_COMM_WORLD,"Searching closest eigenvalues to the right of %s.\n\n",str);
47: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48: Compute the operator matrix that defines the eigensystem, Ax=kx
49: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51: MatCreate(PETSC_COMM_WORLD,&A);
52: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
53: MatSetFromOptions(A);
54: MatSetUp(A);
55: MatMarkovModel(m,A);
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Create the eigensolver and set various options
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: /*
62: Create eigensolver context
63: */
64: EPSCreate(PETSC_COMM_WORLD,&eps);
66: /*
67: Set operators. In this case, it is a standard eigenvalue problem
68: */
69: EPSSetOperators(eps,A,NULL);
70: EPSSetProblemType(eps,EPS_NHEP);
72: /*
73: Set the custom comparing routine in order to obtain the eigenvalues
74: closest to the target on the right only
75: */
76: EPSSetEigenvalueComparison(eps,MyEigenSort,&target);
78: /*
79: Set solver parameters at runtime
80: */
81: EPSSetFromOptions(eps);
83: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84: Solve the eigensystem
85: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87: EPSSolve(eps);
89: /*
90: Optional: Get some information from the solver and display it
91: */
92: EPSGetType(eps,&type);
93: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
94: EPSGetDimensions(eps,&nev,NULL,NULL);
95: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
97: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98: Display solution and clean up
99: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
101: /* show detailed info unless -terse option is given by user */
102: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
103: if (terse) {
104: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
105: } else {
106: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
107: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
108: EPSReasonView(eps,viewer);
109: EPSErrorView(eps,EPS_ERROR_RELATIVE,viewer);
110: PetscViewerPopFormat(viewer);
111: }
112: EPSDestroy(&eps);
113: MatDestroy(&A);
114: SlepcFinalize();
115: return ierr;
116: }
118: /*
119: Matrix generator for a Markov model of a random walk on a triangular grid.
121: This subroutine generates a test matrix that models a random walk on a
122: triangular grid. This test example was used by G. W. Stewart ["{SRRIT} - a
123: FORTRAN subroutine to calculate the dominant invariant subspaces of a real
124: matrix", Tech. report. TR-514, University of Maryland (1978).] and in a few
125: papers on eigenvalue problems by Y. Saad [see e.g. LAA, vol. 34, pp. 269-295
126: (1980) ]. These matrices provide reasonably easy test problems for eigenvalue
127: algorithms. The transpose of the matrix is stochastic and so it is known
128: that one is an exact eigenvalue. One seeks the eigenvector of the transpose
129: associated with the eigenvalue unity. The problem is to calculate the steady
130: state probability distribution of the system, which is the eigevector
131: associated with the eigenvalue one and scaled in such a way that the sum all
132: the components is equal to one.
134: Note: the code will actually compute the transpose of the stochastic matrix
135: that contains the transition probabilities.
136: */
137: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)138: {
139: const PetscReal cst = 0.5/(PetscReal)(m-1);
140: PetscReal pd,pu;
141: PetscInt Istart,Iend,i,j,jmax,ix=0;
142: PetscErrorCode ierr;
145: MatGetOwnershipRange(A,&Istart,&Iend);
146: for (i=1;i<=m;i++) {
147: jmax = m-i+1;
148: for (j=1;j<=jmax;j++) {
149: ix = ix + 1;
150: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
151: if (j!=jmax) {
152: pd = cst*(PetscReal)(i+j-1);
153: /* north */
154: if (i==1) {
155: MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
156: } else {
157: MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
158: }
159: /* east */
160: if (j==1) {
161: MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
162: } else {
163: MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
164: }
165: }
166: /* south */
167: pu = 0.5 - cst*(PetscReal)(i+j-3);
168: if (j>1) {
169: MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
170: }
171: /* west */
172: if (i>1) {
173: MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
174: }
175: }
176: }
177: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
178: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
179: return(0);
180: }
182: /*
183: Function for user-defined eigenvalue ordering criterion.
185: Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
186: one of them as the preferred one according to the criterion.
187: In this example, the preferred value is the one closest to the target,
188: but on the right side.
189: */
190: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)191: {
192: PetscScalar target = *(PetscScalar*)ctx;
193: PetscReal da,db;
194: PetscBool aisright,bisright;
197: if (PetscRealPart(target) < PetscRealPart(ar)) aisright = PETSC_TRUE;
198: else aisright = PETSC_FALSE;
199: if (PetscRealPart(target) < PetscRealPart(br)) bisright = PETSC_TRUE;
200: else bisright = PETSC_FALSE;
201: if (aisright == bisright) {
202: /* both are on the same side of the target */
203: da = SlepcAbsEigenvalue(ar-target,ai);
204: db = SlepcAbsEigenvalue(br-target,bi);
205: if (da < db) *r = -1;
206: else if (da > db) *r = 1;
207: else *r = 0;
208: } else if (aisright && !bisright) *r = -1; /* 'a' is on the right */
209: else *r = 1; /* 'b' is on the right */
210: return(0);
211: }
213: /*TEST
215: test:
216: suffix: 1
217: args: -eps_nev 4 -terse
218: requires: !single
220: TEST*/