Actual source code: test1.c
slepc-3.21.2 2024-09-25
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[] = "Tests B-orthonormality of eigenvectors in a GHEP problem.\n\n";
13: #include <slepceps.h>
15: int main(int argc,char **argv)
16: {
17: Mat A,B; /* matrices */
18: EPS eps; /* eigenproblem solver context */
19: ST st;
20: Vec *X,v;
21: PetscReal lev=0.0,tol=PETSC_SMALL;
22: PetscInt N,n=45,m,Istart,Iend,II,i,j,nconv;
23: PetscBool flag,skiporth=PETSC_FALSE;
24: EPSPowerShiftType variant;
26: PetscFunctionBeginUser;
27: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
28: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
29: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag));
30: if (!flag) m=n;
31: N = n*m;
32: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized Symmetric Eigenproblem, N=%" PetscInt_FMT " (%" PetscInt_FMT "x%" PetscInt_FMT " grid)\n\n",N,n,m));
33: PetscCall(PetscOptionsGetBool(NULL,NULL,"-skiporth",&skiporth,NULL));
35: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36: Compute the matrices that define the eigensystem, Ax=kBx
37: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
39: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
40: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
41: PetscCall(MatSetFromOptions(A));
43: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
44: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N));
45: PetscCall(MatSetFromOptions(B));
47: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
48: for (II=Istart;II<Iend;II++) {
49: i = II/n; j = II-i*n;
50: if (i>0) PetscCall(MatSetValue(A,II,II-n,-1.0,INSERT_VALUES));
51: if (i<m-1) PetscCall(MatSetValue(A,II,II+n,-1.0,INSERT_VALUES));
52: if (j>0) PetscCall(MatSetValue(A,II,II-1,-1.0,INSERT_VALUES));
53: if (j<n-1) PetscCall(MatSetValue(A,II,II+1,-1.0,INSERT_VALUES));
54: PetscCall(MatSetValue(A,II,II,4.0,INSERT_VALUES));
55: PetscCall(MatSetValue(B,II,II,2.0/PetscLogScalar(II+2),INSERT_VALUES));
56: }
58: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
59: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
60: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
61: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
62: PetscCall(MatCreateVecs(B,&v,NULL));
64: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
65: Create the eigensolver and set various options
66: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68: PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
69: PetscCall(EPSSetOperators(eps,A,B));
70: PetscCall(EPSSetProblemType(eps,EPS_GHEP));
71: PetscCall(EPSSetTolerances(eps,tol,PETSC_DEFAULT));
72: PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_NORM));
73: PetscCall(EPSSetFromOptions(eps));
75: /* illustrate how to extract parameters from specific solver types */
76: PetscCall(PetscObjectTypeCompare((PetscObject)eps,EPSPOWER,&flag));
77: if (flag) {
78: PetscCall(EPSGetST(eps,&st));
79: PetscCall(PetscObjectTypeCompare((PetscObject)st,STSHIFT,&flag));
80: if (flag) {
81: PetscCall(EPSPowerGetShiftType(eps,&variant));
82: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Type of shifts used during power iteration: %s\n",EPSPowerShiftTypes[variant]));
83: }
84: }
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Solve the eigensystem
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90: PetscCall(EPSSolve(eps));
92: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93: Display solution and clean up
94: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: PetscCall(EPSGetTolerances(eps,&tol,NULL));
97: PetscCall(EPSErrorView(eps,EPS_ERROR_BACKWARD,NULL));
98: PetscCall(EPSGetConverged(eps,&nconv));
99: if (nconv>1) {
100: PetscCall(VecDuplicateVecs(v,nconv,&X));
101: for (i=0;i<nconv;i++) PetscCall(EPSGetEigenvector(eps,i,X[i],NULL));
102: if (!skiporth) PetscCall(VecCheckOrthonormality(X,nconv,NULL,nconv,B,NULL,&lev));
103: if (lev<10*tol) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality below the tolerance\n"));
104: else PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Level of orthogonality: %g\n",(double)lev));
105: PetscCall(VecDestroyVecs(nconv,&X));
106: }
108: PetscCall(EPSDestroy(&eps));
109: PetscCall(MatDestroy(&A));
110: PetscCall(MatDestroy(&B));
111: PetscCall(VecDestroy(&v));
112: PetscCall(SlepcFinalize());
113: return 0;
114: }
116: /*TEST
118: testset:
119: args: -n 18 -eps_nev 4 -eps_max_it 1500
120: requires: !single
121: output_file: output/test1_1.out
122: test:
123: suffix: 1
124: args: -eps_type {{krylovschur arnoldi gd jd lapack}}
125: test:
126: suffix: 1_subspace
127: args: -eps_type subspace -eps_conv_rel
128: test:
129: suffix: 1_ks_nopurify
130: args: -eps_purify 0
131: test:
132: suffix: 1_ks_trueres
133: args: -eps_true_residual
134: test:
135: suffix: 1_ks_sinvert
136: args: -st_type sinvert -eps_target 22
137: test:
138: suffix: 1_ks_cayley
139: args: -st_type cayley -eps_target 22
140: test:
141: suffix: 1_lanczos
142: args: -eps_type lanczos -eps_lanczos_reorthog full
143: test:
144: suffix: 1_gd2
145: args: -eps_type gd -eps_gd_double_expansion
146: test:
147: suffix: 1_gd_borth
148: args: -eps_type gd -eps_gd_borth
149: test:
150: suffix: 1_jd_borth
151: args: -eps_type jd -eps_jd_borth
152: test:
153: suffix: 1_lobpcg
154: args: -eps_type lobpcg -st_shift 22 -eps_largest_real
155: test:
156: suffix: 1_hpddm
157: requires: hpddm
158: args: -eps_type lobpcg -st_shift 22 -eps_largest_real -st_pc_type lu -st_ksp_type hpddm
159: test:
160: suffix: 1_cholesky
161: args: -mat_type sbaij
162: test:
163: suffix: 1_scalapack
164: nsize: {{1 2 3}}
165: requires: scalapack
166: args: -eps_type scalapack
167: test:
168: suffix: 1_elpa
169: nsize: {{1 2 3}}
170: requires: elpa
171: args: -eps_type elpa
172: filter: grep -v "Buffering level"
173: test:
174: suffix: 1_elemental
175: nsize: {{1 2}}
176: requires: elemental
177: args: -eps_type elemental
179: testset:
180: args: -n 18 -eps_type ciss -rg_interval_endpoints 20.8,22
181: requires: !single
182: output_file: output/test1_1_ciss.out
183: test:
184: suffix: 1_ciss
185: args: -eps_ciss_extraction {{ritz hankel}}
186: test:
187: suffix: 1_ciss_ksps
188: args: -eps_ciss_usest 0 -eps_ciss_integration_points 12
189: test:
190: suffix: 1_ciss_gnhep
191: args: -eps_gen_non_hermitian -skiporth
192: test:
193: suffix: 1_ciss_trapezoidal
194: args: -eps_ciss_quadrule trapezoidal -eps_ciss_integration_points 24 -eps_ciss_extraction hankel -eps_ciss_delta 1e-10 -eps_tol 5e-11 -skiporth
195: test:
196: suffix: 1_ciss_cuda
197: args: -mat_type aijcusparse -st_pc_factor_mat_solver_type cusparse
198: requires: cuda
200: testset:
201: requires: !single
202: args: -eps_tol 1e-10 -st_type sinvert -st_ksp_type preonly -st_pc_type cholesky
203: test:
204: suffix: 2
205: args: -eps_interval .1,1.1
206: test:
207: suffix: 2_open
208: args: -eps_interval -inf,1.1
209: test:
210: suffix: 2_parallel
211: requires: mumps !complex
212: nsize: 3
213: args: -eps_interval .1,1.1 -eps_krylovschur_partitions 2 -st_pc_factor_mat_solver_type mumps -st_mat_mumps_icntl_13 1
214: output_file: output/test1_2.out
216: test:
217: suffix: 3
218: requires: !single
219: args: -n 18 -eps_type power -eps_conv_rel -eps_nev 3
221: test:
222: suffix: 4
223: requires: !single
224: args: -n 18 -eps_type power -eps_conv_rel -eps_nev 3 -st_type sinvert -eps_target 1.149 -eps_power_shift_type {{constant rayleigh wilkinson}}
226: testset:
227: args: -n 18 -eps_nev 3 -eps_smallest_real -eps_max_it 500 -st_pc_type icc
228: output_file: output/test1_5.out
229: test:
230: suffix: 5_rqcg
231: args: -eps_type rqcg
232: test:
233: suffix: 5_lobpcg
234: args: -eps_type lobpcg -eps_lobpcg_blocksize 3
235: test:
236: suffix: 5_hpddm
237: args: -eps_type lobpcg -eps_lobpcg_blocksize 3 -st_pc_type lu -st_ksp_type hpddm
238: requires: hpddm
239: test:
240: suffix: 5_blopex
241: args: -eps_type blopex -eps_conv_abs -st_shift 0.1
242: requires: blopex
244: testset:
245: args: -n 18 -eps_nev 12 -eps_mpd 8 -eps_max_it 3000
246: requires: !single
247: output_file: output/test1_6.out
248: test:
249: suffix: 6
250: args: -eps_type {{krylovschur arnoldi gd}}
251: test:
252: suffix: 6_lanczos
253: args: -eps_type lanczos -eps_lanczos_reorthog full
254: test:
255: suffix: 6_subspace
256: args: -eps_type subspace -eps_conv_rel
258: testset:
259: args: -n 18 -eps_nev 4 -eps_max_it 1500 -mat_type aijcusparse
260: requires: cuda !single
261: output_file: output/test1_1.out
262: test:
263: suffix: 7
264: args: -eps_type {{krylovschur arnoldi gd jd}}
265: test:
266: suffix: 7_subspace
267: args: -eps_type subspace -eps_conv_rel
268: test:
269: suffix: 7_ks_sinvert
270: args: -st_type sinvert -eps_target 22
271: test:
272: suffix: 7_lanczos
273: args: -eps_type lanczos -eps_lanczos_reorthog full
274: test:
275: suffix: 7_ciss
276: args: -eps_type ciss -rg_interval_endpoints 20.8,22 -st_pc_factor_mat_solver_type cusparse
277: output_file: output/test1_1_ciss.out
279: testset:
280: args: -n 18 -eps_nev 3 -eps_smallest_real -eps_max_it 500 -st_pc_type sor -mat_type aijcusparse
281: requires: cuda
282: output_file: output/test1_5.out
283: test:
284: suffix: 8_rqcg
285: args: -eps_type rqcg
286: test:
287: suffix: 8_lobpcg
288: args: -eps_type lobpcg -eps_lobpcg_blocksize 3
290: testset:
291: nsize: 2
292: args: -n 18 -eps_nev 7 -eps_ncv 32 -ds_parallel synchronized
293: filter: grep -v "orthogonality" | sed -e "s/[+-]0\.0*i//g" | sed -e "s/0.61338/0.61339/g"
294: output_file: output/test1_9.out
295: test:
296: suffix: 9_ks_ghep
297: args: -eps_gen_hermitian -st_pc_type redundant -st_type sinvert
298: test:
299: suffix: 9_ks_gnhep
300: args: -eps_gen_non_hermitian -st_pc_type redundant -st_type sinvert
301: test:
302: suffix: 9_ks_ghiep
303: args: -eps_gen_indefinite -st_pc_type redundant -st_type sinvert
304: requires: !single
305: test:
306: suffix: 9_lobpcg_ghep
307: args: -eps_gen_hermitian -eps_type lobpcg -eps_max_it 200 -eps_lobpcg_blocksize 6
308: requires: !single
309: timeoutfactor: 2
310: test:
311: suffix: 9_jd_gnhep
312: args: -eps_gen_non_hermitian -eps_type jd -eps_target 0 -eps_ncv 64
313: requires: !single
314: timeoutfactor: 2
316: test:
317: suffix: 10_feast
318: args: -n 25 -eps_type feast -eps_interval .95,1.1 -eps_conv_rel -eps_tol 1e-6
319: requires: feast
321: TEST*/