24. 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.81124098894189, 54.33837125702939, 59.251287160581754, 128.12681732983054, 140.082992539563, 170.82416393349706, 199.3098182114319, 322.82217832368667, 356.67867106663, 461.5903135038212, 517.3870889633408, 625.720500621145]
1 : [19.739301603250876, 49.38806350192768, 49.399648518730075, 80.2694729622346, 99.9135418055376, 101.16717245415786, 130.9479928050399, 150.7375727811963, 180.39373250563062, 199.29706697899493, 243.37890518670434, 275.0069022158316]
2 : [19.739209736107792, 49.348179164989475, 49.348239248633305, 78.97556274229115, 98.71753152603894, 98.77541632451648, 128.51404525146347, 129.49669541839515, 168.9692403750798, 173.53824643624918, 184.96483075753622, 204.37620613281814]
3 : [19.73920970375884, 49.348050576168056, 49.34805663276713, 78.95714216946529, 98.69661976740994, 98.69738490103255, 128.31140649384972, 128.32786041915895, 167.83585991923258, 168.23479011439838, 177.96647880522897, 197.809337936424]
4 : [19.73920970375293, 49.348050362042535, 49.34805634696672, 78.95704924286461, 98.69648901545305, 98.69653276278318, 128.30635539909952, 128.30659054351074, 167.788481496967, 167.82142238404, 177.68669578321638, 197.4380207156559]
5 : [19.73920970375296, 49.34805036166317, 49.34805634657541, 78.95704877857021, 98.6964879939008, 98.69652648656154, 128.30613642309666, 128.30637471269003, 167.7871194967777, 167.79291728483514, 177.6654036258026, 197.40227381537346]
6 : [19.739209703752937, 49.34805036166336, 49.34805634657455, 78.95704877665513, 98.69648798064448, 98.69652639235382, 128.30613002138617, 128.30636997448585, 167.78703931305026, 167.78762793666692, 177.65937008679148, 197.39974063584378]
7 : [19.739209703752945, 49.3480503616632, 49.34805634657484, 78.9570487766454, 98.6964879804718, 98.69652639140085, 128.30612985593953, 128.30636985368795, 167.78693710277517, 167.78714727793906, 177.6581200276498, 197.39944830515714]
8 : [19.73920970375294, 49.34805036166223, 49.348056346574396, 78.95704877664528, 98.69648798047037, 98.69652639139214, 128.3061298514663, 128.3063698512005, 167.7869075913768, 167.78713024637753, 177.65792227057244, 197.39942301251472]
9 : [19.739209703752927, 49.34805036166259, 49.348056346575234, 78.95704877664548, 98.69648798047007, 98.69652639139208, 128.30612985137904, 128.30636985114066, 167.78690569838363, 167.78712940474762, 177.6579062790221, 197.3994216174728]
10 : [19.73920970375293, 49.3480503616631, 49.34805634657517, 78.9570487766451, 98.69648798047011, 98.69652639139198, 128.30612985137583, 128.3063698511374, 167.78690552972964, 167.78712932968446, 177.6579047088812, 197.39942140912845]
11 : [19.739209703752927, 49.34805036166317, 49.34805634657462, 78.95704877664504, 98.6964879804699, 98.69652639139204, 128.30612985137597, 128.30636985113733, 167.78690551926223, 167.7871293250783, 177.65790458110382, 197.3994213735342]
12 : [19.739209703752937, 49.34805036166288, 49.34805634657482, 78.95704877664485, 98.69648798047015, 98.69652639139207, 128.30612985137572, 128.3063698511375, 167.78690551838676, 167.78712932468125, 177.65790456966144, 197.3994213417151]
13 : [19.739209703752923, 49.348050361663184, 49.34805634657505, 78.95704877664527, 98.69648798047018, 98.69652639139214, 128.30612985137577, 128.30636985113762, 167.7869055182979, 167.78712932464043, 177.65790456854643, 197.39942120701645]
14 : [19.739209703752937, 49.34805036166306, 49.34805634657515, 78.95704877664498, 98.69648798046997, 98.69652639139181, 128.30612985137563, 128.30636985113702, 167.786905518292, 167.7871293246376, 177.65790456847202, 197.3994203184489]
15 : [19.73920970375293, 49.34805036166232, 49.34805634657487, 78.95704877664502, 98.69648798046951, 98.6965263913916, 128.30612985137566, 128.30636985113748, 167.78690551829158, 167.7871293246376, 177.65790456846656, 197.39941443715892]
16 : [19.739209703752913, 49.348050361663375, 49.34805634657478, 78.95704877664507, 98.69648798047004, 98.6965263913922, 128.30612985137537, 128.30636985113725, 167.78690551829212, 167.7871293246379, 177.65790456846594, 197.3993649105956]
17 : [19.739209703752927, 49.3480503616631, 49.34805634657452, 78.95704877664518, 98.69648798046984, 98.69652639139223, 128.30612985137606, 128.30636985113762, 167.78690551829203, 167.78712932463762, 177.65790456846642, 197.39912773111075]
18 : [19.739209703752923, 49.34805036166307, 49.34805634657454, 78.95704877664583, 98.69648798046963, 98.69652639139196, 128.30612985137583, 128.30636985113745, 167.7869055182917, 167.78712932463773, 177.65790456846642, 197.39905497579815]
19
: [19.739209703752937, 49.34805036166266, 49.34805634657472, 78.95704877664502, 98.69648798047035, 98.69652639139267, 128.3061298513756, 128.30636985113748, 167.78690551829158, 167.7871293246378, 177.65790456846608, 197.3990467596755]
2.0000000913485674
5.000002873085355
5.00000347948373
8.000021638955472
10.00004498352302
10.000048875363133
13.000128945108864
13.000153262168812
17.0003678667984
17.000390543122304
18.00050917429447
20.000705067559505
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