Actual source code: ex12.c

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

  6:    This file is part of SLEPc.
  7:       
  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY 
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS 
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for 
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

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

 28:  #include slepceps.h

 30: /* 
 31:    User-defined routines
 32: */
 33: PetscErrorCode MatMarkovModel( PetscInt m, Mat A );

 37: int main( int argc, char **argv )
 38: {
 40:   Vec                  v0,temp;          /* initial vector */
 41:   Vec                  *X,*Y;           /* right and left eigenvectors */
 42:   Mat                  A;                  /* operator matrix */
 43:   EPS                  eps;                  /* eigenproblem solver context */
 44:   const EPSType  type;
 45:   PetscReal            error1, error2, tol, re, im;
 46:   PetscScalar          kr, ki;
 47:   PetscInt             nev, maxit, i, its, nconv, N, m=15;

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

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

 55:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 56:      Compute the operator matrix that defines the eigensystem, Ax=kx
 57:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 59:   MatCreate(PETSC_COMM_WORLD,&A);
 60:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 61:   MatSetFromOptions(A);
 62:   MatMarkovModel( m, A );

 64:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
 65:                 Create the eigensolver and set various options
 66:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 68:   /* 
 69:      Create eigensolver context
 70:   */
 71:   EPSCreate(PETSC_COMM_WORLD,&eps);

 73:   /* 
 74:      Set operators. In this case, it is a standard eigenvalue problem
 75:   */
 76:   EPSSetOperators(eps,A,PETSC_NULL);
 77:   EPSSetProblemType(eps,EPS_NHEP);

 79:   /*
 80:      Select a two-sided version of the eigensolver so that left eigenvectors
 81:      are also computed
 82:   */
 83:   EPSSetClass(eps,EPS_TWO_SIDE);

 85:   /*
 86:      Set solver parameters at runtime
 87:   */
 88:   EPSSetFromOptions(eps);

 90:   /*
 91:      Set the initial vector. This is optional, if not done the initial
 92:      vector is set to random values
 93:   */
 94:   MatGetVecs(A,&v0,&temp);
 95:   VecSet(v0,1.0);
 96:   MatMult(A,v0,temp);
 97:   EPSSetInitialVector(eps,v0);
 98:   EPSSetLeftInitialVector(eps,temp);

100:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
101:                       Solve the eigensystem
102:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

104:   EPSSolve(eps);
105:   EPSGetIterationNumber(eps, &its);
106:   PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %d\n",its);

108:   /*
109:      Optional: Get some information from the solver and display it
110:   */
111:   EPSGetType(eps,&type);
112:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
113:   EPSGetDimensions(eps,&nev,PETSC_NULL,PETSC_NULL);
114:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %d\n",nev);
115:   EPSGetTolerances(eps,&tol,&maxit);
116:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4g, maxit=%d\n",tol,maxit);

118:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
119:                     Display solution and clean up
120:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

122:   /* 
123:      Get number of converged approximate eigenpairs
124:   */
125:   EPSGetConverged(eps,&nconv);
126:   PetscPrintf(PETSC_COMM_WORLD," Number of converged approximate eigenpairs: %d\n\n",nconv);

128:   if (nconv>0) {
129:     /*
130:        Display eigenvalues and relative errors
131:     */
132:     PetscPrintf(PETSC_COMM_WORLD,
133:          "           k          ||Ax-kx||/||kx||   ||y'A-ky'||/||ky||\n"
134:          "   ----------------- ------------------ --------------------\n" );

136:     for( i=0; i<nconv; i++ ) {
137:       /* 
138:         Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
139:         ki (imaginary part)
140:       */
141:       EPSGetValue(eps,i,&kr,&ki);
142:       /*
143:          Compute the relative errors associated to both right and left eigenvectors
144:       */
145:       EPSComputeRelativeError(eps,i,&error1);
146:       EPSComputeRelativeErrorLeft(eps,i,&error2);

148: #ifdef PETSC_USE_COMPLEX
149:       re = PetscRealPart(kr);
150:       im = PetscImaginaryPart(kr);
151: #else
152:       re = kr;
153:       im = ki;
154: #endif 
155:       if (im!=0.0) {
156:         PetscPrintf(PETSC_COMM_WORLD," %9f%+9f j %12g%12g\n",re,im,error1,error2);
157:       } else {
158:         PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g       %12g\n",re,error1,error2);
159:       }
160:     }
161:     PetscPrintf(PETSC_COMM_WORLD,"\n" );

163:     VecDuplicateVecs(v0,nconv,&X);
164:     VecDuplicateVecs(temp,nconv,&Y);
165:     for (i=0;i<nconv;i++) {
166:       EPSGetRightVector(eps,i,X[i],PETSC_NULL);
167:       EPSGetLeftVector(eps,i,Y[i],PETSC_NULL);
168:     }
169:     PetscPrintf(PETSC_COMM_WORLD,
170:          "                   Bi-orthogonality <x,y>                   \n"
171:          "   ---------------------------------------------------------\n" );

173:     SlepcCheckOrthogonality(X,nconv,Y,nconv,PETSC_NULL,PETSC_NULL);
174:     PetscPrintf(PETSC_COMM_WORLD,"\n" );
175:     VecDestroyVecs(X,nconv);
176:     VecDestroyVecs(Y,nconv);

178:   }
179: 
180:   /* 
181:      Free work space
182:   */
183:   VecDestroy(v0);
184:   VecDestroy(temp);
185:   EPSDestroy(eps);
186:   MatDestroy(A);
187:   SlepcFinalize();
188:   return 0;
189: }

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

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

209:     Note: the code will actually compute the transpose of the stochastic matrix
210:     that contains the transition probabilities.
211: */
212: PetscErrorCode MatMarkovModel( PetscInt m, Mat A )
213: {
214:   const PetscReal cst = 0.5/(PetscReal)(m-1);
215:   PetscReal       pd, pu;
216:   PetscErrorCode  ierr;
217:   PetscInt        i, j, jmax, ix=0, Istart, Iend;

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