Actual source code: ex5.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[] = "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 );

 26: int main( int argc, char **argv )
 27: {
 28:   Vec                  v0;                  /* initial vector */
 29:   Mat                  A;                  /* operator matrix */
 30:   EPS                  eps;                  /* eigenproblem solver context */
 31:   EPSType              type;
 32:   PetscReal            error, tol, re, im;
 33:   PetscScalar          kr, ki;
 34:   PetscInt             N, m=15;
 35:   int                  nev, maxit, i, its, nconv;
 37: 
 38:   SlepcInitialize(&argc,&argv,(char*)0,help);

 40:   PetscOptionsGetInt(PETSC_NULL,"-m",&m,PETSC_NULL);
 41:   N = m*(m+1)/2;
 42:   PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%d (m=%d)\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:   MatMarkovModel( m, A );

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

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

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

 68:   /*
 69:      Set solver parameters at runtime
 70:   */
 71:   EPSSetFromOptions(eps);

 73:   /*
 74:      Set the initial vector. This is optional, if not done the initial
 75:      vector is set to random values
 76:   */
 77:   MatGetVecs(A,&v0,PETSC_NULL);
 78:   VecSet(v0,1.0);
 79:   EPSSetInitialVector(eps,v0);

 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||\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:       EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
123:       /*
124:          Compute the relative error associated to each eigenpair
125:       */
126:       EPSComputeRelativeError(eps,i,&error);

128: #ifdef PETSC_USE_COMPLEX
129:       re = PetscRealPart(kr);
130:       im = PetscImaginaryPart(kr);
131: #else
132:       re = kr;
133:       im = ki;
134: #endif 
135:       if (im!=0.0) {
136:         PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12g\n",re,im,error);
137:       } else {
138:         PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g\n",re,error);
139:       }
140:     }
141:     PetscPrintf(PETSC_COMM_WORLD,"\n" );
142:   }
143: 
144:   /* 
145:      Free work space
146:   */
147:   EPSDestroy(eps);
148:   MatDestroy(A);
149:   VecDestroy(v0);
150:   SlepcFinalize();
151:   return 0;
152: }

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

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

172:     Note: the code will actually compute the transpose of the stochastic matrix
173:     that contains the transition probabilities.
174: */
175: PetscErrorCode MatMarkovModel( PetscInt m, Mat A )
176: {
177:   const PetscReal cst = 0.5/(PetscReal)(m-1);
178:   PetscReal       pd, pu;
179:   PetscErrorCode  ierr;
180:   PetscInt        Istart, Iend, i, j, jmax, ix=0;

183:   MatGetOwnershipRange(A,&Istart,&Iend);
184:   for( i=1; i<=m; i++ ) {
185:     jmax = m-i+1;
186:     for( j=1; j<=jmax; j++ ) {
187:       ix = ix + 1;
188:       if( ix-1<Istart || ix>Iend ) continue;  /* compute only owned rows */
189:       if( j!=jmax ) {
190:         pd = cst*(PetscReal)(i+j-1);
191:         /* north */
192:         if( i==1 ) {
193:           MatSetValue( A, ix-1, ix, 2*pd, INSERT_VALUES );
194:         }
195:         else {
196:           MatSetValue( A, ix-1, ix, pd, INSERT_VALUES );
197:         }
198:         /* east */
199:         if( j==1 ) {
200:           MatSetValue( A, ix-1, ix+jmax-1, 2*pd, INSERT_VALUES );
201:         }
202:         else {
203:           MatSetValue( A, ix-1, ix+jmax-1, pd, INSERT_VALUES );
204:         }
205:       }
206:       /* south */
207:       pu = 0.5 - cst*(PetscReal)(i+j-3);
208:       if( j>1 ) {
209:         MatSetValue( A, ix-1, ix-2, pu, INSERT_VALUES );
210:       }
211:       /* west */
212:       if( i>1 ) {
213:         MatSetValue( A, ix-1, ix-jmax-2, pu, INSERT_VALUES );
214:       }
215:     }
216:   }
217:   MatAssemblyBegin( A, MAT_FINAL_ASSEMBLY );
218:   MatAssemblyEnd( A, MAT_FINAL_ASSEMBLY );
219:   return(0);
220: }