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.229093264910017, 22.636447501310272, 38.759126808515624, 56.95688372587462, 64.68904675401627, 71.61510643643915, 76.89311294243652, 86.36950975535069, 90.92730267554337, 99.90890224263467, 108.0770980116467, 123.5307833446752]
1 : [3.2310126754353625, 6.086234263773771, 6.406899223991599, 12.054533474398818, 12.67242102912911, 13.138599686784541, 14.49003224902289, 17.723681706335324, 18.66258439910937, 19.346120603919253, 23.999160959479802, 25.451852240433542]
2 : [3.2203793481865373, 5.8866238226109, 5.905548429292155, 10.853747390745843, 10.90128797020628, 11.033338300496151, 12.490608609998057, 13.713193967558734, 14.184813087466393, 15.348674545985382, 16.52982807563945, 17.96454034923895]
3 : [3.2202534507267324, 5.8807687086777545, 5.881908488027902, 10.709133420839786, 10.726052065983636, 10.778656219111069, 12.34596039697965, 12.760676042992067, 13.760281355689752, 14.018326042967432, 14.607718515637963, 16.592214150599002]
4 : [3.220251466900598, 5.880494279982013, 5.880556917690762, 10.689894953487492, 10.698983267280372, 10.71828110555352, 12.323621310869129, 12.515134575045614, 13.53091563615175, 13.646352490639108, 14.013638403924357, 15.844042102921314]
5 : [3.2202514327806324, 5.880477114537639, 5.880480023263911, 10.686851036557021, 10.69470089102496, 10.700890407059042, 12.319109793665936, 12.416240501278963, 13.40947788990967, 13.552567887227255, 13.749216311201069, 15.246452468152999]
6 : [3.22025143218498, 5.880473890511284, 5.880477843318866, 10.686161008731158, 10.694032876757856, 10.69582053859151, 12.318002400653254, 12.365108987080287, 13.37572215505905, 13.488497017424878, 13.613981824812424, 14.813746365401355]
7 : [3.220251432174624, 5.880473676484475, 5.880477751399407, 10.685969287196565, 10.693915072917038, 10.694415397844809, 12.317696219797753, 12.338610207929385, 13.365112608788833, 13.452949540241514, 13.531662728933414, 14.54811834615401]
8 : [3.220251432174445, 5.880473664772773, 5.880477746015266, 10.68591773176345, 10.69389210813747, 10.69404345765075, 12.317604665085229, 12.32609806879331, 13.3608783849663, 13.43596044702098, 13.480566199380982, 14.398036525545237]
9 : [3.220251432174438, 5.88047366413512, 5.880477745705945, 10.6859048729153, 10.693887275313681, 10.693949433866047, 12.31757553783788, 12.320797131470359, 13.358877986627162, 13.42825810403486, 13.451525427543544, 14.315735123590821]
10 : [3.220251432174443, 5.880473664100332, 5.880477745688498, 10.68590179219415, 10.693886175010217, 10.693926453073352, 12.31756516099855, 12.318715754506481, 13.357892966554576, 13.424793974979444, 13.436298707022425, 14.270543581043746]
11 : [3.2202514321744418, 5.8804736640984085, 5.8804777456875135, 10.685901066924986, 10.693885916583834, 10.693920952758612, 12.3175588316784, 12.317939833650616, 13.357414569398246, 13.423228034323651, 13.428747507676745, 14.245626494281687]
12 : [3.2202514321744404, 5.8804736640983135, 5.880477745687485, 10.6859008973463, 10.693885856844705, 10.69391964867532, 12.317545562865725, 12.317668127061854, 13.357188433372249, 13.422515875154179, 13.42512008738087, 14.231830185469805]
13 : [3.220251432174438, 5.8804736640983055, 5.880477745687485, 10.685900857788846, 10.693885843152731, 10.69391934103933, 12.317512862167568, 12.317598525517825, 13.357083769364031, 13.422189473178223, 13.423410934293502, 14.224178689637185]
14 : [3.2202514321744395, 5.880473664098301, 5.880477745687493, 10.685900848563827, 10.693885840013678, 10.693919268633321, 12.317488501522073, 12.317586424285075, 13.35703596703179, 13.422036763589485, 13.422616262904016, 14.219930196234545]
15 : [3.22025143217444, 5.880473664098318, 5.880477745687521, 10.685900846411192, 10.693885839291793, 10.693919251608182, 12.317478500591161, 12.317583475116834, 13.35701427695884, 13.4219598346521, 13.422253950524366, 14.217569169664038]
16 : [3.220251432174444, 5.880473664098355, 5.880477745687343, 10.685900845908353, 10.693885839124846, 10.693919247601793, 12.317474816761738, 12.31758255617534, 13.357004470738248, 13.421913862341702, 13.422096586791417, 14.216254837425035]
17 : [3.2202514321744418, 5.880473664098205, 5.880477745687441, 10.685900845791004, 10.693885839086139, 10.693919246665095, 12.317473503135414, 12.317582247009359, 13.35700004012301, 13.421883543586336, 13.422032054887083, 14.215520897393487]
18 : [3.2202514321744444, 5.880473664098386, 5.880477745687485, 10.685900845763667, 10.693885839077748, 10.693919246445633, 12.317473040743927, 12.317582140519214, 13.356998038430767, 13.421865751200825, 13.42200682223182, 14.215115956322808]
19 : [3.2202514321744427, 5.880473664098223, 5.880477745687425, 10.685900845757397, 10.693885839075868, 10.693919246394234, 12.317472878734899, 12.317582103525643, 13.356997136491158, 13.421856638786501, 13.421996495703791, 14.214893645338543]
[6]:
print ("Eigenvalues")
for lam in evals:
print (lam)
Eigenvalues
3.2202514321744427
5.880473664098223
5.880477745687425
10.685900845757397
10.693885839075868
10.693919246394234
12.317472878734899
12.317582103525643
13.356997136491158
13.421856638786501
13.421996495703791
14.214893645338543
[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);
[ ]:
[ ]: