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 = 43192
[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 : [8.21312981438697, 17.79981324288422, 37.219825654290396, 44.01930032500969, 62.18233337557714, 65.53333262068308, 73.98716307852285, 81.08231276265101, 88.06001749544146, 91.25055697995714, 113.61980588223211, 122.21673389966705]
1 : [3.238097539155719, 5.9986506305549785, 6.161390071943701, 11.485860551459016, 12.037309681895579, 14.413520502500129, 15.234814500779297, 16.218269370508033, 17.34672695201068, 19.89228392636556, 22.638813209689285, 24.80756093365987]
2 : [3.220542688575668, 5.884876708899267, 5.88822609551715, 10.785471632148738, 10.847156783408213, 11.381986377803349, 13.343238318155038, 13.37270650615316, 14.112393285300925, 14.988639563390914, 15.344080015364577, 17.688329565510738]
3 : [3.2202954145601517, 5.880713446554104, 5.880839526590821, 10.707466526318392, 10.72006917841996, 10.944773310604969, 12.527973954549319, 13.075036549488201, 13.540697915816448, 14.030695748094939, 14.367494218076756, 15.906741569390611]
4 : [3.2202898941021947, 5.880502253290996, 5.880518097321735, 10.69511436308039, 10.699129958356712, 10.787751811345048, 12.35409837066285, 12.871053036661346, 13.445683487771273, 13.649403154219073, 14.213859958921322, 14.607665437391077]
5 : [3.2202897664212125, 5.8804933752546695, 5.8804959631163545, 10.692776367606898, 10.694999509695847, 10.721906317148337, 12.32420257887225, 12.640658943020448, 13.417157954178274, 13.509762214867187, 13.816759271693694, 14.260863114299163]
6 : [3.220289763807026, 5.880492978117859, 5.880494510444945, 10.691097181798142, 10.694161153653802, 10.698643594307601, 12.318883281763688, 12.460746603564504, 13.402394959766507, 13.454795875732096, 13.563293261472136, 14.232709515946997]
7 : [3.2202897637595163, 5.880492958816486, 5.88049442503795, 10.68786162766492, 10.693985793779925, 10.694552792288166, 12.317899571105146, 12.372364425930282, 13.387450842550567, 13.43393941279623, 13.473271157464499, 14.223431358496684]
8 : [3.2202897637586987, 5.880492957834426, 5.880494420312152, 10.68646436809804, 10.693942498545146, 10.694085266047452, 12.317711055329475, 12.337392167381001, 13.37373542749864, 13.426640065230014, 13.441852510393783, 14.219267826273315]
9 : [3.220289763758687, 5.8804929577828275, 5.880494420056971, 10.686108698958849, 10.693925371312476, 10.69400489399861, 12.317673740213694, 12.324680164772042, 13.365194188406006, 13.423921247790018, 13.430320729520341, 14.217191023865118]
10 : [3.220289763758686, 5.880492957780057, 5.880494420043334, 10.686025216520271, 10.693919077435641, 10.693989757202965, 12.317666113412333, 12.320210206487413, 13.360871916858914, 13.422816823149185, 13.425650479552912, 14.21610359485821]
11 : [3.220289763758686, 5.880492957779928, 5.88049442004259, 10.68600601101789, 10.693917414576262, 10.693986570955634, 12.31766449253041, 12.31865920108561, 13.35884929149572, 13.422344037503395, 13.423640048290112, 14.215520284287884]
12 : [3.2202897637586845, 5.8804929577799205, 5.880494420042561, 10.686001612135907, 10.693917021614125, 10.69398585920919, 12.317664119925215, 12.318123649773188, 13.357930686392498, 13.422136539821203, 13.422748122569207, 14.215203524469747]
13 : [3.2202897637586876, 5.8804929577799, 5.880494420042568, 10.686000605526424, 10.693916931073238, 10.69398569752759, 12.317664014900839, 12.31793901225933, 13.357518459926846, 13.422044245774648, 13.422347218840969, 14.215030369576104]
14 : [3.2202897637586867, 5.880492957779877, 5.880494420042552, 10.686000374980662, 10.693916910324093, 10.693985660611734, 12.317663977189008, 12.317875397756389, 13.357334494861478, 13.42200299802874, 13.422165823276192, 14.21493538911811]
15 : [3.220289763758687, 5.880492957780066, 5.880494420043104, 10.686000322143956, 10.693916905560336, 10.693985652167362, 12.317663963000292, 12.317853444062056, 13.357252453180658, 13.421984371458896, 13.42208361024183, 14.214883183911088]
16 : [3.2202897637586854, 5.880492957779983, 5.880494420042415, 10.68600031000016, 10.693916904464935, 10.693985650228816, 12.317663957906365, 12.317845854551532, 13.357215843271762, 13.42197592445536, 13.422046240568335, 14.214854416686027]
17 : [3.2202897637586845, 5.880492957780348, 5.880494420042721, 10.686000307203287, 10.693916904212651, 10.69398564978306, 12.317663956130755, 12.317843221804242, 13.35719944364148, 13.421972073984133, 13.422029220074242, 14.2148384818532]
18 : [3.220289763758689, 5.880492957779712, 5.880494420042257, 10.686000306537897, 10.693916904151333, 10.693985649677796, 12.31766395551302, 12.317842262046225, 13.357191794648443, 13.421970168096665, 13.422021005795505, 14.214829709476758]
19 : [3.220289763758688, 5.880492957779606, 5.880494420042267, 10.686000306405504, 10.693916904140213, 10.693985649655854, 12.317663955307653, 12.317841968695548, 13.357188615779025, 13.421969408342136, 13.422017592825718, 14.214824860494982]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.220289763758688
5.880492957779606
5.880494420042267
10.686000306405504
10.693916904140213
10.693985649655854
12.317663955307653
12.317841968695548
13.357188615779025
13.421969408342136
13.422017592825718
14.214824860494982
[7]:
gfu = GridFunction(fes, multidim=len(evecs))
for i in range(len(evecs)):
    gfu.vecs[i].data = evecs[i]
Draw (Norm(gfu), mesh, order=2, min=0, max=2);
[ ]:

[ ]: