This page was generated from unit-2.4-Maxwell/Maxwellevp.ipynb.
2.4.1 Maxwell eigenvalue problem¶
We solve the Maxwell eigenvalue problem
\[\int \operatorname{curl} u \, \operatorname{curl} v
= \lambda \int u v\]
for \(u, v \; \bot \; \nabla H^1\) using a PINVIT solver from the ngsolve solvers module.
The orthogonality to gradient fields is important to eliminate the huge number of zero eigenvalues. The orthogonal sub-space is implemented using a Poisson projection:
\[P u = u - \nabla \Delta^{-1} \operatorname{div} u\]
The algorithm and example is take form the Phd thesis of Sabine Zaglmayr, p 145-150.
[1]:
from ngsolve import *
from ngsolve.webgui import Draw
from netgen.occ import *
[2]:
from netgen.occ import *
cube1 = Box( (-1,-1,-1), (1,1,1) )
cube2 = Box( (0,0,0), (2,2,2) )
cube2.edges.hpref=1 # mark edges for geometric refinement
fichera = cube1-cube2
Draw (fichera);
[3]:
mesh = Mesh(OCCGeometry(fichera).GenerateMesh(maxh=0.4))
mesh.RefineHP(levels=2, factor=0.2)
Draw (mesh);
[4]:
# SetHeapSize(100*1000*1000)
fes = HCurl(mesh, order=3)
print ("ndof =", fes.ndof)
u,v = fes.TnT()
a = BilinearForm(curl(u)*curl(v)*dx)
m = BilinearForm(u*v*dx)
apre = BilinearForm(curl(u)*curl(v)*dx + u*v*dx)
pre = Preconditioner(apre, "direct", inverse="sparsecholesky")
ndof = 42562
[5]:
with TaskManager():
a.Assemble()
m.Assemble()
apre.Assemble()
# build gradient matrix as sparse matrix (and corresponding scalar FESpace)
gradmat, fesh1 = fes.CreateGradient()
gradmattrans = gradmat.CreateTranspose() # transpose sparse matrix
math1 = gradmattrans @ m.mat @ gradmat # multiply matrices
math1[0,0] += 1 # fix the 1-dim kernel
invh1 = math1.Inverse(inverse="sparsecholesky")
# build the Poisson projector with operator Algebra:
proj = IdentityMatrix() - gradmat @ invh1 @ gradmattrans @ m.mat
projpre = proj @ pre.mat
evals, evecs = solvers.PINVIT(a.mat, m.mat, pre=projpre, num=12, maxit=20)
0 : [6.229093262820175, 22.63644749025295, 38.759126810977094, 56.956883733946135, 64.68904675899891, 71.61510644516521, 76.8931129213988, 86.36950976036896, 90.92730267643832, 99.90890224286441, 108.07709801664532, 123.5307833431483]
1 : [3.2310126754198487, 6.086234263856218, 6.406899224032846, 12.054533475641554, 12.672421029068799, 13.138599689521422, 14.490032245605784, 17.723681706289682, 18.662584394825764, 19.346120603728224, 23.99916096075235, 25.451852239054226]
2 : [3.2203793481863734, 5.886623822612676, 5.9055484292866, 10.853747390931511, 10.901287970101235, 11.03333830073588, 12.490608609845165, 13.713193966177878, 14.184813087343043, 15.348674545512887, 16.529828074222593, 17.964540349584034]
3 : [3.2202534507267297, 5.880768708677691, 5.881908488027309, 10.709133420835567, 10.726052065996036, 10.778656219151694, 12.345960396951515, 12.760676042579911, 13.760281355570093, 14.018326042190127, 14.607718515686285, 16.59221415218176]
4 : [3.2202514669005926, 5.880494279981996, 5.880556917690725, 10.68989495348676, 10.698983267281065, 10.718281105568547, 12.32362131086267, 12.51513457489618, 13.53091563582842, 13.64635249072207, 14.01363840395026, 15.844042105143586]
5 : [3.2202514327806324, 5.880477114537649, 5.8804800232639245, 10.68685103655696, 10.694700891024869, 10.700890407065597, 12.319109793664493, 12.416240501238002, 13.409477889799083, 13.552567887324328, 13.749216311207435, 15.246452470342732]
6 : [3.2202514321849836, 5.880473890511281, 5.88047784331886, 10.686161008731254, 10.694032876757808, 10.695820538594054, 12.318002400652963, 12.365108987082506, 13.375722155010333, 13.48849701749803, 13.613981824873118, 14.813746367028653]
7 : [3.2202514321746287, 5.880473676484508, 5.880477751399401, 10.6859692871966, 10.69391507291705, 10.694415397845678, 12.317696219797678, 12.338610207940281, 13.365112608764258, 13.45294954028565, 13.531662729010746, 14.548118347152705]
8 : [3.220251432174445, 5.8804736647727625, 5.8804777460152735, 10.685917731763448, 10.693892108137472, 10.694043457650977, 12.317604665085224, 12.326098068802114, 13.360878384953496, 13.435960447043813, 13.480566199467477, 14.398036526143148]
9 : [3.22025143217444, 5.880473664135107, 5.880477745705942, 10.685904872914465, 10.693887275313205, 10.693949433865097, 12.317575537839357, 12.320797131231341, 13.358877986643131, 13.428258104934052, 13.451525423009931, 14.315735109093776]
10 : [3.2202514321744435, 5.880473664100336, 5.880477745688503, 10.685901792193802, 10.693886175009986, 10.693926453072823, 12.317565160998619, 12.318715754371285, 13.357892966556575, 13.42479397575723, 13.436298702451683, 14.270543565593396]
11 : [3.220251432174442, 5.880473664098424, 5.880477745687534, 10.685901066925137, 10.693885916583495, 10.693920952760298, 12.31755883168381, 12.31793983414668, 13.35741456947277, 13.423228034794613, 13.428747513149716, 14.245626510428336]
12 : [3.220251432174441, 5.880473664098303, 5.880477745687484, 10.68590089734647, 10.693885856844588, 10.69391964867598, 12.317545562932844, 12.317668127502062, 13.357188433424898, 13.422515875542418, 13.425120104899394, 14.231830230973388]
13 : [3.2202514321744418, 5.880473664098319, 5.880477745687472, 10.685900857786802, 10.693885843152405, 10.693919341030396, 12.317512858384818, 12.317598522896327, 13.357083765940855, 13.422189471707993, 13.423410749052282, 14.224177813833341]
14 : [3.220251432174442, 5.8804736640983, 5.880477745687467, 10.685900848563241, 10.69388584001389, 10.693919268630554, 12.317488500406272, 12.31758642406262, 13.35703596428454, 13.422036767202846, 13.422616346397255, 14.219929933744343]
15 : [3.2202514321744413, 5.880473664098272, 5.88047774568743, 10.685900846411386, 10.693885839291926, 10.693919251609506, 12.317478503796883, 12.317583476125803, 13.357014276911487, 13.421959836489455, 13.42225423137519, 14.21757015538465]
16 : [3.2202514321744458, 5.880473664098191, 5.880477745687658, 10.685900845908813, 10.693885839125228, 10.693919247607438, 12.317474827481078, 12.317582559138748, 13.35700446843974, 13.421913945141908, 13.422096950346706, 14.21625822982107]
17 : [3.2202514321744435, 5.880473664098339, 5.880477745687491, 10.685900845790146, 10.693885839086436, 10.693919246654698, 12.317473464414137, 12.3175822369226, 13.357000009176513, 13.42188255653977, 13.422029680552118, 14.215491631766673]
18 : [3.2202514321744427, 5.880473664098336, 5.880477745687551, 10.685900845763218, 10.693885839077595, 10.693919246440744, 12.31747301499527, 12.317582133884411, 13.356998010575111, 13.421864630681075, 13.422004997876344, 14.21508518159039]
19 : [3.220251432174443, 5.880473664097934, 5.880477745687488, 10.685900845756787, 10.693885839075797, 10.693919246389822, 12.317472855628834, 12.317582098730286, 13.356996932138415, 13.421854446007396, 13.421995192230339, 14.2148274528146]
[6]:
print ("Eigenvalues")
for lam in evals:
print (lam)
Eigenvalues
3.220251432174443
5.880473664097934
5.880477745687488
10.685900845756787
10.693885839075797
10.693919246389822
12.317472855628834
12.317582098730286
13.356996932138415
13.421854446007396
13.421995192230339
14.2148274528146
[7]:
gfu = GridFunction(fes, multidim=len(evecs))
for i in range(len(evecs)):
gfu.vecs[i].data = evecs[i]
Draw (Norm(gfu), mesh, order=4, min=0, max=2);
[ ]:
[ ]: