Eigenvalue problems

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.811240988942085, 54.33837125704101, 59.25128716057947, 128.12681732979485, 140.0829925396113, 170.82416393357428, 199.30981821128336, 322.82217832358367, 356.6786710663352, 461.5903135039409, 517.3870889635418, 625.7205006219798]
1 : [19.739301603250873, 49.38806350193087, 49.399648518732555, 80.2694729622326, 99.91354180553904, 101.16717245414851, 130.94799280502548, 150.7375727812236, 180.3937325056378, 199.29706697891055, 243.3789051868941, 275.0069022158784]
2 : [19.739209736107778, 49.348179164989475, 49.34823924863315, 78.97556274229092, 98.71753152603921, 98.77541632451575, 128.5140452514619, 129.49669541839452, 168.96924037508117, 173.53824643623028, 184.96483075753085, 204.37620613282508]
3 : [19.73920970375886, 49.34805057616807, 49.34805663276679, 78.9571421694649, 98.69661976741001, 98.69738490103244, 128.31140649384932, 128.32786041915892, 167.83585991923263, 168.23479011439665, 177.96647880522957, 197.8093379364212]
4 : [19.73920970375293, 49.348050362042656, 49.34805634696643, 78.95704924286433, 98.69648901545278, 98.69653276278275, 128.30635539909957, 128.30659054351045, 167.78848149696717, 167.82142238404, 177.68669578321598, 197.43802071565725]
5 : [19.73920970375297, 49.348050361663425, 49.34805634657529, 78.9570487785702, 98.69648799390052, 98.69652648656172, 128.30613642309658, 128.30637471269, 167.78711949677762, 167.79291728485182, 177.66540362582137, 197.40227381534297]
6 : [19.739209703752934, 49.348050361663354, 49.34805634657509, 78.95704877665491, 98.69648798064448, 98.69652639235399, 128.30613002138608, 128.3063699744855, 167.7870393130533, 167.78762793655832, 177.65937008666378, 197.39974063564642]
7 : [19.739209703752937, 49.348050361663326, 49.34805634657493, 78.95704877664548, 98.69648798047194, 98.69652639140092, 128.30612985593953, 128.30636985368804, 167.78693710288493, 167.7871472780485, 177.65812002687525, 197.39944830395518]
8 : [19.739209703752934, 49.34805036166354, 49.34805634657498, 78.95704877664537, 98.69648798047017, 98.69652639139217, 128.3061298514664, 128.30636985120066, 167.78690759693745, 167.7871302503458, 177.65792227865532, 197.39942301864343]
9 : [19.739209703752934, 49.34805036166282, 49.348056346574594, 78.95704877664555, 98.69648798046997, 98.69652639139213, 128.30612985137873, 128.3063698511398, 167.78690569939013, 167.78712940546137, 177.65790629694, 197.3994216262459]
10 : [19.73920970375294, 49.34805036166286, 49.34805634657479, 78.9570487766446, 98.69648798046991, 98.6965263913923, 128.30612985137608, 128.3063698511377, 167.78690552874374, 167.7871293292202, 177.65790470851587, 197.3994214044584]
11 : [19.739209703752923, 49.34805036166304, 49.34805634657514, 78.95704877664542, 98.69648798046975, 98.6965263913925, 128.30612985137563, 128.3063698511371, 167.78690551906624, 167.78712932498513, 177.6579045804347, 197.3994213649057]
12 : [19.739209703752927, 49.34805036166239, 49.34805634657526, 78.95704877664538, 98.69648798047021, 98.6965263913919, 128.30612985137515, 128.30636985113756, 167.7869055183641, 167.78712932467013, 177.65790456946863, 197.39942133390767]
13 : [19.739209703752937, 49.348050361663056, 49.348056346574744, 78.95704877664544, 98.69648798046975, 98.69652639139213, 128.30612985137586, 128.30636985113767, 167.78690551829519, 167.787129324639, 177.65790456851948, 197.39942119187288]
14 : [19.73920970375293, 49.34805036166297, 49.34805634657543, 78.9570487766445, 98.6964879804703, 98.69652639139187, 128.3061298513756, 128.3063698511375, 167.7869055182923, 167.78712932463785, 177.6579045684699, 197.39941991884996]
15 : [19.739209703752927, 49.34805036166297, 49.348056346575184, 78.95704877664495, 98.69648798046961, 98.69652639139207, 128.30612985137583, 128.30636985113748, 167.78690551829175, 167.78712932463748, 177.65790456846625, 197.39940600408127]
16 : [19.739209703752927, 49.34805036166286, 49.348056346574445, 78.95704877664471, 98.69648798047037, 98.69652639139221, 128.30612985137589, 128.30636985113708, 167.7869055182917, 167.78712932463773, 177.65790456846622, 197.39924614263873]
17 : [19.73920970375294, 49.34805036166257, 49.348056346574445, 78.95704877664492, 98.69648798047035, 98.69652639139173, 128.30612985137572, 128.30636985113688, 167.78690551829246, 167.78712932463765, 177.65790456846676, 197.3990725580613]
18 : [19.73920970375293, 49.34805036166333, 49.34805634657417, 78.95704877664517, 98.69648798047004, 98.69652639139184, 128.3061298513754, 128.30636985113708, 167.78690551829192, 167.7871293246379, 177.65790456846608, 197.39904841452565]
19 : [19.739209703752934, 49.34805036166291, 49.34805634657483, 78.95704877664541, 98.69648798047021, 98.69652639139261, 128.30612985137554, 128.3063698511372, 167.78690551829195, 167.78712932463742, 177.65790456846625, 197.3990458189554]
2.000000091348567
5.00000287308538
5.00000347948374
8.000021638955511
10.000044983523006
10.000048875363127
13.000128945108859
13.000153262168782
17.000367866798435
17.000390543122265
18.00050917429449
20.000704972244627
gfu = GridFunction(fes, multidim=0)
for evec in evecs:
    gfu.AddMultiDimComponent(evec)

Draw (gfu, min=-2, max=2);