Actual source code: ex18.c

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, 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:   Mat            A;               /* operator matrix */
 42:   EPS            eps;             /* eigenproblem solver context */
 43:   EPSType        type;
 44:   PetscScalar    target=0.5;
 45:   PetscInt       N,m=15,nev;
 46:   PetscBool      terse;
 47:   PetscViewer    viewer;
 49:   char           str[50];

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

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

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

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

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

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

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

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

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

 96:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 97:                       Solve the eigensystem
 98:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

100:   EPSSolve(eps);

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

110:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
111:                     Display solution and clean up
112:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

114:   /* show detailed info unless -terse option is given by user */
115:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
116:   if (terse) {
117:     EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
118:   } else {
119:     PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
120:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
121:     EPSReasonView(eps,viewer);
122:     EPSErrorView(eps,EPS_ERROR_RELATIVE,viewer);
123:     PetscViewerPopFormat(viewer);
124:   }
125:   EPSDestroy(&eps);
126:   MatDestroy(&A);
127:   SlepcFinalize();
128:   return ierr;
129: }

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

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

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

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

199: /*
200:     Function for user-defined eigenvalue ordering criterion.

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

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