27. Eigenvalue problems#
The Helmholtz equation (with hard boundary conditions) is singular when \(\omega^2\) is an eigenvalue of the negative Laplace operator, i.e.
\[
-\Delta u = \lambda u
\]
The weak form is to find eigen pairs \((\lambda, u) \in {\mathbf R} \times H^1\) such that
\[
\int \nabla u \nabla v \, dx = \lambda \int u v \, dx
\]
Finite element discretization leads to the generalized matrix eigenvalue problem
\[
A u = \lambda M u.
\]
from ngsolve import *
from ngsolve.webgui import Draw
mesh = Mesh(unit_square.GenerateMesh(maxh=0.1))
fes = H1(mesh, order=3, dirichlet=".*")
u,v = fes.TnT()
a = BilinearForm(grad(u)*grad(v)*dx).Assemble()
m = BilinearForm(u*v*dx).Assemble()
Eigenvectorsare stationary points of the Rayleigh quotient
\[
\frac{\left< A u, u \right> }{ \left< M u, u \right> }
\]
This characterization is used in the LOBPCG (Locally Optimal Block Preconditioned Conjugate Gradient) method.
import ngsolve.solvers
pre = a.mat.Inverse(freedofs=fes.FreeDofs())
evals, evecs = solvers.LOBPCG(a.mat, m.mat, pre=pre, num=12, maxit=20, printrates=True)
for lam in evals:
print (lam/pi**2)
0 : [19.811240988944178, 54.33837125704087, 59.25128716058129, 128.126817329797, 140.0829925396108, 170.8241639335785, 199.30981821129515, 322.82217832360783, 356.6786710663663, 461.5903135039543, 517.3870889635471, 625.7205006219388]
1 : [19.73930160325087, 49.388063501926396, 49.39964851873146, 80.2694729622347, 99.91354180553937, 101.16717245415118, 130.94799280502545, 150.73757278123063, 180.39373250563935, 199.2970669789112, 243.37890518688477, 275.00690221587183]
2 : [19.73920973610779, 49.348179164990086, 49.34823924863233, 78.97556274229181, 98.71753152603834, 98.77541632451624, 128.51404525146216, 129.4966954183955, 168.96924037508094, 173.53824643623233, 184.96483075753028, 204.37620613282297]
3 : [19.73920970375886, 49.34805057616805, 49.34805663276696, 78.95714216946521, 98.69661976741, 98.69738490103214, 128.31140649384918, 128.3278604191592, 167.83585991923235, 168.23479011439653, 177.96647880523034, 197.8093379364245]
4 : [19.739209703752927, 49.348050362042486, 49.348056346966736, 78.95704924286451, 98.69648901545274, 98.69653276278264, 128.30635539909937, 128.30659054351094, 167.78848149696714, 167.8214223840411, 177.686695783217, 197.43802071565932]
5 : [19.739209703752966, 49.34805036166353, 49.34805634657536, 78.95704877857025, 98.69648799390056, 98.69652648656147, 128.30613642309655, 128.3063747126899, 167.78711949677765, 167.79291728482244, 177.66540362580363, 197.402273815319]
6 : [19.739209703752934, 49.34805036166286, 49.348056346574936, 78.95704877665536, 98.69648798064468, 98.69652639235356, 128.30613002138546, 128.30636997448536, 167.78703931304904, 167.78762793668196, 177.659370086895, 197.39974063603145]
7 : [19.739209703752955, 49.34805036166279, 49.348056346574495, 78.95704877664524, 98.69648798047166, 98.69652639140105, 128.30612985593956, 128.30636985368795, 167.78693710372016, 167.7871472786963, 177.65812004000338, 197.3994483030294]
8 : [19.73920970375294, 49.348050361662246, 49.3480563465745, 78.9570487766451, 98.69648798047017, 98.69652639139206, 128.3061298514662, 128.30636985120046, 167.7869075994817, 167.78713025168935, 177.6579222797308, 197.39942301502822]
9 : [19.73920970375293, 49.34805036166305, 49.34805634657479, 78.95704877664524, 98.69648798047012, 98.69652639139179, 128.30612985137847, 128.30636985113986, 167.78690569341936, 167.7871294027892, 177.65790625612752, 197.39942162596853]
10 : [19.739209703752934, 49.34805036166273, 49.34805634657464, 78.95704877664474, 98.69648798047015, 98.69652639139207, 128.30612985137566, 128.30636985113748, 167.78690552979154, 167.78712932967878, 177.65790471773275, 197.3994214103718]
11 : [19.73920970375293, 49.34805036166355, 49.348056346575156, 78.95704877664498, 98.69648798046995, 98.69652639139238, 128.30612985137572, 128.3063698511373, 167.78690551928403, 167.78712932508785, 177.65790458142067, 197.3994213730823]
12 : [19.739209703752923, 49.348050361663454, 49.348056346574644, 78.95704877664505, 98.69648798046971, 98.6965263913915, 128.3061298513753, 128.30636985113713, 167.78690551839375, 167.78712932468534, 177.6579045694638, 197.39942134379686]
13 : [19.73920970375294, 49.34805036166347, 49.348056346575326, 78.95704877664538, 98.69648798046943, 98.69652639139201, 128.30612985137526, 128.306369851137, 167.78690551829774, 167.7871293246405, 177.657904568533, 197.39942120117325]
14 : [19.73920970375293, 49.34805036166325, 49.3480563465749, 78.95704877664552, 98.69648798047004, 98.69652639139237, 128.3061298513759, 128.3063698511371, 167.78690551829223, 167.787129324638, 177.65790456847213, 197.39942046944293]
15 : [19.739209703752937, 49.34805036166296, 49.34805634657475, 78.95704877664559, 98.69648798047014, 98.69652639139193, 128.30612985137589, 128.30636985113736, 167.7869055182922, 167.7871293246374, 177.6579045684667, 197.3994144761106]
16 : [19.739209703752934, 49.34805036166277, 49.34805634657462, 78.95704877664525, 98.69648798046975, 98.69652639139164, 128.30612985137506, 128.3063698511375, 167.7869055182923, 167.78712932463776, 177.65790456846656, 197.39936865280853]
17 : [19.739209703752934, 49.348050361662814, 49.34805634657466, 78.9570487766451, 98.69648798047, 98.69652639139218, 128.30612985137574, 128.30636985113773, 167.78690551829195, 167.7871293246372, 177.65790456846634, 197.39913322818526]
18 : [19.739209703752937, 49.34805036166338, 49.34805634657509, 78.9570487766453, 98.69648798047001, 98.69652639139186, 128.30612985137608, 128.30636985113767, 167.78690551829143, 167.78712932463807, 177.65790456846588, 197.399054358519]
19 : [19.739209703752934, 49.348050361663326, 49.34805634657516, 78.95704877664511, 98.6964879804695, 98.69652639139171, 128.30612985137608, 128.30636985113702, 167.78690551829192, 167.78712932463787, 177.65790456846682, 197.39904646877545]
2.000000091348567
5.000002873085423
5.000003479483774
8.00002163895548
10.000044983522933
10.000048875363037
13.000128945108914
13.000153262168766
17.000367866798435
17.00039054312231
18.000509174294546
20.000705038085165
gfu = GridFunction(fes, multidim=0)
for evec in evecs:
gfu.AddMultiDimComponent(evec)
Draw (gfu, min=-2, max=2);
compute eigenfunctions on more complex domains.
compute eigenfunctions for the elasticity equation (e.g. a vibrating beam)
Eigenvalues of a tuning fork: https://www.youtube.com/watch?v=35Yh8kGbMLA