Actual source code: ex41.c

slepc-main 2024-11-09
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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[] = "Illustrates the computation of left eigenvectors.\n\n"
 12:   "The problem is the Markov model as in ex5.c.\n"
 13:   "The command line options are:\n"
 14:   "  -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";

 16: #include <slepceps.h>

 18: /*
 19:    User-defined routines
 20: */
 21: PetscErrorCode MatMarkovModel(PetscInt,Mat);
 22: PetscErrorCode ComputeResidualNorm(Mat,PetscBool,PetscScalar,PetscScalar,Vec,Vec,Vec,PetscReal*);

 24: int main(int argc,char **argv)
 25: {
 26:   Vec            v0,w0;           /* initial vectors */
 27:   Mat            A;               /* operator matrix */
 28:   EPS            eps;             /* eigenproblem solver context */
 29:   EPSType        type;
 30:   PetscInt       i,N,m=15,nconv;
 31:   PetscBool      twosided;
 32:   PetscReal      nrmr,nrml=0.0,re,im,lev;
 33:   PetscScalar    *kr,*ki;
 34:   Vec            t,*xr,*xi,*yr,*yi;
 35:   PetscMPIInt    rank;

 37:   PetscFunctionBeginUser;
 38:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));

 40:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
 41:   N = m*(m+1)/2;
 42:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n\n",N,m));

 44:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 45:      Compute the operator matrix that defines the eigensystem, Ax=kx
 46:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 48:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 49:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
 50:   PetscCall(MatSetFromOptions(A));
 51:   PetscCall(MatMarkovModel(m,A));
 52:   PetscCall(MatCreateVecs(A,NULL,&t));

 54:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 55:                 Create the eigensolver and set various options
 56:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 58:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
 59:   PetscCall(EPSSetOperators(eps,A,NULL));
 60:   PetscCall(EPSSetProblemType(eps,EPS_NHEP));

 62:   /* use a two-sided algorithm to compute left eigenvectors as well */
 63:   PetscCall(EPSSetTwoSided(eps,PETSC_TRUE));

 65:   /* allow user to change settings at run time */
 66:   PetscCall(EPSSetFromOptions(eps));
 67:   PetscCall(EPSGetTwoSided(eps,&twosided));

 69:   /*
 70:      Set the initial vectors. This is optional, if not done the initial
 71:      vectors are set to random values
 72:   */
 73:   PetscCall(MatCreateVecs(A,&v0,&w0));
 74:   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
 75:   if (!rank) {
 76:     PetscCall(VecSetValue(v0,0,1.0,INSERT_VALUES));
 77:     PetscCall(VecSetValue(v0,1,1.0,INSERT_VALUES));
 78:     PetscCall(VecSetValue(v0,2,1.0,INSERT_VALUES));
 79:     PetscCall(VecSetValue(w0,0,2.0,INSERT_VALUES));
 80:     PetscCall(VecSetValue(w0,2,0.5,INSERT_VALUES));
 81:   }
 82:   PetscCall(VecAssemblyBegin(v0));
 83:   PetscCall(VecAssemblyBegin(w0));
 84:   PetscCall(VecAssemblyEnd(v0));
 85:   PetscCall(VecAssemblyEnd(w0));
 86:   PetscCall(EPSSetInitialSpace(eps,1,&v0));
 87:   PetscCall(EPSSetLeftInitialSpace(eps,1,&w0));

 89:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 90:                       Solve the eigensystem
 91:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 93:   PetscCall(EPSSolve(eps));

 95:   /*
 96:      Optional: Get some information from the solver and display it
 97:   */
 98:   PetscCall(EPSGetType(eps,&type));
 99:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type));

101:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102:                     Display solution and clean up
103:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

105:   /*
106:      Get number of converged approximate eigenpairs
107:   */
108:   PetscCall(EPSGetConverged(eps,&nconv));
109:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of converged eigenpairs: %" PetscInt_FMT "\n\n",nconv));
110:   PetscCall(PetscMalloc2(nconv,&kr,nconv,&ki));
111:   PetscCall(VecDuplicateVecs(t,nconv,&xr));
112:   PetscCall(VecDuplicateVecs(t,nconv,&xi));
113:   if (twosided) {
114:     PetscCall(VecDuplicateVecs(t,nconv,&yr));
115:     PetscCall(VecDuplicateVecs(t,nconv,&yi));
116:   }

118:   if (nconv>0) {
119:     /*
120:        Display eigenvalues and relative errors
121:     */
122:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,
123:          "           k            ||Ax-kx||         ||y'A-y'k||\n"
124:          "   ---------------- ------------------ ------------------\n"));

126:     for (i=0;i<nconv;i++) {
127:       /*
128:         Get converged eigenpairs: i-th eigenvalue is stored in kr (real part) and
129:         ki (imaginary part)
130:       */
131:       PetscCall(EPSGetEigenpair(eps,i,&kr[i],&ki[i],xr[i],xi[i]));
132:       if (twosided) PetscCall(EPSGetLeftEigenvector(eps,i,yr[i],yi[i]));
133:       /*
134:          Compute the residual norms associated to each eigenpair
135:       */
136:       PetscCall(ComputeResidualNorm(A,PETSC_FALSE,kr[i],ki[i],xr[i],xi[i],t,&nrmr));
137:       if (twosided) PetscCall(ComputeResidualNorm(A,PETSC_TRUE,kr[i],ki[i],yr[i],yi[i],t,&nrml));

139: #if defined(PETSC_USE_COMPLEX)
140:       re = PetscRealPart(kr[i]);
141:       im = PetscImaginaryPart(kr[i]);
142: #else
143:       re = kr[i];
144:       im = ki[i];
145: #endif
146:       if (im!=0.0) PetscCall(PetscPrintf(PETSC_COMM_WORLD," %8f%+8fi %12g %12g\n",(double)re,(double)im,(double)nrmr,(double)nrml));
147:       else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"   %12f       %12g       %12g\n",(double)re,(double)nrmr,(double)nrml));
148:     }
149:     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n"));
150:     /*
151:        Check bi-orthogonality of eigenvectors
152:     */
153:     if (twosided) {
154:       PetscCall(VecCheckOrthogonality(xr,nconv,yr,nconv,NULL,NULL,&lev));
155:       if (lev<100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"  Level of bi-orthogonality of eigenvectors < 100*eps\n\n"));
156:       else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"  Level of bi-orthogonality of eigenvectors: %g\n\n",(double)lev));
157:     }
158:   }

160:   PetscCall(EPSDestroy(&eps));
161:   PetscCall(MatDestroy(&A));
162:   PetscCall(VecDestroy(&v0));
163:   PetscCall(VecDestroy(&w0));
164:   PetscCall(VecDestroy(&t));
165:   PetscCall(PetscFree2(kr,ki));
166:   PetscCall(VecDestroyVecs(nconv,&xr));
167:   PetscCall(VecDestroyVecs(nconv,&xi));
168:   if (twosided) {
169:     PetscCall(VecDestroyVecs(nconv,&yr));
170:     PetscCall(VecDestroyVecs(nconv,&yi));
171:   }
172:   PetscCall(SlepcFinalize());
173:   return 0;
174: }

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

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

192:     Note: the code will actually compute the transpose of the stochastic matrix
193:     that contains the transition probabilities.
194: */
195: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
196: {
197:   const PetscReal cst = 0.5/(PetscReal)(m-1);
198:   PetscReal       pd,pu;
199:   PetscInt        Istart,Iend,i,j,jmax,ix=0;

201:   PetscFunctionBeginUser;
202:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
203:   for (i=1;i<=m;i++) {
204:     jmax = m-i+1;
205:     for (j=1;j<=jmax;j++) {
206:       ix = ix + 1;
207:       if (ix-1<Istart || ix>Iend) continue;  /* compute only owned rows */
208:       if (j!=jmax) {
209:         pd = cst*(PetscReal)(i+j-1);
210:         /* north */
211:         if (i==1) PetscCall(MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES));
212:         else PetscCall(MatSetValue(A,ix-1,ix,pd,INSERT_VALUES));
213:         /* east */
214:         if (j==1) PetscCall(MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES));
215:         else PetscCall(MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES));
216:       }
217:       /* south */
218:       pu = 0.5 - cst*(PetscReal)(i+j-3);
219:       if (j>1) PetscCall(MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES));
220:       /* west */
221:       if (i>1) PetscCall(MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES));
222:     }
223:   }
224:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
225:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
226:   PetscFunctionReturn(PETSC_SUCCESS);
227: }

229: /*
230:    ComputeResidualNorm - Computes the norm of the residual vector
231:    associated with an eigenpair.

233:    Input Parameters:
234:      trans - whether A' must be used instead of A
235:      kr,ki - eigenvalue
236:      xr,xi - eigenvector
237:      u     - work vector
238: */
239: PetscErrorCode ComputeResidualNorm(Mat A,PetscBool trans,PetscScalar kr,PetscScalar ki,Vec xr,Vec xi,Vec u,PetscReal *norm)
240: {
241: #if !defined(PETSC_USE_COMPLEX)
242:   PetscReal      ni,nr;
243: #endif
244:   PetscErrorCode (*matmult)(Mat,Vec,Vec) = trans? MatMultTranspose: MatMult;

246:   PetscFunctionBegin;
247: #if !defined(PETSC_USE_COMPLEX)
248:   if (ki == 0 || PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
249: #endif
250:     PetscCall((*matmult)(A,xr,u));
251:     if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) PetscCall(VecAXPY(u,-kr,xr));
252:     PetscCall(VecNorm(u,NORM_2,norm));
253: #if !defined(PETSC_USE_COMPLEX)
254:   } else {
255:     PetscCall((*matmult)(A,xr,u));
256:     if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
257:       PetscCall(VecAXPY(u,-kr,xr));
258:       PetscCall(VecAXPY(u,ki,xi));
259:     }
260:     PetscCall(VecNorm(u,NORM_2,&nr));
261:     PetscCall((*matmult)(A,xi,u));
262:     if (SlepcAbsEigenvalue(kr,ki) > PETSC_MACHINE_EPSILON) {
263:       PetscCall(VecAXPY(u,-kr,xi));
264:       PetscCall(VecAXPY(u,-ki,xr));
265:     }
266:     PetscCall(VecNorm(u,NORM_2,&ni));
267:     *norm = SlepcAbsEigenvalue(nr,ni);
268:   }
269: #endif
270:   PetscFunctionReturn(PETSC_SUCCESS);
271: }

273: /*TEST

275:    testset:
276:       args: -st_type sinvert -eps_target 1.1 -eps_nev 4
277:       filter: grep -v method | sed -e "s/[+-]0\.0*i//g" | sed -e "s/[0-9]\.[0-9]*e[+-]\([0-9]*\)/removed/g"
278:       requires: !single
279:       output_file: output/ex41_1.out
280:       test:
281:          suffix: 1
282:          args: -eps_type {{power krylovschur}}
283:       test:
284:          suffix: 1_balance
285:          args: -eps_balance {{oneside twoside}} -eps_ncv 17 -eps_krylovschur_locking 0
286:          requires: !__float128

288: TEST*/