Actual source code: ex18.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2014, 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[] = "Solves the same problem as in ex5, but with a user-defined sorting criterion."
 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 a custom spectrum selection.\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: */

 34: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx);
 35: PetscErrorCode MatMarkovModel(PetscInt m,Mat A);

 39: int main(int argc,char **argv)
 40: {
 41:   Vec            v0;              /* initial vector */
 42:   Mat            A;               /* operator matrix */
 43:   EPS            eps;             /* eigenproblem solver context */
 44:   EPSType        type;
 45:   PetscScalar    target=0.5;
 46:   PetscInt       N,m=15,nev;
 48:   char           str[50];

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

 52:   PetscOptionsGetInt(NULL,"-m",&m,NULL);
 53:   N = m*(m+1)/2;
 54:   PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n",N,m);
 55:   PetscOptionsGetScalar(NULL,"-target",&target,NULL);
 56:   SlepcSNPrintfScalar(str,50,target,PETSC_FALSE);
 57:   PetscPrintf(PETSC_COMM_WORLD,"Searching closest eigenvalues to the right of %s.\n\n",str);

 59:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 60:      Compute the operator matrix that defines the eigensystem, Ax=kx
 61:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 63:   MatCreate(PETSC_COMM_WORLD,&A);
 64:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 65:   MatSetFromOptions(A);
 66:   MatSetUp(A);
 67:   MatMarkovModel(m,A);

 69:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 70:                 Create the eigensolver and set various options
 71:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 73:   /*
 74:      Create eigensolver context
 75:   */
 76:   EPSCreate(PETSC_COMM_WORLD,&eps);

 78:   /*
 79:      Set operators. In this case, it is a standard eigenvalue problem
 80:   */
 81:   EPSSetOperators(eps,A,NULL);
 82:   EPSSetProblemType(eps,EPS_NHEP);

 84:   /*
 85:      Set the custom comparing routine in order to obtain the eigenvalues
 86:      closest to the target on the right only
 87:   */
 88:   EPSSetEigenvalueComparison(eps,MyEigenSort,&target);

 90:   /*
 91:      Set solver parameters at runtime
 92:   */
 93:   EPSSetFromOptions(eps);

 95:   /*
 96:      Set the initial vector. This is optional, if not done the initial
 97:      vector is set to random values
 98:   */
 99:   MatGetVecs(A,&v0,NULL);
100:   VecSet(v0,1.0);
101:   EPSSetInitialSpace(eps,1,&v0);

103:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104:                       Solve the eigensystem
105:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

107:   EPSSolve(eps);

109:   /*
110:      Optional: Get some information from the solver and display it
111:   */
112:   EPSGetType(eps,&type);
113:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
114:   EPSGetDimensions(eps,&nev,NULL,NULL);
115:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

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

121:   EPSPrintSolution(eps,NULL);
122:   EPSDestroy(&eps);
123:   MatDestroy(&A);
124:   VecDestroy(&v0);
125:   SlepcFinalize();
126:   return 0;
127: }

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

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

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

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

197: /*
198:     Function for user-defined eigenvalue ordering criterion.

200:     Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
201:     one of them as the preferred one according to the criterion.
202:     In this example, the preferred value is the one closest to the target,
203:     but on the right side.
204: */
205: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
206: {
207:   PetscScalar target = *(PetscScalar*)ctx;
208:   PetscReal   da,db;
209:   PetscBool   aisright,bisright;

212:   if (PetscRealPart(target) < PetscRealPart(ar)) aisright = PETSC_TRUE;
213:   else aisright = PETSC_FALSE;
214:   if (PetscRealPart(target) < PetscRealPart(br)) bisright = PETSC_TRUE;
215:   else bisright = PETSC_FALSE;
216:   if (aisright == bisright) {
217:     /* both are on the same side of the target */
218:     da = SlepcAbsEigenvalue(ar-target,ai);
219:     db = SlepcAbsEigenvalue(br-target,bi);
220:     if (da < db) *r = -1;
221:     else if (da > db) *r = 1;
222:     else *r = 0;
223:   } else if (aisright && !bisright) *r = -1; /* 'a' is on the right */
224:   else *r = 1;  /* 'b' is on the right */
225:   return(0);
226: }