Actual source code: ex5.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2015, 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;
 44:   PetscBool      terse;

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

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

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

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

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

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

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

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

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

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

 95:   EPSSolve(eps);

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

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

109:   /* show detailed info unless -terse option is given by user */
110:   PetscOptionsHasName(NULL,"-terse",&terse);
111:   if (terse) {
112:     EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
113:   } else {
114:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
115:     EPSReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
116:     EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
117:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
118:   }
119:   EPSDestroy(&eps);
120:   MatDestroy(&A);
121:   VecDestroy(&v0);
122:   SlepcFinalize();
123:   return 0;
124: }

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

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

144:     Note: the code will actually compute the transpose of the stochastic matrix
145:     that contains the transition probabilities.
146: */
147: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
148: {
149:   const PetscReal cst = 0.5/(PetscReal)(m-1);
150:   PetscReal       pd,pu;
151:   PetscInt        Istart,Iend,i,j,jmax,ix=0;
152:   PetscErrorCode  ierr;

155:   MatGetOwnershipRange(A,&Istart,&Iend);
156:   for (i=1;i<=m;i++) {
157:     jmax = m-i+1;
158:     for (j=1;j<=jmax;j++) {
159:       ix = ix + 1;
160:       if (ix-1<Istart || ix>Iend) continue;  /* compute only owned rows */
161:       if (j!=jmax) {
162:         pd = cst*(PetscReal)(i+j-1);
163:         /* north */
164:         if (i==1) {
165:           MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
166:         } else {
167:           MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
168:         }
169:         /* east */
170:         if (j==1) {
171:           MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
172:         } else {
173:           MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
174:         }
175:       }
176:       /* south */
177:       pu = 0.5 - cst*(PetscReal)(i+j-3);
178:       if (j>1) {
179:         MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
180:       }
181:       /* west */
182:       if (i>1) {
183:         MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
184:       }
185:     }
186:   }
187:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
188:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
189:   return(0);
190: }