Actual source code: ex18.c
slepc-3.18.3 2023-03-24
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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;
35: char str[50];
38: SlepcInitialize(&argc,&argv,(char*)0,help);
40: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
41: N = m*(m+1)/2;
42: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n",N,m);
43: PetscOptionsGetScalar(NULL,NULL,"-target",&target,NULL);
44: SlepcSNPrintfScalar(str,sizeof(str),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: %" PetscInt_FMT "\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) EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
104: else {
105: PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
106: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
107: EPSConvergedReasonView(eps,viewer);
108: EPSErrorView(eps,EPS_ERROR_RELATIVE,viewer);
109: PetscViewerPopFormat(viewer);
110: }
111: EPSDestroy(&eps);
112: MatDestroy(&A);
113: SlepcFinalize();
114: return 0;
115: }
117: /*
118: Matrix generator for a Markov model of a random walk on a triangular grid.
120: This subroutine generates a test matrix that models a random walk on a
121: triangular grid. This test example was used by G. W. Stewart ["{SRRIT} - a
122: FORTRAN subroutine to calculate the dominant invariant subspaces of a real
123: matrix", Tech. report. TR-514, University of Maryland (1978).] and in a few
124: papers on eigenvalue problems by Y. Saad [see e.g. LAA, vol. 34, pp. 269-295
125: (1980) ]. These matrices provide reasonably easy test problems for eigenvalue
126: algorithms. The transpose of the matrix is stochastic and so it is known
127: that one is an exact eigenvalue. One seeks the eigenvector of the transpose
128: associated with the eigenvalue unity. The problem is to calculate the steady
129: state probability distribution of the system, which is the eigevector
130: associated with the eigenvalue one and scaled in such a way that the sum all
131: the components is equal to one.
133: Note: the code will actually compute the transpose of the stochastic matrix
134: that contains the transition probabilities.
135: */
136: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
137: {
138: const PetscReal cst = 0.5/(PetscReal)(m-1);
139: PetscReal pd,pu;
140: PetscInt Istart,Iend,i,j,jmax,ix=0;
143: MatGetOwnershipRange(A,&Istart,&Iend);
144: for (i=1;i<=m;i++) {
145: jmax = m-i+1;
146: for (j=1;j<=jmax;j++) {
147: ix = ix + 1;
148: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
149: if (j!=jmax) {
150: pd = cst*(PetscReal)(i+j-1);
151: /* north */
152: if (i==1) MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
153: else MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
154: /* east */
155: if (j==1) MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
156: else MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
157: }
158: /* south */
159: pu = 0.5 - cst*(PetscReal)(i+j-3);
160: if (j>1) MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
161: /* west */
162: if (i>1) MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
163: }
164: }
165: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
166: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
167: return 0;
168: }
170: /*
171: Function for user-defined eigenvalue ordering criterion.
173: Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
174: one of them as the preferred one according to the criterion.
175: In this example, the preferred value is the one closest to the target,
176: but on the right side.
177: */
178: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
179: {
180: PetscScalar target = *(PetscScalar*)ctx;
181: PetscReal da,db;
182: PetscBool aisright,bisright;
185: if (PetscRealPart(target) < PetscRealPart(ar)) aisright = PETSC_TRUE;
186: else aisright = PETSC_FALSE;
187: if (PetscRealPart(target) < PetscRealPart(br)) bisright = PETSC_TRUE;
188: else bisright = PETSC_FALSE;
189: if (aisright == bisright) {
190: /* both are on the same side of the target */
191: da = SlepcAbsEigenvalue(ar-target,ai);
192: db = SlepcAbsEigenvalue(br-target,bi);
193: if (da < db) *r = -1;
194: else if (da > db) *r = 1;
195: else *r = 0;
196: } else if (aisright && !bisright) *r = -1; /* 'a' is on the right */
197: else *r = 1; /* 'b' is on the right */
198: return 0;
199: }
201: /*TEST
203: test:
204: suffix: 1
205: args: -eps_nev 4 -terse
206: requires: !single
207: filter: sed -e "s/[+-]0\.0*i//g"
209: TEST*/