Actual source code: test9.c
slepc-3.14.2 2021-02-01
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: {
31: PetscInt i,nconv;
32: Mat A;
33: Vec xr,xi;
34: PetscReal error=0.0,normr;
35: #if !defined(PETSC_USE_COMPLEX)
36: PetscReal normi;
37: #endif
40: EPSGetConverged(eps,&nconv);
41: if (nconv>0) {
42: EPSGetOperators(eps,&A,NULL);
43: MatCreateVecs(A,&xr,&xi);
44: for (i=0;i<nconv;i++) {
45: EPSGetEigenvector(eps,i,xr,xi);
46: #if defined(PETSC_USE_COMPLEX)
47: VecNorm(xr,NORM_2,&normr);
48: error = PetscMax(error,PetscAbsReal(normr-PetscRealConstant(1.0)));
49: #else
50: VecNormBegin(xr,NORM_2,&normr);
51: VecNormBegin(xi,NORM_2,&normi);
52: VecNormEnd(xr,NORM_2,&normr);
53: VecNormEnd(xi,NORM_2,&normi);
54: error = PetscMax(error,PetscAbsReal(SlepcAbsEigenvalue(normr,normi)-PetscRealConstant(1.0)));
55: #endif
56: }
57: VecDestroy(&xr);
58: VecDestroy(&xi);
59: if (error>100*PETSC_MACHINE_EPSILON) {
60: PetscPrintf(PETSC_COMM_WORLD,"Vectors are not normalized. Error=%g\n",(double)error);
61: }
62: }
63: return(0);
64: }
66: int main(int argc,char **argv)
67: {
68: Vec v0; /* initial vector */
69: Mat A; /* operator matrix */
70: EPS eps; /* eigenproblem solver context */
71: PetscReal tol=1000*PETSC_MACHINE_EPSILON;
72: PetscInt N,m=15,nev;
73: PetscScalar origin=0.0;
74: PetscBool flg,delay,skipnorm=PETSC_FALSE;
77: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
79: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
80: N = m*(m+1)/2;
81: PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n\n",N,m);
82: PetscOptionsGetBool(NULL,NULL,"-skipnorm",&skipnorm,NULL);
84: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85: Compute the operator matrix that defines the eigensystem, Ax=kx
86: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88: MatCreate(PETSC_COMM_WORLD,&A);
89: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
90: MatSetFromOptions(A);
91: MatSetUp(A);
92: MatMarkovModel(m,A);
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Create the eigensolver and set various options
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: /*
99: Create eigensolver context
100: */
101: EPSCreate(PETSC_COMM_WORLD,&eps);
103: /*
104: Set operators. In this case, it is a standard eigenvalue problem
105: */
106: EPSSetOperators(eps,A,NULL);
107: EPSSetProblemType(eps,EPS_NHEP);
108: EPSSetTolerances(eps,tol,PETSC_DEFAULT);
110: /*
111: Set the custom comparing routine in order to obtain the eigenvalues
112: closest to the target on the right only
113: */
114: EPSSetEigenvalueComparison(eps,MyEigenSort,&origin);
117: /*
118: Set solver parameters at runtime
119: */
120: EPSSetFromOptions(eps);
121: PetscObjectTypeCompare((PetscObject)eps,EPSARNOLDI,&flg);
122: if (flg) {
123: EPSArnoldiGetDelayed(eps,&delay);
124: if (delay) {
125: PetscPrintf(PETSC_COMM_WORLD," Warning: delayed reorthogonalization may be unstable\n");
126: }
127: }
129: /*
130: Set the initial vector. This is optional, if not done the initial
131: vector is set to random values
132: */
133: MatCreateVecs(A,&v0,NULL);
134: VecSetValue(v0,0,-1.5,INSERT_VALUES);
135: VecSetValue(v0,1,2.1,INSERT_VALUES);
136: VecAssemblyBegin(v0);
137: VecAssemblyEnd(v0);
138: EPSSetInitialSpace(eps,1,&v0);
140: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141: Solve the eigensystem
142: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144: EPSSolve(eps);
145: EPSGetDimensions(eps,&nev,NULL,NULL);
146: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
148: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149: Display solution and clean up
150: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152: EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
153: if (!skipnorm) { CheckNormalizedVectors(eps); }
154: EPSDestroy(&eps);
155: MatDestroy(&A);
156: VecDestroy(&v0);
157: SlepcFinalize();
158: return ierr;
159: }
161: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
162: {
163: const PetscReal cst = 0.5/(PetscReal)(m-1);
164: PetscReal pd,pu;
165: PetscInt Istart,Iend,i,j,jmax,ix=0;
166: PetscErrorCode ierr;
169: MatGetOwnershipRange(A,&Istart,&Iend);
170: for (i=1;i<=m;i++) {
171: jmax = m-i+1;
172: for (j=1;j<=jmax;j++) {
173: ix = ix + 1;
174: if (ix-1<Istart || ix>Iend) continue; /* compute only owned rows */
175: if (j!=jmax) {
176: pd = cst*(PetscReal)(i+j-1);
177: /* north */
178: if (i==1) {
179: MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
180: } else {
181: MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
182: }
183: /* east */
184: if (j==1) {
185: MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
186: } else {
187: MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
188: }
189: }
190: /* south */
191: pu = 0.5 - cst*(PetscReal)(i+j-3);
192: if (j>1) {
193: MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
194: }
195: /* west */
196: if (i>1) {
197: MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
198: }
199: }
200: }
201: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
202: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
203: return(0);
204: }
206: /*
207: Function for user-defined eigenvalue ordering criterion.
209: Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
210: one of them as the preferred one according to the criterion.
211: In this example, the preferred value is the one furthest away from the origin.
212: */
213: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
214: {
215: PetscScalar origin = *(PetscScalar*)ctx;
216: PetscReal d;
219: d = (SlepcAbsEigenvalue(br-origin,bi) - SlepcAbsEigenvalue(ar-origin,ai))/PetscMax(SlepcAbsEigenvalue(ar-origin,ai),SlepcAbsEigenvalue(br-origin,bi));
220: *r = d > PETSC_SQRT_MACHINE_EPSILON ? 1 : (d < -PETSC_SQRT_MACHINE_EPSILON ? -1 : PetscSign(PetscRealPart(br)));
221: return(0);
222: }
224: /*TEST
226: testset:
227: args: -eps_nev 4
228: requires: !single
229: output_file: output/test9_1.out
230: test:
231: suffix: 1
232: args: -eps_type {{krylovschur arnoldi lapack}} -eps_ncv 7 -eps_max_it 300
233: test:
234: suffix: 1_gd
235: args: -eps_type gd -st_pc_type none
236: test:
237: suffix: 1_gd2
238: args: -eps_type gd -eps_gd_double_expansion -st_pc_type none
240: test:
241: suffix: 2
242: args: -eps_balance {{none oneside twoside}} -eps_krylovschur_locking {{0 1}} -eps_nev 4 -eps_ncv 7 -eps_max_it 1500
243: requires: double
244: output_file: output/test9_1.out
246: test:
247: suffix: 3
248: nsize: 2
249: args: -eps_type arnoldi -eps_arnoldi_delayed -eps_largest_real -eps_nev 3 -eps_tol 1e-7 -bv_orthog_refine {{never ifneeded}} -skipnorm
250: requires: !single
251: output_file: output/test9_3.out
253: test:
254: suffix: 4
255: args: -eps_nev 4 -eps_true_residual
256: output_file: output/test9_1.out
258: test:
259: suffix: 5
260: args: -eps_type jd -eps_nev 3 -eps_target .5 -eps_harmonic -st_ksp_type bicg -st_pc_type lu -eps_jd_minv 2
261: requires: !single
263: test:
264: suffix: 5_arpack
265: args: -eps_nev 3 -st_type sinvert -eps_target .5 -eps_type arpack -eps_ncv 10
266: requires: arpack
267: output_file: output/test9_5.out
269: testset:
270: 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
271: requires: double
272: output_file: output/test9_6.out
273: test:
274: suffix: 6
275: test:
276: suffix: 6_hankel
277: args: -eps_ciss_extraction hankel -eps_ciss_spurious_threshold 1e-6 -eps_ncv 64
278: test:
279: suffix: 6_cheby
280: args: -eps_ciss_quadrule chebyshev
281: test:
282: suffix: 6_hankel_cheby
283: args: -eps_ciss_extraction hankel -eps_ciss_quadrule chebyshev -eps_ncv 64
284: test:
285: suffix: 6_refine
286: args: -eps_ciss_refine_inner 1 -eps_ciss_refine_blocksize 1
287: test:
288: suffix: 6_bcgs
289: args: -eps_ciss_realmats -eps_ciss_ksp_type bcgs -eps_ciss_pc_type sor -eps_ciss_integration_points 12
291: test:
292: suffix: 7
293: args: -eps_nev 4 -eps_two_sided -eps_view_vectors ::ascii_info -eps_view_values
294: requires: !single
295: filter: sed -e "s/\(0x[0-9a-fA-F]*\)/objectid/"
297: test:
298: suffix: 8
299: args: -eps_nev 4 -eps_view_values draw -eps_monitor_lg
300: requires: x
301: output_file: output/test9_1.out
303: test:
304: suffix: 5_ksphpddm
305: args: -eps_nev 3 -st_type sinvert -eps_target .5 -st_ksp_type hpddm -st_ksp_hpddm_type gcrodr -eps_ncv 10
306: requires: hpddm
307: output_file: output/test9_5.out
309: test:
310: suffix: 5_pchpddm
311: args: -eps_nev 3 -st_type sinvert -eps_target .5 -st_pc_type hpddm -st_pc_hpddm_coarse_pc_type lu -eps_ncv 10
312: requires: hpddm
313: output_file: output/test9_5.out
315: TEST*/