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