Actual source code: test9.c

slepc-3.22.1 2024-10-28
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[] = "Eigenvalue problem associated with a Markov model of a random walk on a triangular grid. "
 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 the initial vector.\n\n"
 14:   "The command line options are:\n"
 15:   "  -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";

 17: #include <slepceps.h>

 19: /*
 20:    User-defined routines
 21: */
 22: PetscErrorCode MatMarkovModel(PetscInt m,Mat A);
 23: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx);

 25: /*
 26:    Check if computed eigenvectors have unit norm
 27: */
 28: PetscErrorCode CheckNormalizedVectors(EPS eps)
 29: {
 30:   PetscInt       i,nconv;
 31:   Mat            A;
 32:   Vec            xr,xi;
 33:   PetscReal      error=0.0,normr;
 34: #if !defined(PETSC_USE_COMPLEX)
 35:   PetscReal      normi;
 36: #endif

 38:   PetscFunctionBeginUser;
 39:   PetscCall(EPSGetConverged(eps,&nconv));
 40:   if (nconv>0) {
 41:     PetscCall(EPSGetOperators(eps,&A,NULL));
 42:     PetscCall(MatCreateVecs(A,&xr,&xi));
 43:     for (i=0;i<nconv;i++) {
 44:       PetscCall(EPSGetEigenvector(eps,i,xr,xi));
 45: #if defined(PETSC_USE_COMPLEX)
 46:       PetscCall(VecNorm(xr,NORM_2,&normr));
 47:       error = PetscMax(error,PetscAbsReal(normr-PetscRealConstant(1.0)));
 48: #else
 49:       PetscCall(VecNormBegin(xr,NORM_2,&normr));
 50:       PetscCall(VecNormBegin(xi,NORM_2,&normi));
 51:       PetscCall(VecNormEnd(xr,NORM_2,&normr));
 52:       PetscCall(VecNormEnd(xi,NORM_2,&normi));
 53:       error = PetscMax(error,PetscAbsReal(SlepcAbsEigenvalue(normr,normi)-PetscRealConstant(1.0)));
 54: #endif
 55:     }
 56:     PetscCall(VecDestroy(&xr));
 57:     PetscCall(VecDestroy(&xi));
 58:     if (error>100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Vectors are not normalized. Error=%g\n",(double)error));
 59:   }
 60:   PetscFunctionReturn(PETSC_SUCCESS);
 61: }

 63: int main(int argc,char **argv)
 64: {
 65:   Vec            v0;              /* initial vector */
 66:   Mat            A;               /* operator matrix */
 67:   EPS            eps;             /* eigenproblem solver context */
 68:   PetscReal      tol=0.5*PETSC_SMALL;
 69:   PetscInt       N,m=15,nev;
 70:   PetscScalar    origin=0.0;
 71:   PetscBool      flg,delay,skipnorm=PETSC_FALSE;

 73:   PetscFunctionBeginUser;
 74:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));

 76:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
 77:   N = m*(m+1)/2;
 78:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n\n",N,m));
 79:   PetscCall(PetscOptionsGetBool(NULL,NULL,"-skipnorm",&skipnorm,NULL));

 81:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 82:      Compute the operator matrix that defines the eigensystem, Ax=kx
 83:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 85:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 86:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
 87:   PetscCall(MatSetFromOptions(A));
 88:   PetscCall(MatMarkovModel(m,A));

 90:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 91:                 Create the eigensolver and set various options
 92:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 94:   /*
 95:      Create eigensolver context
 96:   */
 97:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));

 99:   /*
100:      Set operators. In this case, it is a standard eigenvalue problem
101:   */
102:   PetscCall(EPSSetOperators(eps,A,NULL));
103:   PetscCall(EPSSetProblemType(eps,EPS_NHEP));
104:   PetscCall(EPSSetTolerances(eps,tol,PETSC_CURRENT));

106:   /*
107:      Set the custom comparing routine in order to obtain the eigenvalues
108:      closest to the target on the right only
109:   */
110:   PetscCall(EPSSetEigenvalueComparison(eps,MyEigenSort,&origin));

112:   /*
113:      Set solver parameters at runtime
114:   */
115:   PetscCall(EPSSetFromOptions(eps));
116:   PetscCall(PetscObjectTypeCompare((PetscObject)eps,EPSARNOLDI,&flg));
117:   if (flg) {
118:     PetscCall(EPSArnoldiGetDelayed(eps,&delay));
119:     if (delay) PetscCall(PetscPrintf(PETSC_COMM_WORLD," Warning: delayed reorthogonalization may be unstable\n"));
120:   }

122:   /*
123:      Set the initial vector. This is optional, if not done the initial
124:      vector is set to random values
125:   */
126:   PetscCall(MatCreateVecs(A,&v0,NULL));
127:   PetscCall(VecSetValue(v0,0,-1.5,INSERT_VALUES));
128:   PetscCall(VecSetValue(v0,1,2.1,INSERT_VALUES));
129:   PetscCall(VecAssemblyBegin(v0));
130:   PetscCall(VecAssemblyEnd(v0));
131:   PetscCall(EPSSetInitialSpace(eps,1,&v0));

133:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134:                       Solve the eigensystem
135:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

137:   PetscCall(EPSSolve(eps));
138:   PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
139:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));

141:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142:                     Display solution and clean up
143:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

145:   PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
146:   if (!skipnorm) PetscCall(CheckNormalizedVectors(eps));
147:   PetscCall(EPSDestroy(&eps));
148:   PetscCall(MatDestroy(&A));
149:   PetscCall(VecDestroy(&v0));
150:   PetscCall(SlepcFinalize());
151:   return 0;
152: }

154: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
155: {
156:   const PetscReal cst = 0.5/(PetscReal)(m-1);
157:   PetscReal       pd,pu;
158:   PetscInt        Istart,Iend,i,j,jmax,ix=0;

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

188: /*
189:     Function for user-defined eigenvalue ordering criterion.

191:     Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
192:     one of them as the preferred one according to the criterion.
193:     In this example, the preferred value is the one furthest away from the origin.
194: */
195: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
196: {
197:   PetscScalar origin = *(PetscScalar*)ctx;
198:   PetscReal   d;

200:   PetscFunctionBeginUser;
201:   d = (SlepcAbsEigenvalue(br-origin,bi) - SlepcAbsEigenvalue(ar-origin,ai))/PetscMax(SlepcAbsEigenvalue(ar-origin,ai),SlepcAbsEigenvalue(br-origin,bi));
202:   *r = d > PETSC_SQRT_MACHINE_EPSILON ? 1 : (d < -PETSC_SQRT_MACHINE_EPSILON ? -1 : PetscSign(PetscRealPart(br)));
203:   PetscFunctionReturn(PETSC_SUCCESS);
204: }

206: /*TEST

208:    testset:
209:       args: -eps_nev 4
210:       output_file: output/test9_1.out
211:       filter: sed -e "s/97136/97137/g"
212:       test:
213:          suffix: 1
214:          args: -eps_type {{krylovschur arnoldi lapack}} -eps_ncv 8 -eps_max_it 300
215:       test:
216:          suffix: 1_gd
217:          args: -eps_type gd -st_pc_type none
218:       test:
219:          suffix: 1_gd2
220:          args: -eps_type gd -eps_gd_double_expansion -st_pc_type none

222:    test:
223:       suffix: 2
224:       args: -eps_balance {{none oneside twoside}} -eps_krylovschur_locking {{0 1}} -eps_nev 4 -eps_max_it 1500
225:       requires: double
226:       output_file: output/test9_1.out

228:    test:
229:       suffix: 3
230:       nsize: 2
231:       args: -eps_type arnoldi -eps_arnoldi_delayed -eps_largest_real -eps_nev 3 -eps_tol 1e-7 -bv_orthog_refine {{never ifneeded}} -skipnorm
232:       requires: !single
233:       output_file: output/test9_3.out

235:    test:
236:       suffix: 4
237:       args: -eps_nev 4 -eps_true_residual
238:       requires: !single
239:       output_file: output/test9_1.out

241:    test:
242:       suffix: 5
243:       args: -eps_type jd -eps_nev 3 -eps_target .5 -eps_harmonic -st_ksp_type bicg -st_pc_type lu -eps_jd_minv 2
244:       filter: sed -e "s/[+-]0\.0*i//g"
245:       requires: !single

247:    test:
248:       suffix: 5_arpack
249:       args: -eps_nev 3 -st_type sinvert -eps_target .5 -eps_type arpack -eps_ncv 10
250:       requires: arpack !single
251:       output_file: output/test9_5.out

253:    testset:
254:       args: -eps_type ciss -eps_tol 1e-9 -rg_type ellipse -rg_ellipse_center 0.55 -rg_ellipse_radius 0.05 -rg_ellipse_vscale 0.1 -eps_ciss_usest 0 -eps_all
255:       requires: !single
256:       output_file: output/test9_6.out
257:       test:
258:          suffix: 6
259:       test:
260:          suffix: 6_hankel
261:          args: -eps_ciss_extraction hankel -eps_ciss_spurious_threshold 1e-6 -eps_ncv 64
262:       test:
263:          suffix: 6_cheby
264:          args: -eps_ciss_quadrule chebyshev
265:       test:
266:          suffix: 6_hankel_cheby
267:          args: -eps_ciss_extraction hankel -eps_ciss_quadrule chebyshev -eps_ncv 64
268:       test:
269:          suffix: 6_refine
270:          args: -eps_ciss_moments 4 -eps_ciss_blocksize 5 -eps_ciss_refine_inner 1 -eps_ciss_refine_blocksize 2
271:       test:
272:          suffix: 6_bcgs
273:          args: -eps_ciss_realmats -eps_ciss_ksp_type bcgs -eps_ciss_pc_type ilu -eps_ciss_integration_points 8

275:    test:
276:       suffix: 6_cheby_interval
277:       args: -eps_type ciss -eps_tol 1e-9 -rg_type interval -rg_interval_endpoints 0.5,0.6 -eps_ciss_quadrule chebyshev -eps_ciss_usest 0 -eps_all
278:       requires: !single
279:       output_file: output/test9_6.out

281:    testset:
282:       args: -eps_nev 4 -eps_two_sided -eps_view_vectors ::ascii_info -eps_view_values
283:       filter: sed -e "s/\(0x[0-9a-fA-F]*\)/objectid/"
284:       test:
285:          suffix: 7_real
286:          requires: !single !complex
287:       test:
288:          suffix: 7
289:          requires: !single complex

291:    test:
292:       suffix: 8
293:       args: -eps_nev 4 -eps_ncv 7 -eps_view_values draw -eps_monitor draw::draw_lg
294:       requires: x
295:       output_file: output/test9_1.out

297:    test:
298:       suffix: 5_ksphpddm
299:       args: -eps_nev 3 -st_type sinvert -eps_target .5 -st_ksp_type hpddm -st_ksp_hpddm_type gcrodr -eps_ncv 10
300:       requires: hpddm
301:       output_file: output/test9_5.out

303:    test:
304:       suffix: 5_pchpddm
305:       args: -eps_nev 3 -st_type sinvert -eps_target .5 -st_pc_type hpddm -st_pc_hpddm_coarse_pc_type lu -eps_ncv 10
306:       requires: hpddm
307:       output_file: output/test9_5.out

309: TEST*/