Actual source code: ex18.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[] = "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;
 47:   PetscBool      terse;
 48:   PetscViewer    viewer;
 50:   char           str[50];

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

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

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

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

 71:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 72:                 Create the eigensolver and set various options
 73:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

 92:   /*
 93:      Set solver parameters at runtime
 94:   */
 95:   EPSSetFromOptions(eps);

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

105:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106:                       Solve the eigensystem
107:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

109:   EPSSolve(eps);

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

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

123:   /* show detailed info unless -terse option is given by user */
124:   PetscOptionsHasName(NULL,"-terse",&terse);
125:   if (terse) {
126:     EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
127:   } else {
128:     PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
129:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
130:     EPSReasonView(eps,viewer);
131:     EPSErrorView(eps,EPS_ERROR_RELATIVE,viewer);
132:     PetscViewerPopFormat(viewer);
133:   }
134:   EPSDestroy(&eps);
135:   MatDestroy(&A);
136:   VecDestroy(&v0);
137:   SlepcFinalize();
138:   return 0;
139: }

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

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

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

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

209: /*
210:     Function for user-defined eigenvalue ordering criterion.

212:     Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
213:     one of them as the preferred one according to the criterion.
214:     In this example, the preferred value is the one closest to the target,
215:     but on the right side.
216: */
217: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
218: {
219:   PetscScalar target = *(PetscScalar*)ctx;
220:   PetscReal   da,db;
221:   PetscBool   aisright,bisright;

224:   if (PetscRealPart(target) < PetscRealPart(ar)) aisright = PETSC_TRUE;
225:   else aisright = PETSC_FALSE;
226:   if (PetscRealPart(target) < PetscRealPart(br)) bisright = PETSC_TRUE;
227:   else bisright = PETSC_FALSE;
228:   if (aisright == bisright) {
229:     /* both are on the same side of the target */
230:     da = SlepcAbsEigenvalue(ar-target,ai);
231:     db = SlepcAbsEigenvalue(br-target,bi);
232:     if (da < db) *r = -1;
233:     else if (da > db) *r = 1;
234:     else *r = 0;
235:   } else if (aisright && !bisright) *r = -1; /* 'a' is on the right */
236:   else *r = 1;  /* 'b' is on the right */
237:   return(0);
238: }