Actual source code: acoustic_wave_1d.c
slepc-3.5.2 2014-10-10
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2014, 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];
53: SlepcInitialize(&argc,&argv,(char*)0,help);
55: PetscOptionsGetInt(NULL,"-n",&n,NULL);
56: PetscOptionsGetScalar(NULL,"-z",&z,NULL);
57: SlepcSNPrintfScalar(str,50,z,PETSC_FALSE);
58: PetscPrintf(PETSC_COMM_WORLD,"\nAcoustic wave 1-D, n=%D z=%s\n\n",n,str);
60: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
62: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64: /* K is a tridiagonal */
65: MatCreate(PETSC_COMM_WORLD,&K);
66: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
67: MatSetFromOptions(K);
68: MatSetUp(K);
69:
70: MatGetOwnershipRange(K,&Istart,&Iend);
71: for (i=Istart;i<Iend;i++) {
72: if (i>0) {
73: MatSetValue(K,i,i-1,-1.0*n,INSERT_VALUES);
74: }
75: if (i<n-1) {
76: MatSetValue(K,i,i,2.0*n,INSERT_VALUES);
77: MatSetValue(K,i,i+1,-1.0*n,INSERT_VALUES);
78: } else {
79: MatSetValue(K,i,i,1.0*n,INSERT_VALUES);
80: }
81: }
83: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
84: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
86: /* C is the zero matrix but one element*/
87: MatCreate(PETSC_COMM_WORLD,&C);
88: MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,n,n);
89: MatSetFromOptions(C);
90: MatSetUp(C);
92: MatGetOwnershipRange(C,&Istart,&Iend);
93: if (n-1>=Istart && n-1<Iend) {
94: MatSetValue(C,n-1,n-1,-2*PETSC_PI/z,INSERT_VALUES);
95: }
96: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
97: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
98:
99: /* M is a diagonal matrix */
100: MatCreate(PETSC_COMM_WORLD,&M);
101: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
102: MatSetFromOptions(M);
103: MatSetUp(M);
105: MatGetOwnershipRange(M,&Istart,&Iend);
106: for (i=Istart;i<Iend;i++) {
107: if (i<n-1) {
108: MatSetValue(M,i,i,4*PETSC_PI*PETSC_PI/n,INSERT_VALUES);
109: } else {
110: MatSetValue(M,i,i,2*PETSC_PI*PETSC_PI/n,INSERT_VALUES);
111: }
112: }
113: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
114: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
115:
116: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117: Create the eigensolver and solve the problem
118: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120: PEPCreate(PETSC_COMM_WORLD,&pep);
121: A[0] = K; A[1] = C; A[2] = M;
122: PEPSetOperators(pep,3,A);
123: PEPSetFromOptions(pep);
124: PEPSolve(pep);
126: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127: Display solution and clean up
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129:
130: PEPPrintSolution(pep,NULL);
131: PEPDestroy(&pep);
132: MatDestroy(&M);
133: MatDestroy(&C);
134: MatDestroy(&K);
135: SlepcFinalize();
136: return 0;
137: }