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)
/usr/lib/python3/dist-packages/scipy/__init__.py:146: UserWarning: A NumPy version >=1.17.3 and <1.25.0 is required for this version of SciPy (detected version 1.26.4
  warnings.warn(f"A NumPy version >={np_minversion} and <{np_maxversion}"
0 : [8.213129812122666, 17.799813245531915, 37.219825664830545, 44.01930032912656, 62.182333370725836, 65.53333262600384, 73.98716308030491, 81.08231276008247, 88.06001748808589, 91.2505569707406, 113.61980588275604, 122.21673393053193]
1 : [3.2380975391571063, 5.998650630527236, 6.161390072147758, 11.485860551519867, 12.037309681264139, 14.41352050347481, 15.234814500453556, 16.21826937055013, 17.346726952703, 19.892283925625584, 22.638813216124397, 24.807560938691854]
2 : [3.2205426885757444, 5.884876708898261, 5.888226095522729, 10.785471632163086, 10.847156783298555, 11.381986377970442, 13.343238319030267, 13.37270650611903, 14.112393285428636, 14.988639563656623, 15.344080016682849, 17.68832956934841]
3 : [3.220295414560154, 5.880713446554281, 5.880839526590846, 10.707466526320449, 10.720069178417065, 10.944773310689657, 12.527973954780718, 13.075036549871008, 13.540697915757788, 14.030695748994336, 14.367494218023662, 15.90674157675796]
4 : [3.220289894102194, 5.880502253290988, 5.880518097321764, 10.695114363079439, 10.699129958359945, 10.787751811465204, 12.354098370707604, 12.871053037645979, 13.4456834877476, 13.649403154706432, 14.213859959482843, 14.607665443364922]
5 : [3.2202897664212085, 5.880493375254653, 5.8804959631163785, 10.692776367606877, 10.694999509696592, 10.721906317240387, 12.324202578881227, 12.640658944353879, 13.41715795418364, 13.509762215139625, 13.816759274356434, 14.260863114643659]
6 : [3.2202897638070245, 5.880492978117823, 5.880494510444942, 10.69109718180432, 10.694161153653965, 10.698643594339481, 12.318883281765658, 12.460746604406305, 13.402394959824125, 13.454795875882386, 13.563293262452582, 14.232709516019119]
7 : [3.220289763759515, 5.880492958816482, 5.880494425037949, 10.687861627672339, 10.693985793780039, 10.694552792292066, 12.31789957110563, 12.37236442629531, 13.387450842660245, 13.433939412862845, 13.473271157777504, 14.223431358525424]
8 : [3.220289763758694, 5.880492957834421, 5.880494420312165, 10.686464368099834, 10.693942498545209, 10.694085266048008, 12.317711055329603, 12.33739216751629, 13.373735427594351, 13.426640065258166, 13.441852510475414, 14.219267826290833]
9 : [3.220289763758685, 5.880492957782827, 5.880494420056993, 10.68610869895507, 10.69392537131241, 10.6940048939979, 12.317673740213772, 12.32468016463063, 13.36519418828596, 13.423921247810375, 13.430320729314223, 14.217191023837174]
10 : [3.220289763758685, 5.880492957780049, 5.880494420043325, 10.686025216519289, 10.69391907743567, 10.693989757202713, 12.317666113412328, 12.320210206425758, 13.360871916742985, 13.42281682316593, 13.425650479449894, 14.216103594811019]
11 : [3.220289763758681, 5.880492957779937, 5.880494420042617, 10.686006011051642, 10.693917414577932, 10.693986570960563, 12.317664492531788, 12.318659204387547, 13.358849296248938, 13.422344036927779, 13.423640054482668, 14.215520284444674]
12 : [3.220289763758688, 5.880492957779877, 5.880494420042572, 10.686001612157106, 10.69391702161566, 10.69398585921227, 12.317664119925036, 12.318123653136919, 13.35793069207255, 13.422136539367523, 13.422748129556258, 14.215203524225455]
13 : [3.2202897637586863, 5.880492957779889, 5.880494420042553, 10.686000605536007, 10.693916931075739, 10.69398569752855, 12.317664014900231, 12.317939018840637, 13.357518473304975, 13.42204424153143, 13.422347216123875, 14.215030360336803]
14 : [3.2202897637586863, 5.880492957779907, 5.8804944200425675, 10.686000374943365, 10.693916910320558, 10.693985660604064, 12.317663977185925, 12.317875376220421, 13.357334414725615, 13.422002986409387, 13.422165791758907, 14.21493535333446]
15 : [3.2202897637586867, 5.880492957779877, 5.880494420042704, 10.686000322114799, 10.69391690555671, 10.693985652161881, 12.317663962996349, 12.317853421919759, 13.357252357509015, 13.421984355470304, 13.422083554521095, 14.214883148928891]
16 : [3.2202897637586867, 5.880492957779741, 5.880494420042556, 10.686000309990654, 10.693916904464098, 10.69398565022707, 12.317663957903555, 12.317845844115817, 13.357215783084866, 13.421975921576848, 13.422046203846026, 14.214854388932642]
17 : [3.2202897637586836, 5.880492957779618, 5.880494420042574, 10.686000307202386, 10.693916904212163, 10.693985649782874, 12.317663956134945, 12.317843218924912, 13.35719944027143, 13.421972067148085, 13.422029221579058, 14.214838532959808]
18 : [3.220289763758684, 5.880492957780035, 5.880494420042671, 10.686000306557714, 10.693916904154475, 10.693985649679774, 12.317663955514492, 12.317842303080296, 13.35719201358384, 13.421970322050447, 13.4220214169237, 14.214829385321448]
19 : [3.220289763758684, 5.88049295777995, 5.880494420042747, 10.686000306410342, 10.693916904141057, 10.693985649656964, 12.317663955308186, 12.317841990282167, 13.357188794210597, 13.421969521537594, 13.422017699744412, 14.214824204303405]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.220289763758684
5.88049295777995
5.880494420042747
10.686000306410342
10.693916904141057
10.693985649656964
12.317663955308186
12.317841990282167
13.357188794210597
13.421969521537594
13.422017699744412
14.214824204303405
[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);
[ ]:

[ ]: