Actual source code: ex5.c

  2: static char help[] = "Eigenvalue problem associated with a Markov model of a random walk on a triangular grid. "
  3:   "It is a standard symmetric eigenproblem and the rightmost eigenvalue is known to be 1.\n\n"
  4:   "This example illustrates how the user can set the initial vector.\n\n"
  5:   "The command line options are:\n\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: extern int MatMarkovModel( int m, Mat A );

 17: int main( int argc, char **argv )
 18: {
 19:   Vec         v0;              /* initial vector */
 20:   Mat         A;               /* operator matrix */
 21:   EPS         eps;             /* eigenproblem solver context */
 22:   EPSType     type;
 23:   PetscReal   error, tol, re, im;
 24:   PetscScalar kr, ki;
 25:   int         N, m=15, nev, ierr, maxit, i, its, nconv;
 26:   PetscScalar alpha;

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

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

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

 38:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,N,N,&A);
 39:   MatSetFromOptions(A);
 40:   MatMarkovModel( m, A );

 42:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 43:                 Create the eigensolver and set various options
 44:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 46:   /* 
 47:      Create eigensolver context
 48:   */
 49:   EPSCreate(PETSC_COMM_WORLD,&eps);

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

 57:   /*
 58:      Set solver parameters at runtime
 59:   */
 60:   EPSSetFromOptions(eps);

 62:   /*
 63:      Set the initial vector. This is optional, if not done the initial
 64:      vector is set to random values
 65:   */
 66:   MatGetVecs(A,&v0,PETSC_NULL);
 67:   alpha = 1.0;
 68:   VecSet(&alpha,v0);
 69:   EPSSetInitialVector(eps,v0);

 71:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 72:                       Solve the eigensystem
 73:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 75:   EPSSolve(eps);
 76:   EPSGetIterationNumber(eps, &its);
 77:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);

 79:   /*
 80:      Optional: Get some information from the solver and display it
 81:   */
 82:   EPSGetType(eps,&type);
 83:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
 84:   EPSGetDimensions(eps,&nev,PETSC_NULL);
 85:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
 86:   EPSGetTolerances(eps,&tol,&maxit);
 87:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);

 89:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 90:                     Display solution and clean up
 91:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 93:   /* 
 94:      Get number of converged approximate eigenpairs
 95:   */
 96:   EPSGetConverged(eps,&nconv);
 97:   PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);
 98: 

100:   if (nconv>0) {
101:     /*
102:        Display eigenvalues and relative errors
103:     */
104:     PetscPrintf(PETSC_COMM_WORLD,
105:          "           k          ||Ax-kx||/||kx||\n"
106:          "   ----------------- ------------------\n" );

108:     for( i=0; i<nconv; i++ ) {
109:       /* 
110:         Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
111:         ki (imaginary part)
112:       */
113:       EPSGetEigenpair(eps,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
114:       /*
115:          Compute the relative error associated to each eigenpair
116:       */
117:       EPSComputeRelativeError(eps,i,&error);

119: #ifdef PETSC_USE_COMPLEX
120:       re = PetscRealPart(kr);
121:       im = PetscImaginaryPart(kr);
122: #else
123:       re = kr;
124:       im = ki;
125: #endif 
126:       if (im!=0.0) {
127:         PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12f\n",re,im,error);
128:       } else {
129:         PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12f\n",re,error);
130:       }
131:     }
132:     PetscPrintf(PETSC_COMM_WORLD,"\n" );
133:   }
134: 
135:   /* 
136:      Free work space
137:   */
138:   EPSDestroy(eps);
139:   MatDestroy(A);
140:   VecDestroy(v0);
141:   SlepcFinalize();
142:   return 0;
143: }

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

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

163:     Note: the code will actually compute the transpose of the stochastic matrix
164:     that contains the transition probabilities.
165: */
166: int MatMarkovModel( int m, Mat A )
167: {
168:   const PetscReal cst = 0.5/(PetscReal)(m-1);
169:   PetscReal pd, pu;
170:   int ierr, i, j, jmax, ix=0, Istart, Iend;

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