Actual source code: acoustic_wave_1d.c
slepc-3.6.1 2015-09-03
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
8: SLEPc is free software: you can redistribute it and/or modify it under the
9: terms of version 3 of the GNU Lesser General Public License as published by
10: the Free Software Foundation.
12: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
13: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
14: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
15: more details.
17: You should have received a copy of the GNU Lesser General Public License
18: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20: */
21: /*
22: This example implements one of the problems found at
23: NLEVP: A Collection of Nonlinear Eigenvalue Problems,
24: The University of Manchester.
25: The details of the collection can be found at:
26: [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
27: Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.
29: The acoustic_wave_1d problem is a QEP from an acoustics application.
30: Here we solve it with the eigenvalue scaled by the imaginary unit, to be
31: able to use real arithmetic, so the computed eigenvalues should be scaled
32: back.
33: */
35: static char help[] = "NLEVP problem: acoustic_wave_1d.\n\n"
36: "The command line options are:\n"
37: " -n <n>, where <n> = dimension of the matrices.\n"
38: " -z <z>, where <z> = impedance (default 1.0).\n\n";
40: #include <slepcpep.h>
44: int main(int argc,char **argv)
45: {
46: Mat M,C,K,A[3]; /* problem matrices */
47: PEP pep; /* polynomial eigenproblem solver context */
48: PetscInt n=10,Istart,Iend,i;
49: PetscScalar z=1.0;
50: char str[50];
51: PetscBool terse;
54: SlepcInitialize(&argc,&argv,(char*)0,help);
56: PetscOptionsGetInt(NULL,"-n",&n,NULL);
57: PetscOptionsGetScalar(NULL,"-z",&z,NULL);
58: SlepcSNPrintfScalar(str,50,z,PETSC_FALSE);
59: PetscPrintf(PETSC_COMM_WORLD,"\nAcoustic wave 1-D, n=%D z=%s\n\n",n,str);
61: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
63: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
65: /* K is a tridiagonal */
66: MatCreate(PETSC_COMM_WORLD,&K);
67: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
68: MatSetFromOptions(K);
69: MatSetUp(K);
70:
71: MatGetOwnershipRange(K,&Istart,&Iend);
72: for (i=Istart;i<Iend;i++) {
73: if (i>0) {
74: MatSetValue(K,i,i-1,-1.0*n,INSERT_VALUES);
75: }
76: if (i<n-1) {
77: MatSetValue(K,i,i,2.0*n,INSERT_VALUES);
78: MatSetValue(K,i,i+1,-1.0*n,INSERT_VALUES);
79: } else {
80: MatSetValue(K,i,i,1.0*n,INSERT_VALUES);
81: }
82: }
84: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
85: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
87: /* C is the zero matrix but one element*/
88: MatCreate(PETSC_COMM_WORLD,&C);
89: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
90: MatSetFromOptions(C);
91: MatSetUp(C);
93: MatGetOwnershipRange(C,&Istart,&Iend);
94: if (n-1>=Istart && n-1<Iend) {
95: MatSetValue(C,n-1,n-1,-2*PETSC_PI/z,INSERT_VALUES);
96: }
97: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
98: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
99:
100: /* M is a diagonal matrix */
101: MatCreate(PETSC_COMM_WORLD,&M);
102: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
103: MatSetFromOptions(M);
104: MatSetUp(M);
106: MatGetOwnershipRange(M,&Istart,&Iend);
107: for (i=Istart;i<Iend;i++) {
108: if (i<n-1) {
109: MatSetValue(M,i,i,4*PETSC_PI*PETSC_PI/n,INSERT_VALUES);
110: } else {
111: MatSetValue(M,i,i,2*PETSC_PI*PETSC_PI/n,INSERT_VALUES);
112: }
113: }
114: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
115: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
116:
117: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118: Create the eigensolver and solve the problem
119: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
121: PEPCreate(PETSC_COMM_WORLD,&pep);
122: A[0] = K; A[1] = C; A[2] = M;
123: PEPSetOperators(pep,3,A);
124: PEPSetFromOptions(pep);
125: PEPSolve(pep);
127: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: Display solution and clean up
129: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130:
131: /* show detailed info unless -terse option is given by user */
132: PetscOptionsHasName(NULL,"-terse",&terse);
133: if (terse) {
134: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
135: } else {
136: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
137: PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
138: PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
139: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
140: }
141: PEPDestroy(&pep);
142: MatDestroy(&M);
143: MatDestroy(&C);
144: MatDestroy(&K);
145: SlepcFinalize();
146: return 0;
147: }