Actual source code: acoustic_wave_1d.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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;
41: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
43: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
44: PetscOptionsGetScalar(NULL,NULL,"-z",&z,NULL);
45: SlepcSNPrintfScalar(str,50,z,PETSC_FALSE);
46: PetscPrintf(PETSC_COMM_WORLD,"\nAcoustic wave 1-D, n=%D z=%s\n\n",n,str);
48: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
50: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
52: /* K is a tridiagonal */
53: MatCreate(PETSC_COMM_WORLD,&K);
54: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
55: MatSetFromOptions(K);
56: MatSetUp(K);
58: MatGetOwnershipRange(K,&Istart,&Iend);
59: for (i=Istart;i<Iend;i++) {
60: if (i>0) {
61: MatSetValue(K,i,i-1,-1.0*n,INSERT_VALUES);
62: }
63: if (i<n-1) {
64: MatSetValue(K,i,i,2.0*n,INSERT_VALUES);
65: MatSetValue(K,i,i+1,-1.0*n,INSERT_VALUES);
66: } else {
67: MatSetValue(K,i,i,1.0*n,INSERT_VALUES);
68: }
69: }
71: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
72: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
74: /* C is the zero matrix but one element*/
75: MatCreate(PETSC_COMM_WORLD,&C);
76: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
77: MatSetFromOptions(C);
78: MatSetUp(C);
80: MatGetOwnershipRange(C,&Istart,&Iend);
81: if (n-1>=Istart && n-1<Iend) {
82: MatSetValue(C,n-1,n-1,-2*PETSC_PI/z,INSERT_VALUES);
83: }
84: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
85: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
87: /* M is a diagonal matrix */
88: MatCreate(PETSC_COMM_WORLD,&M);
89: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
90: MatSetFromOptions(M);
91: MatSetUp(M);
93: MatGetOwnershipRange(M,&Istart,&Iend);
94: for (i=Istart;i<Iend;i++) {
95: if (i<n-1) {
96: MatSetValue(M,i,i,4*PETSC_PI*PETSC_PI/n,INSERT_VALUES);
97: } else {
98: MatSetValue(M,i,i,2*PETSC_PI*PETSC_PI/n,INSERT_VALUES);
99: }
100: }
101: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
102: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
104: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105: Create the eigensolver and solve the problem
106: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108: PEPCreate(PETSC_COMM_WORLD,&pep);
109: A[0] = K; A[1] = C; A[2] = M;
110: PEPSetOperators(pep,3,A);
111: PEPSetFromOptions(pep);
112: PEPSolve(pep);
114: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115: Display solution and clean up
116: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118: /* show detailed info unless -terse option is given by user */
119: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
120: if (terse) {
121: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
122: } else {
123: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
124: PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
125: PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
126: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
127: }
128: PEPDestroy(&pep);
129: MatDestroy(&M);
130: MatDestroy(&C);
131: MatDestroy(&K);
132: SlepcFinalize();
133: return ierr;
134: }
136: /*TEST
138: testset:
139: args: -pep_nev 4 -pep_tol 1e-7 -n 24 -terse
140: output_file: output/acoustic_wave_1d_1.out
141: requires: !complex !single
142: test:
143: suffix: 1
144: args: -st_type sinvert -st_transform -pep_type {{toar qarnoldi linear}}
145: test:
146: suffix: 1_stoar
147: args: -st_type sinvert -st_transform -pep_type stoar -pep_hermitian -pep_stoar_locking 0
148: test:
149: suffix: 2
150: args: -st_type sinvert -st_transform -pep_type toar -pep_extract {{none norm residual}}
151: test:
152: suffix: 3
153: args: -st_type sinvert -pep_type linear -pep_extract {{none norm residual}}
154: test:
155: suffix: 4
156: args: -pep_type jd
157: TEST*/