Actual source code: acoustic_wave_1d.c
slepc-3.17.2 2022-08-09
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: This example implements one of the problems found at
12: NLEVP: A Collection of Nonlinear Eigenvalue Problems,
13: The University of Manchester.
14: The details of the collection can be found at:
15: [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
16: Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.
18: The acoustic_wave_1d problem is a QEP from an acoustics application.
19: Here we solve it with the eigenvalue scaled by the imaginary unit, to be
20: able to use real arithmetic, so the computed eigenvalues should be scaled
21: back.
22: */
24: static char help[] = "Quadratic eigenproblem from an acoustics application (1-D).\n\n"
25: "The command line options are:\n"
26: " -n <n>, where <n> = dimension of the matrices.\n"
27: " -z <z>, where <z> = impedance (default 1.0).\n\n";
29: #include <slepcpep.h>
31: int main(int argc,char **argv)
32: {
33: Mat M,C,K,A[3]; /* problem matrices */
34: PEP pep; /* polynomial eigenproblem solver context */
35: PetscInt n=10,Istart,Iend,i;
36: PetscScalar z=1.0;
37: char str[50];
38: PetscBool terse;
40: SlepcInitialize(&argc,&argv,(char*)0,help);
42: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
43: PetscOptionsGetScalar(NULL,NULL,"-z",&z,NULL);
44: SlepcSNPrintfScalar(str,sizeof(str),z,PETSC_FALSE);
45: PetscPrintf(PETSC_COMM_WORLD,"\nAcoustic wave 1-D, n=%" PetscInt_FMT " z=%s\n\n",n,str);
47: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
49: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
51: /* K is a tridiagonal */
52: MatCreate(PETSC_COMM_WORLD,&K);
53: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
54: MatSetFromOptions(K);
55: MatSetUp(K);
57: MatGetOwnershipRange(K,&Istart,&Iend);
58: for (i=Istart;i<Iend;i++) {
59: if (i>0) MatSetValue(K,i,i-1,-1.0*n,INSERT_VALUES);
60: if (i<n-1) {
61: MatSetValue(K,i,i,2.0*n,INSERT_VALUES);
62: MatSetValue(K,i,i+1,-1.0*n,INSERT_VALUES);
63: } else MatSetValue(K,i,i,1.0*n,INSERT_VALUES);
64: }
66: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
67: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
69: /* C is the zero matrix but one element*/
70: MatCreate(PETSC_COMM_WORLD,&C);
71: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
72: MatSetFromOptions(C);
73: MatSetUp(C);
75: MatGetOwnershipRange(C,&Istart,&Iend);
76: if (n-1>=Istart && n-1<Iend) MatSetValue(C,n-1,n-1,-2*PETSC_PI/z,INSERT_VALUES);
77: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
78: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
80: /* M is a diagonal matrix */
81: MatCreate(PETSC_COMM_WORLD,&M);
82: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
83: MatSetFromOptions(M);
84: MatSetUp(M);
86: MatGetOwnershipRange(M,&Istart,&Iend);
87: for (i=Istart;i<Iend;i++) {
88: if (i<n-1) MatSetValue(M,i,i,4*PETSC_PI*PETSC_PI/n,INSERT_VALUES);
89: else MatSetValue(M,i,i,2*PETSC_PI*PETSC_PI/n,INSERT_VALUES);
90: }
91: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
92: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Create the eigensolver and solve the problem
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: PEPCreate(PETSC_COMM_WORLD,&pep);
99: A[0] = K; A[1] = C; A[2] = M;
100: PEPSetOperators(pep,3,A);
101: PEPSetFromOptions(pep);
102: PEPSolve(pep);
104: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: Display solution and clean up
106: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108: /* show detailed info unless -terse option is given by user */
109: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
110: if (terse) PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
111: else {
112: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
113: PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
114: PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
115: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
116: }
117: PEPDestroy(&pep);
118: MatDestroy(&M);
119: MatDestroy(&C);
120: MatDestroy(&K);
121: SlepcFinalize();
122: return 0;
123: }
125: /*TEST
127: testset:
128: args: -pep_nev 4 -pep_tol 1e-7 -n 24 -terse
129: output_file: output/acoustic_wave_1d_1.out
130: requires: !single
131: test:
132: suffix: 1
133: args: -st_type sinvert -st_transform -pep_type {{toar qarnoldi linear}}
134: test:
135: suffix: 1_stoar
136: args: -st_type sinvert -st_transform -pep_type stoar -pep_hermitian -pep_stoar_locking 0 -pep_stoar_nev 11 -pep_ncv 10
137: test:
138: suffix: 2
139: args: -st_type sinvert -st_transform -pep_type toar -pep_extract {{none norm residual}}
140: test:
141: suffix: 3
142: args: -st_type sinvert -pep_type linear -pep_extract {{none norm residual}}
143: test:
144: suffix: 4
145: args: -pep_type jd
147: TEST*/