Actual source code: ex5.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[] = "Eigenvalue problem associated with a Markov model of a random walk on a triangular grid. "
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 the initial vector.\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: */
22: PetscErrorCode MatMarkovModel(PetscInt m,Mat A);
24: int main(int argc,char **argv)
25: {
26: Vec v0; /* initial vector */
27: Mat A; /* operator matrix */
28: EPS eps; /* eigenproblem solver context */
29: EPSType type;
30: PetscInt N,m=15,nev;
31: PetscMPIInt rank;
32: PetscBool terse;
35: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
37: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
38: N = m*(m+1)/2;
39: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n\n",N,m);
41: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
42: Compute the operator matrix that defines the eigensystem, Ax=kx
43: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
45: MatCreate(PETSC_COMM_WORLD,&A);
46: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
47: MatSetFromOptions(A);
48: MatSetUp(A);
49: MatMarkovModel(m,A);
51: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52: Create the eigensolver and set various options
53: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55: /*
56: Create eigensolver context
57: */
58: EPSCreate(PETSC_COMM_WORLD,&eps);
60: /*
61: Set operators. In this case, it is a standard eigenvalue problem
62: */
63: EPSSetOperators(eps,A,NULL);
64: EPSSetProblemType(eps,EPS_NHEP);
66: /*
67: Set solver parameters at runtime
68: */
69: EPSSetFromOptions(eps);
71: /*
72: Set the initial vector. This is optional, if not done the initial
73: vector is set to random values
74: */
75: MatCreateVecs(A,&v0,NULL);
76: MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
77: if (!rank) {
78: VecSetValue(v0,0,1.0,INSERT_VALUES);
79: VecSetValue(v0,1,1.0,INSERT_VALUES);
80: VecSetValue(v0,2,1.0,INSERT_VALUES);
81: }
82: VecAssemblyBegin(v0);
83: VecAssemblyEnd(v0);
84: EPSSetInitialSpace(eps,1,&v0);
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Solve the eigensystem
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90: EPSSolve(eps);
92: /*
93: Optional: Get some information from the solver and display it
94: */
95: EPSGetType(eps,&type);
96: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
97: EPSGetDimensions(eps,&nev,NULL,NULL);
98: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
100: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101: Display solution and clean up
102: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104: /* show detailed info unless -terse option is given by user */
105: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
106: if (terse) {
107: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
108: } else {
109: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
110: EPSReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
111: EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
112: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
113: }
114: EPSDestroy(&eps);
115: MatDestroy(&A);
116: VecDestroy(&v0);
117: SlepcFinalize();
118: return ierr;
119: }
121: /*
122: Matrix generator for a Markov model of a random walk on a triangular grid.
124: This subroutine generates a test matrix that models a random walk on a
125: triangular grid. This test example was used by G. W. Stewart ["{SRRIT} - a
126: FORTRAN subroutine to calculate the dominant invariant subspaces of a real
127: matrix", Tech. report. TR-514, University of Maryland (1978).] and in a few
128: papers on eigenvalue problems by Y. Saad [see e.g. LAA, vol. 34, pp. 269-295
129: (1980) ]. These matrices provide reasonably easy test problems for eigenvalue
130: algorithms. The transpose of the matrix is stochastic and so it is known
131: that one is an exact eigenvalue. One seeks the eigenvector of the transpose
132: associated with the eigenvalue unity. The problem is to calculate the steady
133: state probability distribution of the system, which is the eigevector
134: associated with the eigenvalue one and scaled in such a way that the sum all
135: the components is equal to one.
137: Note: the code will actually compute the transpose of the stochastic matrix
138: that contains the transition probabilities.
139: */
140: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
141: {
142: const PetscReal cst = 0.5/(PetscReal)(m-1);
143: PetscReal pd,pu;
144: PetscInt Istart,Iend,i,j,jmax,ix=0;
145: PetscErrorCode ierr;
148: MatGetOwnershipRange(A,&Istart,&Iend);
149: for (i=1;i<=m;i++) {
150: jmax = m-i+1;
151: for (j=1;j<=jmax;j++) {
152: ix = ix + 1;
153: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
154: if (j!=jmax) {
155: pd = cst*(PetscReal)(i+j-1);
156: /* north */
157: if (i==1) {
158: MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
159: } else {
160: MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
161: }
162: /* east */
163: if (j==1) {
164: MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
165: } else {
166: MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
167: }
168: }
169: /* south */
170: pu = 0.5 - cst*(PetscReal)(i+j-3);
171: if (j>1) {
172: MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
173: }
174: /* west */
175: if (i>1) {
176: MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
177: }
178: }
179: }
180: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
181: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
182: return(0);
183: }
185: /*TEST
187: test:
188: suffix: 1
189: args: -eps_largest_real -eps_nev 4 -terse
191: TEST*/