Actual source code: ex5.c

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

  6:    This file is part of SLEPc.

  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[] = "Eigenvalue problem associated with a Markov model of a random walk on a triangular grid. "
 23:   "It is a standard nonsymmetric eigenproblem with real eigenvalues and the rightmost eigenvalue is known to be 1.\n"
 24:   "This example illustrates how the user can set the initial vector.\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: {
 39:   Vec            v0;              /* initial vector */
 40:   Mat            A;               /* operator matrix */
 41:   EPS            eps;             /* eigenproblem solver context */
 42:   EPSType        type;
 43:   PetscInt       N,m=15,nev;

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

 48:   PetscOptionsGetInt(NULL,"-m",&m,NULL);
 49:   N = m*(m+1)/2;
 50:   PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n\n",N,m);

 52:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 53:      Compute the operator matrix that defines the eigensystem, Ax=kx
 54:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 56:   MatCreate(PETSC_COMM_WORLD,&A);
 57:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 58:   MatSetFromOptions(A);
 59:   MatSetUp(A);
 60:   MatMarkovModel(m,A);

 62:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 63:                 Create the eigensolver and set various options
 64:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

 77:   /*
 78:      Set solver parameters at runtime
 79:   */
 80:   EPSSetFromOptions(eps);

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

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

 94:   EPSSolve(eps);

 96:   /*
 97:      Optional: Get some information from the solver and display it
 98:   */
 99:   EPSGetType(eps,&type);
100:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
101:   EPSGetDimensions(eps,&nev,NULL,NULL);
102:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

104:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105:                     Display solution and clean up
106:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

108:   EPSPrintSolution(eps,NULL);
109:   EPSDestroy(&eps);
110:   MatDestroy(&A);
111:   VecDestroy(&v0);
112:   SlepcFinalize();
113:   return 0;
114: }

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

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

134:     Note: the code will actually compute the transpose of the stochastic matrix
135:     that contains the transition probabilities.
136: */
137: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
138: {
139:   const PetscReal cst = 0.5/(PetscReal)(m-1);
140:   PetscReal       pd,pu;
141:   PetscInt        Istart,Iend,i,j,jmax,ix=0;
142:   PetscErrorCode  ierr;

145:   MatGetOwnershipRange(A,&Istart,&Iend);
146:   for (i=1;i<=m;i++) {
147:     jmax = m-i+1;
148:     for (j=1;j<=jmax;j++) {
149:       ix = ix + 1;
150:       if (ix-1<Istart || ix>Iend) continue;  /* compute only owned rows */
151:       if (j!=jmax) {
152:         pd = cst*(PetscReal)(i+j-1);
153:         /* north */
154:         if (i==1) {
155:           MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
156:         } else {
157:           MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
158:         }
159:         /* east */
160:         if (j==1) {
161:           MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
162:         } else {
163:           MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
164:         }
165:       }
166:       /* south */
167:       pu = 0.5 - cst*(PetscReal)(i+j-3);
168:       if (j>1) {
169:         MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
170:       }
171:       /* west */
172:       if (i>1) {
173:         MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
174:       }
175:     }
176:   }
177:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
178:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
179:   return(0);
180: }