Actual source code: ex41.c
slepc-3.11.2 2019-07-30
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[] = "Illustrates the computation of left eigenvectors.\n\n"
12: "The problem is the Markov model as in ex5.c.\n"
13: "The command line options are:\n"
14: " -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";
16: #include <slepceps.h>
18: /*
19: User-defined routines
20: */
21: PetscErrorCode MatMarkovModel(PetscInt,Mat);
22: PetscErrorCode ComputeResidualNorm(Mat,PetscBool,PetscScalar,PetscScalar,Vec,Vec,Vec,PetscReal*);
24: int main(int argc,char **argv)
25: {
26: Mat A; /* operator matrix */
27: EPS eps; /* eigenproblem solver context */
28: EPSType type;
29: PetscInt i,N,m=15,nconv;
30: PetscBool twosided;
31: PetscReal nrmr,nrml=0.0,re,im,lev;
32: PetscScalar *kr,*ki;
33: Vec t,*xr,*xi,*yr,*yi;
36: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
38: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
39: N = m*(m+1)/2;
40: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n\n",N,m);
42: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43: Compute the operator matrix that defines the eigensystem, Ax=kx
44: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
46: MatCreate(PETSC_COMM_WORLD,&A);
47: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
48: MatSetFromOptions(A);
49: MatSetUp(A);
50: MatMarkovModel(m,A);
51: MatCreateVecs(A,NULL,&t);
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Create the eigensolver and set various options
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57: EPSCreate(PETSC_COMM_WORLD,&eps);
58: EPSSetOperators(eps,A,NULL);
59: EPSSetProblemType(eps,EPS_NHEP);
61: /* use a two-sided algorithm to compute left eigenvectors as well */
62: EPSSetTwoSided(eps,PETSC_TRUE);
64: /* allow user to change settings at run time */
65: EPSSetFromOptions(eps);
66: EPSGetTwoSided(eps,&twosided);
68: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69: Solve the eigensystem
70: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
72: EPSSolve(eps);
74: /*
75: Optional: Get some information from the solver and display it
76: */
77: EPSGetType(eps,&type);
78: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
80: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81: Display solution and clean up
82: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
84: /*
85: Get number of converged approximate eigenpairs
86: */
87: EPSGetConverged(eps,&nconv);
88: PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %D\n\n",nconv);
89: PetscMalloc2(nconv,&kr,nconv,&ki);
90: VecDuplicateVecs(t,nconv,&xr);
91: VecDuplicateVecs(t,nconv,&xi);
92: if (twosided) {
93: VecDuplicateVecs(t,nconv,&yr);
94: VecDuplicateVecs(t,nconv,&yi);
95: }
97: if (nconv>0) {
98: /*
99: Display eigenvalues and relative errors
100: */
101: PetscPrintf(PETSC_COMM_WORLD,
102: " k ||Ax-kx|| ||y'A-y'k||\n"
103: " ---------------- ------------------ ------------------\n");
105: for (i=0;i<nconv;i++) {
106: /*
107: Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
108: ki (imaginary part)
109: */
110: EPSGetEigenpair(eps,i,&kr[i],&ki[i],xr[i],xi[i]);
111: if (twosided) {
112: EPSGetLeftEigenvector(eps,i,yr[i],yi[i]);
113: }
114: /*
115: Compute the residual norms associated to each eigenpair
116: */
117: ComputeResidualNorm(A,PETSC_FALSE,kr[i],ki[i],xr[i],xi[i],t,&nrmr);
118: if (twosided) {
119: ComputeResidualNorm(A,PETSC_TRUE,kr[i],ki[i],yr[i],yi[i],t,&nrml);
120: }
122: #if defined(PETSC_USE_COMPLEX)
123: re = PetscRealPart(kr[i]);
124: im = PetscImaginaryPart(kr[i]);
125: #else
126: re = kr[i];
127: im = ki[i];
128: #endif
129: if (im!=0.0) {
130: PetscPrintf(PETSC_COMM_WORLD," %8f%+8fi %12g %12g\n",(double)re,(double)im,(double)nrmr,(double)nrml);
131: } else {
132: PetscPrintf(PETSC_COMM_WORLD," %12f %12g %12g\n",(double)re,(double)nrmr,(double)nrml);
133: }
134: }
135: PetscPrintf(PETSC_COMM_WORLD,"\n");
136: /*
137: Check bi-orthogonality of eigenvectors
138: */
139: if (twosided) {
140: VecCheckOrthogonality(xr,nconv,yr,nconv,NULL,NULL,&lev);
141: if (lev<100*PETSC_MACHINE_EPSILON) {
142: PetscPrintf(PETSC_COMM_WORLD," Level of bi-orthogonality of eigenvectors < 100*eps\n\n");
143: } else {
144: PetscPrintf(PETSC_COMM_WORLD," Level of bi-orthogonality of eigenvectors: %g\n\n",(double)lev);
145: }
146: }
147: }
149: EPSDestroy(&eps);
150: MatDestroy(&A);
151: VecDestroy(&t);
152: PetscFree2(kr,ki);
153: VecDestroyVecs(nconv,&xr);
154: VecDestroyVecs(nconv,&xi);
155: if (twosided) {
156: VecDestroyVecs(nconv,&yr);
157: VecDestroyVecs(nconv,&yi);
158: }
159: SlepcFinalize();
160: return ierr;
161: }
163: /*
164: Matrix generator for a Markov model of a random walk on a triangular grid.
166: This subroutine generates a test matrix that models a random walk on a
167: triangular grid. This test example was used by G. W. Stewart ["{SRRIT} - a
168: FORTRAN subroutine to calculate the dominant invariant subspaces of a real
169: matrix", Tech. report. TR-514, University of Maryland (1978).] and in a few
170: papers on eigenvalue problems by Y. Saad [see e.g. LAA, vol. 34, pp. 269-295
171: (1980) ]. These matrices provide reasonably easy test problems for eigenvalue
172: algorithms. The transpose of the matrix is stochastic and so it is known
173: that one is an exact eigenvalue. One seeks the eigenvector of the transpose
174: associated with the eigenvalue unity. The problem is to calculate the steady
175: state probability distribution of the system, which is the eigevector
176: associated with the eigenvalue one and scaled in such a way that the sum all
177: the components is equal to one.
179: Note: the code will actually compute the transpose of the stochastic matrix
180: that contains the transition probabilities.
181: */
182: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
183: {
184: const PetscReal cst = 0.5/(PetscReal)(m-1);
185: PetscReal pd,pu;
186: PetscInt Istart,Iend,i,j,jmax,ix=0;
187: PetscErrorCode ierr;
190: MatGetOwnershipRange(A,&Istart,&Iend);
191: for (i=1;i<=m;i++) {
192: jmax = m-i+1;
193: for (j=1;j<=jmax;j++) {
194: ix = ix + 1;
195: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
196: if (j!=jmax) {
197: pd = cst*(PetscReal)(i+j-1);
198: /* north */
199: if (i==1) {
200: MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
201: } else {
202: MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
203: }
204: /* east */
205: if (j==1) {
206: MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
207: } else {
208: MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
209: }
210: }
211: /* south */
212: pu = 0.5 - cst*(PetscReal)(i+j-3);
213: if (j>1) {
214: MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
215: }
216: /* west */
217: if (i>1) {
218: MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
219: }
220: }
221: }
222: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
223: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
224: return(0);
225: }
227: /*
228: ComputeResidualNorm - Computes the norm of the residual vector
229: associated with an eigenpair.
231: Input Parameters:
232: trans - whether A' must be used instead of A
233: kr,ki - eigenvalue
234: xr,xi - eigenvector
235: u - work vector
236: */
237: PetscErrorCode ComputeResidualNorm(Mat A,PetscBool trans,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,Vec u,PetscReal *norm)
238: {
240: #if !defined(PETSC_USE_COMPLEX)
241: PetscReal ni,nr;
242: #endif
243: PetscErrorCode (*matmult)(Mat,Vec,Vec) = trans? MatMultTranspose: MatMult;
246: #if !defined(PETSC_USE_COMPLEX)
247: if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
248: #endif
249: (*matmult)(A,xr,u);
250: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
251: VecAXPY(u,-kr,xr);
252: }
253: VecNorm(u,NORM_2,norm);
254: #if !defined(PETSC_USE_COMPLEX)
255: } else {
256: (*matmult)(A,xr,u);
257: if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
258: VecAXPY(u,-kr,xr);
259: VecAXPY(u,ki,xi);
260: }
261: VecNorm(u,NORM_2,&nr);
262: (*matmult)(A,xi,u);
263: if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
264: VecAXPY(u,-kr,xi);
265: VecAXPY(u,-ki,xr);
266: }
267: VecNorm(u,NORM_2,&ni);
268: *norm = SlepcAbsEigenvalue(nr,ni);
269: }
270: #endif
271: return(0);
272: }
274: /*TEST
276: testset:
277: args: -st_type sinvert -eps_target 1.1 -eps_nev 4
278: filter: grep -v method | sed -e "s/[+-]0.0*i//" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
279: requires: !single
280: output_file: output/ex41_1.out
281: test:
282: suffix: 1
283: args: -eps_type {{power krylovschur}}
284: test:
285: suffix: 1_balance
286: args: -eps_balance {{oneside twoside}} -eps_ncv 16
288: TEST*/