Actual source code: ex29.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: static char help[] = "Solves the same problem as in ex5, with a user-defined stopping test."
 12:   "It is a standard nonsymmetric eigenproblem with real eigenvalues and the rightmost eigenvalue is known to be 1.\n"
 13:   "This example illustrates how the user can set a custom stopping test function.\n\n"
 14:   "The command line options are:\n"
 15:   "  -m <m>, where <m> = number of grid subdivisions in each dimension.\n"
 16:   "  -seconds <s>, where <s> = maximum time in seconds allowed for computation.\n\n";

 18: #include <slepceps.h>
 19: #include <petsctime.h>

 21: /*
 22:    User-defined routines
 23: */

 25: PetscErrorCode MyStoppingTest(EPS,PetscInt,PetscInt,PetscInt,PetscInt,EPSConvergedReason*,void*);
 26: PetscErrorCode MatMarkovModel(PetscInt,Mat);

 28: int main(int argc,char **argv)
 29: {
 30:   Mat                A;               /* operator matrix */
 31:   EPS                eps;             /* eigenproblem solver context */
 32:   PetscReal          seconds=2.5;     /* maximum time allowed for computation */
 33:   PetscLogDouble     deadline;        /* time to abort computation */
 34:   PetscInt           N,m=15,nconv;
 35:   PetscBool          terse;
 36:   PetscViewer        viewer;
 37:   EPSConvergedReason reason;
 38:   PetscErrorCode     ierr;

 40:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 42:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 43:   N = m*(m+1)/2;
 44:   PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n",N,m);
 45:   PetscOptionsGetReal(NULL,NULL,"-seconds",&seconds,NULL);
 46:   PetscPrintf(PETSC_COMM_WORLD,"Maximum time for computation is set to %g seconds.\n\n",(double)seconds);
 47:   deadline = seconds;
 48:   PetscTimeAdd(&deadline);

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

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

 60:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 61:                 Create the eigensolver and set various options
 62:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 64:   EPSCreate(PETSC_COMM_WORLD,&eps);
 65:   EPSSetOperators(eps,A,NULL);
 66:   EPSSetProblemType(eps,EPS_NHEP);
 67:   EPSSetStoppingTestFunction(eps,MyStoppingTest,&deadline,NULL);
 68:   EPSSetFromOptions(eps);

 70:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 71:                       Solve the eigensystem
 72:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 74:   EPSSolve(eps);

 76:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 77:                     Display solution and clean up
 78:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 80:   /* show detailed info unless -terse option is given by user */
 81:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
 82:   if (terse) {
 83:     EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
 84:   } else {
 85:     PetscViewerASCIIGetStdout(PETSC_COMM_WORLD,&viewer);
 86:     PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO_DETAIL);
 87:     EPSGetConvergedReason(eps,&reason);
 88:     if (reason!=EPS_CONVERGED_USER) {
 89:       EPSReasonView(eps,viewer);
 90:       EPSErrorView(eps,EPS_ERROR_RELATIVE,viewer);
 91:     } else {
 92:       EPSGetConverged(eps,&nconv);
 93:       PetscViewerASCIIPrintf(viewer,"Eigensolve finished with %D converged eigenpairs; reason=%s\n",nconv,EPSConvergedReasons[reason]);
 94:     }
 95:     PetscViewerPopFormat(viewer);
 96:   }
 97:   EPSDestroy(&eps);
 98:   MatDestroy(&A);
 99:   SlepcFinalize();
100:   return ierr;
101: }

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

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

119:     Note: the code will actually compute the transpose of the stochastic matrix
120:     that contains the transition probabilities.
121: */
122: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
123: {
124:   const PetscReal cst = 0.5/(PetscReal)(m-1);
125:   PetscReal       pd,pu;
126:   PetscInt        Istart,Iend,i,j,jmax,ix=0;
127:   PetscErrorCode  ierr;

130:   MatGetOwnershipRange(A,&Istart,&Iend);
131:   for (i=1;i<=m;i++) {
132:     jmax = m-i+1;
133:     for (j=1;j<=jmax;j++) {
134:       ix = ix + 1;
135:       if (ix-1<Istart || ix>Iend) continue;  /* compute only owned rows */
136:       if (j!=jmax) {
137:         pd = cst*(PetscReal)(i+j-1);
138:         /* north */
139:         if (i==1) {
140:           MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
141:         } else {
142:           MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
143:         }
144:         /* east */
145:         if (j==1) {
146:           MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
147:         } else {
148:           MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
149:         }
150:       }
151:       /* south */
152:       pu = 0.5 - cst*(PetscReal)(i+j-3);
153:       if (j>1) {
154:         MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
155:       }
156:       /* west */
157:       if (i>1) {
158:         MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
159:       }
160:     }
161:   }
162:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
163:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
164:   return(0);
165: }

167: /*
168:     Function for user-defined stopping test.

170:     Checks that the computing time has not exceeded the deadline.
171: */
172: PetscErrorCode MyStoppingTest(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
173: {
175:   PetscLogDouble now,deadline = *(PetscLogDouble*)ctx;

178:   /* check if usual termination conditions are met */
179:   EPSStoppingBasic(eps,its,max_it,nconv,nev,reason,NULL);
180:   if (*reason==EPS_CONVERGED_ITERATING) {
181:     /* check if deadline has expired */
182:     PetscTime(&now);
183:     if (now>deadline) *reason = EPS_CONVERGED_USER;
184:   }
185:   return(0);
186: }

188: /*TEST

190:    test:
191:       suffix: 1
192:       args: -m 350

194: TEST*/