Actual source code: ex12.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:       SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:       Copyright (c) 2002-2007, Universidad Politecnica de Valencia, Spain

  6:       This file is part of SLEPc. See the README file for conditions of use
  7:       and additional information.
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: static char help[] = "Solves the same eigenproblem as in example ex5, but computing also left eigenvectors. "
 12:   "It is a Markov model of a random walk on a triangular grid. "
 13:   "A standard nonsymmetric eigenproblem with real eigenvalues. The rightmost eigenvalue is known to be 1.\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 );

 26: int main( int argc, char **argv )
 27: {
 29:   Vec                  v0,temp;          /* initial vector */
 30:   Vec                  *X,*Y;           /* right and left eigenvectors */
 31:   Mat                  A;                  /* operator matrix */
 32:   EPS                  eps;                  /* eigenproblem solver context */
 33:   EPSType              type;
 34:   PetscReal            error1, error2, tol, re, im;
 35:   PetscScalar          kr, ki;
 36:   int                  nev, maxit, i, its, nconv;
 37:   PetscInt             N, m=15;

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

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

 45:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 46:      Compute the operator matrix that defines the eigensystem, Ax=kx
 47:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 49:   MatCreate(PETSC_COMM_WORLD,&A);
 50:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 51:   MatSetFromOptions(A);
 52:   MatMarkovModel( m, A );

 54:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 55:                 Create the eigensolver and set various options
 56:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 58:   /* 
 59:      Create eigensolver context
 60:   */
 61:   EPSCreate(PETSC_COMM_WORLD,&eps);

 63:   /* 
 64:      Set operators. In this case, it is a standard eigenvalue problem
 65:   */
 66:   EPSSetOperators(eps,A,PETSC_NULL);
 67:   EPSSetProblemType(eps,EPS_NHEP);

 69:   /*
 70:      Select a two-sided version of the eigensolver so that left eigenvectors
 71:      are also computed
 72:   */
 73:   EPSSetClass(eps,EPS_TWO_SIDE);

 75:   /*
 76:      Set solver parameters at runtime
 77:   */
 78:   EPSSetFromOptions(eps);

 80:   /*
 81:      Set the initial vector. This is optional, if not done the initial
 82:      vector is set to random values
 83:   */
 84:   MatGetVecs(A,&v0,&temp);
 85:   VecSet(v0,1.0);
 86:   MatMult(A,v0,temp);
 87:   EPSSetInitialVector(eps,v0);
 88:   EPSSetLeftInitialVector(eps,temp);

 90:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 91:                       Solve the eigensystem
 92:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 94:   EPSSolve(eps);
 95:   EPSGetIterationNumber(eps, &its);
 96:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);

 98:   /*
 99:      Optional: Get some information from the solver and display it
100:   */
101:   EPSGetType(eps,&type);
102:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
103:   EPSGetDimensions(eps,&nev,PETSC_NULL);
104:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
105:   EPSGetTolerances(eps,&tol,&maxit);
106:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);

108:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
109:                     Display solution and clean up
110:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

112:   /* 
113:      Get number of converged approximate eigenpairs
114:   */
115:   EPSGetConverged(eps,&nconv);
116:   PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);

118:   if (nconv>0) {
119:     /*
120:        Display eigenvalues and relative errors
121:     */
122:     PetscPrintf(PETSC_COMM_WORLD,
123:          "           k          ||Ax-kx||/||kx||   ||y'A-ky'||/||ky||\n"
124:          "   ----------------- ------------------ --------------------\n" );

126:     for( i=0; i<nconv; i++ ) {
127:       /* 
128:         Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
129:         ki (imaginary part)
130:       */
131:       EPSGetValue(eps,i,&kr,&ki);
132:       /*
133:          Compute the relative errors associated to both right and left eigenvectors
134:       */
135:       EPSComputeRelativeError(eps,i,&error1);
136:       EPSComputeRelativeErrorLeft(eps,i,&error2);

138: #ifdef PETSC_USE_COMPLEX
139:       re = PetscRealPart(kr);
140:       im = PetscImaginaryPart(kr);
141: #else
142:       re = kr;
143:       im = ki;
144: #endif 
145:       if (im!=0.0) {
146:         PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12g%12g\n",re,im,error1,error2);
147:       } else {
148:         PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g       %12g\n",re,error1,error2);
149:       }
150:     }
151:     PetscPrintf(PETSC_COMM_WORLD,"\n" );

153:     VecDuplicateVecs(v0,nconv,&X);
154:     VecDuplicateVecs(temp,nconv,&Y);
155:     for (i=0;i<nconv;i++) {
156:       EPSGetRightVector(eps,i,X[i],PETSC_NULL);
157:       EPSGetLeftVector(eps,i,Y[i],PETSC_NULL);
158:     }
159:     PetscPrintf(PETSC_COMM_WORLD,
160:          "                   Bi-orthogonality <x,y>                   \n"
161:          "   ---------------------------------------------------------\n" );

163:     SlepcCheckOrthogonality(X,nconv,Y,nconv,PETSC_NULL,PETSC_NULL);
164:     PetscPrintf(PETSC_COMM_WORLD,"\n" );
165:     VecDestroyVecs(X,nconv);
166:     VecDestroyVecs(Y,nconv);

168:   }
169: 
170:   /* 
171:      Free work space
172:   */
173:   VecDestroy(v0);
174:   VecDestroy(temp);
175:   EPSDestroy(eps);
176:   MatDestroy(A);
177:   SlepcFinalize();
178:   return 0;
179: }

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

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

199:     Note: the code will actually compute the transpose of the stochastic matrix
200:     that contains the transition probabilities.
201: */
202: PetscErrorCode MatMarkovModel( PetscInt m, Mat A )
203: {
204:   const PetscReal cst = 0.5/(PetscReal)(m-1);
205:   PetscReal       pd, pu;
206:   PetscErrorCode  ierr;
207:   PetscInt        i, j, jmax, ix=0, Istart, Iend;

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