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.229093262376307, 22.636447503549284, 38.759126808803885, 56.95688372999922, 64.68904675073665, 71.61510642452876, 76.89311293866422, 86.36950975302643, 90.92730267567258, 99.90890224551634, 108.07709802279689, 123.53078334165765]
1 : [3.2310126754198696, 6.0862342638011295, 6.406899223986624, 12.05453347334034, 12.672421029383163, 13.138599687809759, 14.490032246971388, 17.72368170651508, 18.662584396240455, 19.346120605793114, 23.99916095892357, 25.451852239100393]
2 : [3.2203793481863547, 5.886623822609558, 5.905548429288102, 10.85374739074184, 10.901287970059592, 11.033338300583145, 12.490608609960082, 13.713193967654698, 14.184813087307006, 15.348674547006746, 16.529828077630487, 17.964540348018893]
3 : [3.220253450726725, 5.88076870867771, 5.8819084880279, 10.709133420840962, 10.726052065986194, 10.778656219117055, 12.345960396985742, 12.760676043254755, 13.760281355651495, 14.018326043969342, 14.607718515444361, 16.592214148636195]
4 : [3.220251466900592, 5.880494279982023, 5.8805569176907895, 10.689894953486398, 10.698983267283454, 10.718281105549773, 12.323621310873046, 12.515134575175228, 13.530915636500483, 13.64635249035674, 14.013638403796556, 15.844042100368293]
5 : [3.2202514327806315, 5.880477114537657, 5.880480023263918, 10.686851036556348, 10.69470089102519, 10.700890407056024, 12.319109793667765, 12.416240501308188, 13.409477890011722, 13.55256788702332, 13.749216311053587, 15.246452465742502]
6 : [3.220251432184984, 5.880473890511269, 5.880477843318862, 10.686161008730789, 10.694032876757806, 10.695820538590045, 12.318002400654079, 12.365108987067046, 13.375722155108082, 13.488497017304994, 13.613981824603512, 14.813746363651548]
7 : [3.220251432174622, 5.880473676484497, 5.880477751399423, 10.685969287196405, 10.69391507291702, 10.69441539784431, 12.317696219798085, 12.338610207911527, 13.365112608820217, 13.45294954017943, 13.531662728717356, 14.548118345054826]
8 : [3.2202514321744453, 5.880473664772771, 5.880477746015269, 10.685917731763421, 10.693892108137469, 10.694043457650611, 12.317604665085383, 12.326098068779558, 13.360878384986124, 13.435960446995658, 13.480566199243615, 14.398036524906468]
9 : [3.220251432174438, 5.880473664135115, 5.880477745705952, 10.6859048729148, 10.693887275314315, 10.693949433864695, 12.317575537839515, 12.320797130997056, 13.358877986652955, 13.428258104747398, 13.451525424425377, 14.315735109706619]
10 : [3.220251432174441, 5.880473664100284, 5.8804777456885, 10.685901792193926, 10.693886175010395, 10.69392645307258, 12.31756516099852, 12.318715754552224, 13.357892966659588, 13.424793975854536, 13.436298709578633, 14.270543577098868]
11 : [3.220251432174444, 5.8804736640984165, 5.880477745687536, 10.68590106692516, 10.69388591658371, 10.693920952759568, 12.31755883168256, 12.317939834104939, 13.35741456950004, 13.423228034879285, 13.428747516711448, 14.245626512741467]
12 : [3.220251432174443, 5.880473664098305, 5.8804777456874815, 10.685900897346338, 10.693885856844707, 10.693919648675887, 12.317545562905138, 12.317668127348785, 13.357188433607563, 13.422515876412755, 13.425120098481754, 14.23183020607282]
13 : [3.2202514321744404, 5.880473664098321, 5.880477745687494, 10.685900857787738, 10.693885843152724, 10.693919341029789, 12.317512858180072, 12.317598522411586, 13.357083768121779, 13.422189473014937, 13.423410785547661, 14.224177976797115]
14 : [3.2202514321744413, 5.880473664098306, 5.880477745687476, 10.68590084856366, 10.693885840013696, 10.693919268631019, 12.317488498901344, 12.31758642323177, 13.35703596571532, 13.422036752390145, 13.422616110413417, 14.219929718063439]
15 : [3.2202514321744418, 5.880473664098243, 5.88047774568741, 10.685900846411213, 10.693885839291639, 10.693919251609486, 12.317478502854087, 12.317583475735077, 13.35701427705286, 13.421959832599617, 13.422254052278541, 14.217569954466171]
16 : [3.220251432174443, 5.880473664098548, 5.880477745687563, 10.685900845908597, 10.693885839123736, 10.693919247601642, 12.317474809974131, 12.317582552892382, 13.35700445485342, 13.421913920608153, 13.4220963657161, 14.216255547683073]
17 : [3.2202514321744395, 5.880473664098458, 5.880477745687365, 10.685900845791313, 10.693885839086361, 10.693919246662515, 12.317473491399639, 12.317582242680416, 13.357000029242656, 13.421883568681329, 13.422032394590767, 14.215521370760248]
18 : [3.220251432174442, 5.880473664098351, 5.8804777456874255, 10.685900845763735, 10.69388583907761, 10.693919246444562, 12.317473032901658, 12.317582137948458, 13.356998038282118, 13.421865641034554, 13.422006646805317, 14.215112693656588]
19 : [3.2202514321744373, 5.88047366409837, 5.8804777456875, 10.685900845757553, 10.69388583907593, 10.693919246393902, 12.317472874550278, 12.317582102255193, 13.356997137112666, 13.4218564709155, 13.421996213436499, 14.214888741188242]
[6]:
print ("Eigenvalues")
for lam in evals:
print (lam)
Eigenvalues
3.2202514321744373
5.88047366409837
5.8804777456875
10.685900845757553
10.69388583907593
10.693919246393902
12.317472874550278
12.317582102255193
13.356997137112666
13.4218564709155
13.421996213436499
14.214888741188242
[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);
[ ]:
[ ]: