Eigenvalue problems

29. 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.796097336039125, 58.533583335792414, 63.491962149323065, 139.22793680443462, 165.03490076786736, 203.19369937104176, 299.6645516236699, 329.9618332109861, 356.25089482426995, 453.7279250935756, 493.3263167078267, 647.6039282736585]
1 : [19.739382255585443, 49.42226035466374, 49.47862028468172, 82.86369754652067, 103.25375531647632, 106.80442334016817, 144.2635523384763, 153.9425900192458, 187.78710846556461, 206.2443988408144, 267.48597833259953, 273.439848948533]
2 : [19.739209836751634, 49.348341666754365, 49.348381247627636, 79.00911877776073, 98.78252061968159, 98.90545072429421, 128.75030179710203, 131.76428226269545, 170.5425837436024, 173.97438798164143, 186.4129018571641, 207.35433351741102]
3 : [19.73920975377069, 49.34805382893231, 49.348059764343425, 78.9574605072455, 98.69749585807088, 98.6984941226605, 128.31571256801416, 128.40753094537672, 167.98632689680142, 168.43932591412707, 178.23363601224906, 198.2867221246065]
4 : [19.73920975374541, 49.34805312850327, 49.34805880174678, 78.95706884504324, 98.69653353689009, 98.6965892107736, 128.3066684302896, 128.3087076202016, 167.7934663916605, 167.8430222862805, 177.69426484606896, 197.47468039631661]
5 : [19.73920975374541, 49.34805312716073, 49.348058800217835, 78.9570662426703, 98.6965246912183, 98.69657087071222, 128.30627128437354, 128.30647571080738, 167.78734565849805, 167.79225804876344, 177.6622002754929, 197.40378555054465]
6 : [19.739209753745385, 49.34805312715688, 49.348058800216414, 78.95706622606022, 98.69652463289026, 98.69657059515876, 128.3062479904469, 128.30643644851403, 167.78712023256648, 167.78805130561358, 177.65859150034083, 197.40037064411828]
7 : [19.739209753745396, 49.34805312715686, 49.34805880021673, 78.95706622595054, 98.6965246322767, 98.69657059190757, 128.30624744050655, 128.30643556744246, 167.7871019363846, 167.7876406746849, 177.65822193799661, 197.4002144768528]
8 : [19.739209753745385, 49.34805312715699, 49.34805880021662, 78.95706622595026, 98.69652463227212, 98.69657059186018, 128.30624742860337, 128.3064355403203, 167.7871007126377, 167.7876134766808, 177.65819956574015, 197.40020467501878]
9 : [19.7392097537454, 49.348053127156994, 49.34805880021657, 78.95706622594992, 98.69652463227204, 98.69657059185927, 128.30624742817358, 128.306435539818, 167.78710066170345, 167.7876128616985, 177.65819876923527, 197.40020350589555]
10 : [19.73920975374539, 49.348053127157094, 49.348058800216194, 78.95706622594987, 98.69652463227186, 98.69657059185953, 128.30624742816315, 128.30643553980676, 167.78710065903613, 167.78761284306384, 177.65819873430388, 197.40020335762122]
11 : [19.739209753745392, 49.34805312715713, 49.348058800216116, 78.95706622595006, 98.69652463227195, 98.69657059185967, 128.30624742816252, 128.30643553980622, 167.7871006588647, 167.78761284210182, 177.6581987286562, 197.40020331710215]
12 : [19.739209753745403, 49.348053127157094, 49.348058800216485, 78.95706622595024, 98.69652463227237, 98.69657059185947, 128.30624742816244, 128.3064355398064, 167.78710065885238, 167.78761284196526, 177.6581987260406, 197.40020324183666]
13 : [19.73920975374539, 49.348053127156845, 49.348058800215696, 78.9570662259503, 98.69652463227227, 98.69657059185968, 128.3062474281623, 128.30643553980545, 167.78710065885136, 167.7876128419476, 177.65819872550517, 197.40020215034963]
14 : [19.739209753745406, 49.348053127156604, 49.34805880021617, 78.95706622595023, 98.69652463227169, 98.6965705918601, 128.30624742816195, 128.30643553980582, 167.7871006588514, 167.78761284194644, 177.65819872547854, 197.40019060953847]
15 : [19.739209753745403, 49.34805312715662, 49.348058800215824, 78.9570662259501, 98.69652463227173, 98.69657059185951, 128.30624742816246, 128.30643553980616, 167.78710065885122, 167.7876128419465, 177.65819872547678, 197.400082269594]
16 : [19.7392097537454, 49.34805312715733, 49.34805880021614, 78.95706622594987, 98.69652463227216, 98.69657059185958, 128.30624742816244, 128.30643553980613, 167.78710065885144, 167.7876128419462, 177.65819872547672, 197.39949982070624]
17
 : [19.739209753745406, 49.348053127157165, 49.34805880021631, 78.95706622594994, 98.69652463227207, 98.69657059185958, 128.30624742816215, 128.30643553980636, 167.7871006588507, 167.78761284194655, 177.65819872547695, 197.3993928187826]
18
 : [19.73920975374539, 49.3480531271568, 49.348058800216016, 78.95706622595033, 98.69652463227175, 98.69657059185963, 128.3062474281621, 128.30643553980633, 167.7871006588514, 167.78761284194644, 177.6581987254767, 197.3993825012752]
19
 : [19.739209753745396, 49.34805312715746, 49.34805880021549, 78.9570662259502, 98.69652463227204, 98.69657059185982, 128.30624742816212, 128.3064355398062, 167.78710065885147, 167.78761284194647, 177.65819872547738, 197.3993816084001]
2.0000000964138622
5.000003153288563
5.000003728089516
8.000023406939727
10.00004869712695
10.000053353806784
13.000140858128043
13.000159917822478
17.000387638670905
17.00043953366833
18.00053897863103
20.00073899482862
gfu = GridFunction(fes, multidim=0)
for evec in evecs:
    gfu.AddMultiDimComponent(evec)

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