Actual source code: ex12.c

  2: static char help[] = "Solves the same eigenproblem as in example ex5, but computing also left eigenvectors. "
  3:   "It is a Markov model of a random walk on a triangular grid. "
  4:   "A standard nonsymmetric eigenproblem with real eigenvalues. The rightmost eigenvalue is known to be 1.\n\n"
  5:   "The command line options are:\n"
  6:   "  -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";

 8:  #include slepceps.h

 10: /* 
 11:    User-defined routines
 12: */
 13: PetscErrorCode MatMarkovModel( PetscInt m, Mat A );

 17: int main( int argc, char **argv )
 18: {
 20:   Vec                  v0,temp;          /* initial vector */
 21:   Vec                  *X,*Y;           /* right and left eigenvectors */
 22:   Mat                  A;                  /* operator matrix */
 23:   EPS                  eps;                  /* eigenproblem solver context */
 24:   EPSType              type;
 25:   PetscReal            error1, error2, tol, re, im;
 26:   PetscScalar          kr, ki;
 27:   int                  nev, maxit, i, its, nconv;
 28:   PetscInt             N, m=15;

 30:   SlepcInitialize(&argc,&argv,(char*)0,help);

 32:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 33:   N = m*(m+1)/2;
 34:   PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%d (m=%d)\n\n",N,m);

 36:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 37:      Compute the operator matrix that defines the eigensystem, Ax=kx
 38:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 40:   MatCreate(PETSC_COMM_WORLD,&A);
 41:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 42:   MatSetFromOptions(A);
 43:   MatMarkovModel( m, A );

 45:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 46:                 Create the eigensolver and set various options
 47:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 49:   /* 
 50:      Create eigensolver context
 51:   */
 52:   EPSCreate(PETSC_COMM_WORLD,&eps);

 54:   /* 
 55:      Set operators. In this case, it is a standard eigenvalue problem
 56:   */
 57:   EPSSetOperators(eps,A,PETSC_NULL);
 58:   EPSSetProblemType(eps,EPS_NHEP);

 60:   /*
 61:      Select a two-sided version of the eigensolver so that left eigenvectors
 62:      are also computed
 63:   */
 64:   EPSSetClass(eps,EPS_TWO_SIDE);

 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:   MatGetVecs(A,&v0,&temp);
 76:   VecSet(v0,1.0);
 77:   MatMult(A,v0,temp);
 78:   EPSSetInitialVector(eps,v0);
 79:   EPSSetLeftInitialVector(eps,temp);

 81:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 82:                       Solve the eigensystem
 83:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 85:   EPSSolve(eps);
 86:   EPSGetIterationNumber(eps, &its);
 87:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);

 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,PETSC_NULL);
 95:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
 96:   EPSGetTolerances(eps,&tol,&maxit);
 97:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);

 99:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
100:                     Display solution and clean up
101:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

103:   /* 
104:      Get number of converged approximate eigenpairs
105:   */
106:   EPSGetConverged(eps,&nconv);
107:   PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);

109:   if (nconv>0) {
110:     /*
111:        Display eigenvalues and relative errors
112:     */
113:     PetscPrintf(PETSC_COMM_WORLD,
114:          "           k          ||Ax-kx||/||kx||   ||y'A-ky'||/||ky||\n"
115:          "   ----------------- ------------------ --------------------\n" );

117:     for( i=0; i<nconv; i++ ) {
118:       /* 
119:         Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
120:         ki (imaginary part)
121:       */
122:       EPSGetValue(eps,i,&kr,&ki);
123:       /*
124:          Compute the relative errors associated to both right and left eigenvectors
125:       */
126:       EPSComputeRelativeError(eps,i,&error1);
127:       EPSComputeRelativeErrorLeft(eps,i,&error2);

129: #ifdef PETSC_USE_COMPLEX
130:       re = PetscRealPart(kr);
131:       im = PetscImaginaryPart(kr);
132: #else
133:       re = kr;
134:       im = ki;
135: #endif 
136:       if (im!=0.0) {
137:         PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12g%12g\n",re,im,error1,error2);
138:       } else {
139:         PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g       %12g\n",re,error1,error2);
140:       }
141:     }
142:     PetscPrintf(PETSC_COMM_WORLD,"\n" );

144:     VecDuplicateVecs(v0,nconv,&X);
145:     VecDuplicateVecs(temp,nconv,&Y);
146:     for (i=0;i<nconv;i++) {
147:       EPSGetRightVector(eps,i,X[i],PETSC_NULL);
148:       EPSGetLeftVector(eps,i,Y[i],PETSC_NULL);
149:     }
150:     PetscPrintf(PETSC_COMM_WORLD,
151:          "                   Bi-orthogonality <x,y>                   \n"
152:          "   ---------------------------------------------------------\n" );

154:     SlepcCheckOrthogonality(X,nconv,Y,nconv,PETSC_NULL,PETSC_NULL);
155:     PetscPrintf(PETSC_COMM_WORLD,"\n" );
156:     VecDestroyVecs(X,nconv);
157:     VecDestroyVecs(Y,nconv);

159:   }
160: 
161:   /* 
162:      Free work space
163:   */
164:   VecDestroy(v0);
165:   VecDestroy(temp);
166:   EPSDestroy(eps);
167:   MatDestroy(A);
168:   SlepcFinalize();
169:   return 0;
170: }

174: /*
175:     Matrix generator for a Markov model of a random walk on a triangular grid.

177:     This subroutine generates a test matrix that models a random walk on a 
178:     triangular grid. This test example was used by G. W. Stewart ["{SRRIT} - a 
179:     FORTRAN subroutine to calculate the dominant invariant subspaces of a real
180:     matrix", Tech. report. TR-514, University of Maryland (1978).] and in a few
181:     papers on eigenvalue problems by Y. Saad [see e.g. LAA, vol. 34, pp. 269-295
182:     (1980) ]. These matrices provide reasonably easy test problems for eigenvalue
183:     algorithms. The transpose of the matrix  is stochastic and so it is known 
184:     that one is an exact eigenvalue. One seeks the eigenvector of the transpose 
185:     associated with the eigenvalue unity. The problem is to calculate the steady
186:     state probability distribution of the system, which is the eigevector 
187:     associated with the eigenvalue one and scaled in such a way that the sum all
188:     the components is equal to one.

190:     Note: the code will actually compute the transpose of the stochastic matrix
191:     that contains the transition probabilities.
192: */
193: PetscErrorCode MatMarkovModel( PetscInt m, Mat A )
194: {
195:   const PetscReal cst = 0.5/(PetscReal)(m-1);
196:   PetscReal       pd, pu;
197:   PetscErrorCode  ierr;
198:   PetscInt        i, j, jmax, ix=0, Istart, Iend;

201:   MatGetOwnershipRange(A,&Istart,&Iend);
202:   for( i=1; i<=m; i++ ) {
203:     jmax = m-i+1;
204:     for( j=1; j<=jmax; j++ ) {
205:       ix = ix + 1;
206:       if( ix-1<Istart || ix>Iend ) continue;  /* compute only owned rows */
207:       if( j!=jmax ) {
208:         pd = cst*(PetscReal)(i+j-1);
209:         /* north */
210:         if( i==1 ) {
211:           MatSetValue( A, ix-1, ix, 2*pd, INSERT_VALUES );
212:         }
213:         else {
214:           MatSetValue( A, ix-1, ix, pd, INSERT_VALUES );
215:         }
216:         /* east */
217:         if( j==1 ) {
218:           MatSetValue( A, ix-1, ix+jmax-1, 2*pd, INSERT_VALUES );
219:         }
220:         else {
221:           MatSetValue( A, ix-1, ix+jmax-1, pd, INSERT_VALUES );
222:         }
223:       }
224:       /* south */
225:       pu = 0.5 - cst*(PetscReal)(i+j-3);
226:       if( j>1 ) {
227:         MatSetValue( A, ix-1, ix-2, pu, INSERT_VALUES );
228:       }
229:       /* west */
230:       if( i>1 ) {
231:         MatSetValue( A, ix-1, ix-jmax-2, pu, INSERT_VALUES );
232:       }
233:     }
234:   }
235:   MatAssemblyBegin( A, MAT_FINAL_ASSEMBLY );
236:   MatAssemblyEnd( A, MAT_FINAL_ASSEMBLY );
237:   return(0);
238: }