Actual source code: test2.c
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: */
10: /*
11: Example based on spring problem in NLEVP collection [1]. See the parameters
12: meaning at Example 2 in [2].
14: [1] T. Betcke, N. J. Higham, V. Mehrmann, C. Schroder, and F. Tisseur,
15: NLEVP: A Collection of Nonlinear Eigenvalue Problems, MIMS EPrint
16: 2010.98, November 2010.
17: [2] F. Tisseur, Backward error and condition of polynomial eigenvalue
18: problems, Linear Algebra and its Applications, 309 (2000), pp. 339--361,
19: April 2000.
20: */
22: static char help[] = "Test the solution of a PEP from a finite element model of damped mass-spring system.\n\n"
23: "Problem from NLEVP collection.\n"
24: "The command line options are:\n"
25: " -n <n> ... number of grid subdivisions.\n"
26: " -mu <value> ... mass (default 1).\n"
27: " -tau <value> ... damping constant of the dampers (default 10).\n"
28: " -kappa <value> ... damping constant of the springs (default 5).\n"
29: " -initv ... set an initial vector.\n\n";
31: #include <slepcpep.h>
33: /*
34: Check if computed eigenvectors have unit norm
35: */
36: PetscErrorCode CheckNormalizedVectors(PEP pep)
37: {
38: PetscInt i,nconv;
39: Mat A;
40: Vec xr,xi;
41: PetscReal error=0.0,normr;
42: #if !defined(PETSC_USE_COMPLEX)
43: PetscReal normi;
44: #endif
46: PetscFunctionBeginUser;
47: PetscCall(PEPGetConverged(pep,&nconv));
48: if (nconv>0) {
49: PetscCall(PEPGetOperators(pep,0,&A));
50: PetscCall(MatCreateVecs(A,&xr,&xi));
51: for (i=0;i<nconv;i++) {
52: PetscCall(PEPGetEigenpair(pep,i,NULL,NULL,xr,xi));
53: #if defined(PETSC_USE_COMPLEX)
54: PetscCall(VecNorm(xr,NORM_2,&normr));
55: error = PetscMax(error,PetscAbsReal(normr-PetscRealConstant(1.0)));
56: #else
57: PetscCall(VecNormBegin(xr,NORM_2,&normr));
58: PetscCall(VecNormBegin(xi,NORM_2,&normi));
59: PetscCall(VecNormEnd(xr,NORM_2,&normr));
60: PetscCall(VecNormEnd(xi,NORM_2,&normi));
61: error = PetscMax(error,PetscAbsReal(SlepcAbsEigenvalue(normr,normi)-PetscRealConstant(1.0)));
62: #endif
63: }
64: PetscCall(VecDestroy(&xr));
65: PetscCall(VecDestroy(&xi));
66: if (error>100*PETSC_MACHINE_EPSILON) PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Vectors are not normalized. Error=%g\n",(double)error));
67: }
68: PetscFunctionReturn(PETSC_SUCCESS);
69: }
71: int main(int argc,char **argv)
72: {
73: Mat M,C,K,A[3]; /* problem matrices */
74: PEP pep; /* polynomial eigenproblem solver context */
75: PetscInt n=30,Istart,Iend,i,nev;
76: PetscReal mu=1.0,tau=10.0,kappa=5.0;
77: PetscBool initv=PETSC_FALSE,skipnorm=PETSC_FALSE;
78: Vec IV[2];
80: PetscFunctionBeginUser;
81: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
83: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
84: PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&mu,NULL));
85: PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
86: PetscCall(PetscOptionsGetReal(NULL,NULL,"-kappa",&kappa,NULL));
87: PetscCall(PetscOptionsGetBool(NULL,NULL,"-initv",&initv,NULL));
88: PetscCall(PetscOptionsGetBool(NULL,NULL,"-skipnorm",&skipnorm,NULL));
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: /* K is a tridiagonal */
95: PetscCall(MatCreate(PETSC_COMM_WORLD,&K));
96: PetscCall(MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n));
97: PetscCall(MatSetFromOptions(K));
99: PetscCall(MatGetOwnershipRange(K,&Istart,&Iend));
100: for (i=Istart;i<Iend;i++) {
101: if (i>0) PetscCall(MatSetValue(K,i,i-1,-kappa,INSERT_VALUES));
102: PetscCall(MatSetValue(K,i,i,kappa*3.0,INSERT_VALUES));
103: if (i<n-1) PetscCall(MatSetValue(K,i,i+1,-kappa,INSERT_VALUES));
104: }
106: PetscCall(MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY));
107: PetscCall(MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY));
109: /* C is a tridiagonal */
110: PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
111: PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n));
112: PetscCall(MatSetFromOptions(C));
114: PetscCall(MatGetOwnershipRange(C,&Istart,&Iend));
115: for (i=Istart;i<Iend;i++) {
116: if (i>0) PetscCall(MatSetValue(C,i,i-1,-tau,INSERT_VALUES));
117: PetscCall(MatSetValue(C,i,i,tau*3.0,INSERT_VALUES));
118: if (i<n-1) PetscCall(MatSetValue(C,i,i+1,-tau,INSERT_VALUES));
119: }
121: PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY));
122: PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY));
124: /* M is a diagonal matrix */
125: PetscCall(MatCreate(PETSC_COMM_WORLD,&M));
126: PetscCall(MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n));
127: PetscCall(MatSetFromOptions(M));
128: PetscCall(MatGetOwnershipRange(M,&Istart,&Iend));
129: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(M,i,i,mu,INSERT_VALUES));
130: PetscCall(MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY));
131: PetscCall(MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY));
133: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134: Create the eigensolver and set various options
135: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
137: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
138: A[0] = K; A[1] = C; A[2] = M;
139: PetscCall(PEPSetOperators(pep,3,A));
140: PetscCall(PEPSetProblemType(pep,PEP_GENERAL));
141: PetscCall(PEPSetTolerances(pep,PETSC_SMALL,PETSC_CURRENT));
142: if (initv) { /* initial vector */
143: PetscCall(MatCreateVecs(K,&IV[0],NULL));
144: PetscCall(VecSetValue(IV[0],0,-1.0,INSERT_VALUES));
145: PetscCall(VecSetValue(IV[0],1,0.5,INSERT_VALUES));
146: PetscCall(VecAssemblyBegin(IV[0]));
147: PetscCall(VecAssemblyEnd(IV[0]));
148: PetscCall(MatCreateVecs(K,&IV[1],NULL));
149: PetscCall(VecSetValue(IV[1],0,4.0,INSERT_VALUES));
150: PetscCall(VecSetValue(IV[1],2,1.5,INSERT_VALUES));
151: PetscCall(VecAssemblyBegin(IV[1]));
152: PetscCall(VecAssemblyEnd(IV[1]));
153: PetscCall(PEPSetInitialSpace(pep,2,IV));
154: PetscCall(VecDestroy(&IV[0]));
155: PetscCall(VecDestroy(&IV[1]));
156: }
157: PetscCall(PEPSetFromOptions(pep));
159: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
160: Solve the eigensystem
161: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163: PetscCall(PEPSolve(pep));
164: PetscCall(PEPGetDimensions(pep,&nev,NULL,NULL));
165: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
167: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168: Display solution and clean up
169: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
171: PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
172: if (!skipnorm) PetscCall(CheckNormalizedVectors(pep));
173: PetscCall(PEPDestroy(&pep));
174: PetscCall(MatDestroy(&M));
175: PetscCall(MatDestroy(&C));
176: PetscCall(MatDestroy(&K));
177: PetscCall(SlepcFinalize());
178: return 0;
179: }
181: /*TEST
183: testset:
184: args: -pep_nev 4 -initv
185: requires: !single
186: output_file: output/test2_1.out
187: test:
188: suffix: 1
189: args: -pep_type {{toar linear}}
190: test:
191: suffix: 1_toar_mgs
192: args: -pep_type toar -bv_orthog_type mgs
193: test:
194: suffix: 1_qarnoldi
195: args: -pep_type qarnoldi -bv_orthog_refine never
196: test:
197: suffix: 1_linear_gd
198: args: -pep_type linear -pep_linear_eps_type gd -pep_linear_explicitmatrix
200: testset:
201: args: -pep_target -0.43 -pep_nev 4 -pep_ncv 20 -st_type sinvert
202: output_file: output/test2_2.out
203: test:
204: suffix: 2
205: args: -pep_type {{toar linear}}
206: test:
207: suffix: 2_toar_scaleboth
208: args: -pep_type toar -pep_scale both
209: test:
210: suffix: 2_toar_transform
211: args: -pep_type toar -st_transform
212: test:
213: suffix: 2_qarnoldi
214: args: -pep_type qarnoldi -bv_orthog_refine always
215: test:
216: suffix: 2_linear_explicit
217: args: -pep_type linear -pep_linear_explicitmatrix -pep_linear_linearization 0,1
218: test:
219: suffix: 2_linear_explicit_her
220: args: -pep_type linear -pep_linear_explicitmatrix -pep_hermitian -pep_linear_linearization {{0,1 1,0 .3,.7}}
221: test:
222: suffix: 2_stoar
223: args: -pep_type stoar -pep_hermitian
224: requires: !single
225: test:
226: suffix: 2_jd
227: args: -pep_type jd -st_type precond -pep_max_it 200 -pep_ncv 24
228: requires: !single
230: test:
231: suffix: 3
232: args: -pep_nev 12 -pep_extract {{none norm residual structured}} -pep_monitor_cancel
233: requires: !single
235: testset:
236: args: -st_type sinvert -pep_target -0.43 -pep_nev 4
237: output_file: output/test2_2.out
238: test:
239: suffix: 4_schur
240: args: -pep_refine simple -pep_refine_scheme schur
241: test:
242: suffix: 4_mbe
243: args: -pep_refine simple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
244: test:
245: suffix: 4_explicit
246: args: -pep_refine simple -pep_refine_scheme explicit
247: test:
248: suffix: 4_multiple_schur
249: args: -pep_refine multiple -pep_refine_scheme schur
250: requires: !single
251: test:
252: suffix: 4_multiple_mbe
253: args: -pep_refine multiple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu -pep_refine_pc_factor_shift_type nonzero
254: test:
255: suffix: 4_multiple_explicit
256: args: -pep_refine multiple -pep_refine_scheme explicit
257: requires: !single
259: test:
260: suffix: 5
261: nsize: 2
262: args: -pep_type linear -pep_linear_explicitmatrix -pep_general -pep_target -0.43 -pep_nev 4 -pep_ncv 20 -st_type sinvert -pep_linear_st_ksp_type bcgs -pep_linear_st_pc_type bjacobi
263: output_file: output/test2_2.out
265: test:
266: suffix: 6
267: args: -pep_type linear -pep_nev 12 -pep_extract {{none norm residual}}
268: requires: !single
269: output_file: output/test2_3.out
271: test:
272: suffix: 7
273: args: -pep_nev 12 -pep_extract {{none norm residual structured}} -pep_refine multiple
274: requires: !single
275: output_file: output/test2_3.out
277: testset:
278: args: -st_type sinvert -pep_target -0.43 -pep_nev 4 -st_transform
279: output_file: output/test2_2.out
280: test:
281: suffix: 8_schur
282: args: -pep_refine simple -pep_refine_scheme schur
283: test:
284: suffix: 8_mbe
285: args: -pep_refine simple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
286: test:
287: suffix: 8_explicit
288: args: -pep_refine simple -pep_refine_scheme explicit
289: test:
290: suffix: 8_multiple_schur
291: args: -pep_refine multiple -pep_refine_scheme schur
292: test:
293: suffix: 8_multiple_mbe
294: args: -pep_refine multiple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
295: test:
296: suffix: 8_multiple_explicit
297: args: -pep_refine multiple -pep_refine_scheme explicit
299: testset:
300: nsize: 2
301: args: -st_type sinvert -pep_target -0.49 -pep_nev 4 -pep_refine_partitions 2 -st_ksp_type bcgs -st_pc_type bjacobi -pep_scale diagonal -pep_scale_its 4
302: output_file: output/test2_2.out
303: test:
304: suffix: 9_mbe
305: args: -pep_refine simple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
306: test:
307: suffix: 9_explicit
308: args: -pep_refine simple -pep_refine_scheme explicit
309: test:
310: suffix: 9_multiple_mbe
311: args: -pep_refine multiple -pep_refine_scheme mbe -pep_refine_ksp_type preonly -pep_refine_pc_type lu
312: requires: !single
313: test:
314: suffix: 9_multiple_explicit
315: args: -pep_refine multiple -pep_refine_scheme explicit
316: requires: !single
318: test:
319: suffix: 10
320: nsize: 4
321: args: -st_type sinvert -pep_target -0.43 -pep_nev 4 -pep_refine simple -pep_refine_scheme explicit -pep_refine_partitions 2 -st_ksp_type bcgs -st_pc_type bjacobi -pep_scale diagonal -pep_scale_its 4
322: output_file: output/test2_2.out
324: testset:
325: args: -pep_nev 4 -initv -mat_type aijcusparse
326: output_file: output/test2_1.out
327: requires: cuda !single
328: test:
329: suffix: 11_cuda
330: args: -pep_type {{toar linear}}
331: test:
332: suffix: 11_cuda_qarnoldi
333: args: -pep_type qarnoldi -bv_orthog_refine never
334: test:
335: suffix: 11_cuda_linear_gd
336: args: -pep_type linear -pep_linear_eps_type gd -pep_linear_explicitmatrix
338: test:
339: suffix: 12
340: nsize: 2
341: args: -pep_type jd -ds_parallel synchronized -pep_target -0.43 -pep_nev 4 -pep_ncv 20
342: requires: !single
344: test:
345: suffix: 13
346: args: -pep_nev 12 -pep_view_values draw -pep_monitor draw::draw_lg
347: requires: x !single
348: output_file: output/test2_3.out
350: test:
351: suffix: 14
352: requires: complex double
353: args: -pep_type ciss -rg_type ellipse -rg_ellipse_center -48.5 -rg_ellipse_radius 1.5 -pep_ciss_delta 1e-10
355: testset:
356: args: -pep_nev 4 -initv -mat_type aijhipsparse
357: output_file: output/test2_1.out
358: requires: hip !single
359: test:
360: suffix: 15_hip
361: args: -pep_type {{toar linear}}
362: test:
363: suffix: 15_hip_qarnoldi
364: args: -pep_type qarnoldi -bv_orthog_refine never
365: test:
366: suffix: 15_hip_linear_gd
367: args: -pep_type linear -pep_linear_eps_type gd -pep_linear_explicitmatrix
369: TEST*/