Actual source code: acoustic_wave_2d.c
slepc-3.7.3 2016-09-29
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2016, 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_2d problem is a 2-D version of acoustic_wave_1d, also
30: scaled for real arithmetic.
31: */
33: static char help[] = "Quadratic eigenproblem from an acoustics application (2-D).\n\n"
34: "The command line options are:\n"
35: " -m <m>, where <m> = grid size, the matrices have dimension m*(m-1).\n"
36: " -z <z>, where <z> = impedance (default 1.0).\n\n";
38: #include <slepcpep.h>
42: int main(int argc,char **argv)
43: {
44: Mat M,C,K,A[3]; /* problem matrices */
45: PEP pep; /* polynomial eigenproblem solver context */
46: PetscInt m=6,n,II,Istart,Iend,i,j;
47: PetscScalar z=1.0;
48: PetscReal h;
49: char str[50];
50: PetscBool terse;
53: SlepcInitialize(&argc,&argv,(char*)0,help);
55: PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
56: if (m<2) SETERRQ(PETSC_COMM_SELF,1,"m must be at least 2");
57: PetscOptionsGetScalar(NULL,NULL,"-z",&z,NULL);
58: h = 1.0/m;
59: n = m*(m-1);
60: SlepcSNPrintfScalar(str,50,z,PETSC_FALSE);
61: PetscPrintf(PETSC_COMM_WORLD,"\nAcoustic wave 2-D, n=%D (m=%D), z=%s\n\n",n,m,str);
63: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64: Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
65: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
67: /* K has a pattern similar to the 2D Laplacian */
68: MatCreate(PETSC_COMM_WORLD,&K);
69: MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,n,n);
70: MatSetFromOptions(K);
71: MatSetUp(K);
72:
73: MatGetOwnershipRange(K,&Istart,&Iend);
74: for (II=Istart;II<Iend;II++) {
75: i = II/m; j = II-i*m;
76: if (i>0) { MatSetValue(K,II,II-m,(j==m-1)?-0.5:-1.0,INSERT_VALUES); }
77: if (i<m-2) { MatSetValue(K,II,II+m,(j==m-1)?-0.5:-1.0,INSERT_VALUES); }
78: if (j>0) { MatSetValue(K,II,II-1,-1.0,INSERT_VALUES); }
79: if (j<m-1) { MatSetValue(K,II,II+1,-1.0,INSERT_VALUES); }
80: MatSetValue(K,II,II,(j==m-1)?2.0:4.0,INSERT_VALUES);
81: }
83: MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
84: MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);
86: /* C is the zero matrix except for a few nonzero elements on the diagonal */
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: for (i=Istart;i<Iend;i++) {
94: if (i%m==m-1) {
95: MatSetValue(C,i,i,-2*PETSC_PI*h/z,INSERT_VALUES);
96: }
97: }
98: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
99: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
100:
101: /* M is a diagonal matrix */
102: MatCreate(PETSC_COMM_WORLD,&M);
103: MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,n,n);
104: MatSetFromOptions(M);
105: MatSetUp(M);
107: MatGetOwnershipRange(M,&Istart,&Iend);
108: for (i=Istart;i<Iend;i++) {
109: if (i%m==m-1) {
110: MatSetValue(M,i,i,2*PETSC_PI*PETSC_PI*h*h,INSERT_VALUES);
111: } else {
112: MatSetValue(M,i,i,4*PETSC_PI*PETSC_PI*h*h,INSERT_VALUES);
113: }
114: }
115: MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
116: MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
117:
118: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119: Create the eigensolver and solve the problem
120: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122: PEPCreate(PETSC_COMM_WORLD,&pep);
123: A[0] = K; A[1] = C; A[2] = M;
124: PEPSetOperators(pep,3,A);
125: PEPSetFromOptions(pep);
126: PEPSolve(pep);
128: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129: Display solution and clean up
130: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131:
132: /* show detailed info unless -terse option is given by user */
133: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
134: if (terse) {
135: PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL);
136: } else {
137: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
138: PEPReasonView(pep,PETSC_VIEWER_STDOUT_WORLD);
139: PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
140: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
141: }
142: PEPDestroy(&pep);
143: MatDestroy(&M);
144: MatDestroy(&C);
145: MatDestroy(&K);
146: SlepcFinalize();
147: return ierr;
148: }