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.229093264702743, 22.636447493332746, 38.75912680796266, 56.95688372815061, 64.68904675602266, 71.6151064215854, 76.89311293757974, 86.36950975000735, 90.92730267719605, 99.90890224816397, 108.07709801799932, 123.53078334588308]
1 : [3.2310126754296573, 6.0862342637162685, 6.406899224048421, 12.054533474443556, 12.672421029209339, 13.138599688893457, 14.490032245974174, 17.72368170716995, 18.66258440032833, 19.346120604313242, 23.999160954696467, 25.451852240424653]
2 : [3.220379348186461, 5.886623822610067, 5.905548429287503, 10.85374739080105, 10.90128797008324, 11.033338300615899, 12.490608609931913, 13.713193968097578, 14.184813087394302, 15.34867454583331, 16.52982807401637, 17.96454034806697]
3 : [3.220253450726732, 5.8807687086776275, 5.8819084880274355, 10.709133420835983, 10.726052065979124, 10.778656219121304, 12.345960396979809, 12.760676043108143, 13.760281355578835, 14.01832604286101, 14.607718515235726, 16.592214148827882]
4 : [3.2202514669005944, 5.880494279981992, 5.880556917690724, 10.689894953485808, 10.698983267279571, 10.718281105553432, 12.323621310872117, 12.51513457506993, 13.53091563611022, 13.646352490555426, 14.013638403691404, 15.844042100616626]
5 : [3.220251432780628, 5.880477114537656, 5.880480023263922, 10.686851036556126, 10.694700891024654, 10.700890407057575, 12.31910979366758, 12.416240501257302, 13.409477889900058, 13.55256788716293, 13.749216310962613, 15.24645246597867]
6 : [3.2202514321849876, 5.880473890511289, 5.880477843318865, 10.68616100873073, 10.694032876757733, 10.695820538590594, 12.318002400653944, 12.365108987045693, 13.375722155055369, 13.488497017394572, 13.61398182455544, 14.813746363843293]
7 : [3.220251432174622, 5.8804736764845, 5.880477751399416, 10.685969287196409, 10.693915072917022, 10.69441539784443, 12.31769621979804, 12.33861020790337, 13.365112608786578, 13.45294954022866, 13.531662728703642, 14.548118345194396]
8 : [3.2202514321744453, 5.880473664772787, 5.880477746015275, 10.68591773176339, 10.693892108137469, 10.694043457650626, 12.317604665085351, 12.326098068775448, 13.360878384964767, 13.435960447007234, 13.480566199148813, 14.39803652491641]
9 : [3.220251432174437, 5.880473664135112, 5.880477745705933, 10.685904872913927, 10.693887275313632, 10.693949433862919, 12.317575537840142, 12.320797130323182, 13.358877986653807, 13.428258104919683, 13.451525412357942, 14.315735088991557]
10 : [3.2202514321744413, 5.8804736641003235, 5.880477745688503, 10.68590179219335, 10.693886175010086, 10.693926453071475, 12.317565160998454, 12.31871575381345, 13.35789296628391, 13.424793977655028, 13.436298694491471, 14.270543537516643]
11 : [3.2202514321744427, 5.880473664098429, 5.880477745687527, 10.685901066924764, 10.69388591658361, 10.693920952758663, 12.317558831672562, 12.31793983359151, 13.35741456916091, 13.423228036492146, 13.428747504480109, 14.24562647101346]
12 : [3.220251432174441, 5.88047366409832, 5.880477745687475, 10.685900897346334, 10.693885856844805, 10.693919648675335, 12.31754556287273, 12.317668127173524, 13.357188433110492, 13.422515876005525, 13.425120095594798, 14.231830193829234]
13 : [3.220251432174439, 5.880473664098283, 5.880477745687478, 10.68590085778872, 10.693885843152536, 10.693919341040505, 12.317512863018015, 12.317598526260495, 13.357083769221198, 13.42218947303075, 13.423410961488505, 14.224178775107049]
14 : [3.2202514321744373, 5.880473664098298, 5.880477745687484, 10.68590084856391, 10.693885840013829, 10.6939192686345, 12.317488502233305, 12.317586424503578, 13.357035967179314, 13.422036765448446, 13.42261629070797, 14.219930340448318]
15 : [3.2202514321744413, 5.8804736640982815, 5.880477745687459, 10.68590084641141, 10.693885839291767, 10.693919251610547, 12.317478504242507, 12.317583476121808, 13.357014278666956, 13.421959842857056, 13.422254132838905, 14.217570323842114]
16 : [3.2202514321744404, 5.88047366409831, 5.8804777456875685, 10.685900845908808, 10.693885839125384, 10.693919247609113, 12.317474830676995, 12.317582560125347, 13.357004471079998, 13.421913956679962, 13.422096929670005, 14.216258604475598]
17 : [3.220251432174441, 5.880473664098474, 5.880477745687609, 10.685900845791057, 10.693885839086484, 10.693919246667255, 12.317473511217447, 12.317582249335695, 13.357000038837713, 13.421883678985806, 13.422032681703506, 14.215525265904184]
18 : [3.220251432174442, 5.880473664098161, 5.880477745687394, 10.6859008457628, 10.693885839077915, 10.69391924644418, 12.31747304282984, 12.317582141683848, 13.35699780954583, 13.421865251083528, 13.422007317745356, 14.215104833891592]
19 : [3.220251432174441, 5.880473664098301, 5.880477745687783, 10.685900845757008, 10.693885839075739, 10.693919246393346, 12.317472877655769, 12.317582103575255, 13.356996964469651, 13.421856122060275, 13.421996511917209, 14.214878697682334]
[6]:
print ("Eigenvalues")
for lam in evals:
print (lam)
Eigenvalues
3.220251432174441
5.880473664098301
5.880477745687783
10.685900845757008
10.693885839075739
10.693919246393346
12.317472877655769
12.317582103575255
13.356996964469651
13.421856122060275
13.421996511917209
14.214878697682334
[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);
[ ]:
[ ]: