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.213129796235323, 17.799813243418246, 37.219825629339184, 44.01930032911768, 62.18233337507144, 65.53333261620158, 73.98716307614761, 81.08231276921508, 88.06001749391193, 91.250556937817, 113.61980593054172, 122.21673394672841]
1 : [3.2380975391162212, 5.998650630547501, 6.16139007178734, 11.485860551430559, 12.037309680994456, 14.41352049969458, 15.234814500126422, 16.218269373059293, 17.34672695269671, 19.892283936385937, 22.63881321085537, 24.807560973195507]
2 : [3.2205426885751627, 5.884876708898056, 5.888226095523973, 10.785471632156723, 10.847156783331846, 11.381986378449936, 13.343238318572071, 13.372706506976664, 14.112393285600902, 14.988639565142904, 15.344080019772324, 17.688329576796352]
3 : [3.220295414560139, 5.880713446554296, 5.880839526591342, 10.707466526320491, 10.72006917842288, 10.944773311188973, 12.527973954550417, 13.07503655062608, 13.540697915979589, 14.030695750150903, 14.367494219084799, 15.906741582284736]
4 : [3.220289894102193, 5.880502253291035, 5.8805180973217945, 10.695114363079934, 10.699129958359517, 10.787751811775168, 12.35409837066802, 12.871053038487599, 13.44568348779125, 13.649403155336389, 14.213859960407099, 14.607665446448593]
5 : [3.220289766421209, 5.880493375254647, 5.880495963116386, 10.69277636760799, 10.694999509696292, 10.721906317381878, 12.324202578874866, 12.640658945147148, 13.417157954167502, 13.509762215342992, 13.816759276079692, 14.260863114739823]
6 : [3.220289763807025, 5.880492978117858, 5.880494510444944, 10.69109718181325, 10.694161153653853, 10.698643594377259, 12.318883281764611, 12.460746604846095, 13.402394959770982, 13.454795875917291, 13.563293263216828, 14.232709516043275]
7 : [3.220289763759515, 5.880492958816489, 5.88049442503792, 10.687861627680963, 10.693985793780035, 10.694552792296019, 12.31789957110542, 12.37236442648047, 13.387450842602162, 13.43393941286365, 13.473271158147854, 14.223431358537056]
8 : [3.2202897637586956, 5.880492957834426, 5.880494420312175, 10.686464368102667, 10.69394249854524, 10.694085266048722, 12.31771105532956, 12.337392167593014, 13.373735427566896, 13.4266400652557, 13.441852510656744, 14.219267826296532]
9 : [3.220289763758685, 5.880492957782818, 5.880494420056987, 10.686108698964775, 10.693925371312737, 10.694004893999631, 12.31767374021444, 12.324680165045981, 13.365194188570747, 13.423921247788515, 13.430320729913966, 14.217191023849617]
10 : [3.220289763758683, 5.880492957780054, 5.880494420043363, 10.686025216524257, 10.693919077435933, 10.69398975720348, 12.317666113412638, 12.320210206732831, 13.360871917041921, 13.42281682314608, 13.42565047998606, 14.216103594785002]
11 : [3.22028976375868, 5.880492957779936, 5.880494420042594, 10.686006011058428, 10.693917414576964, 10.693986570961044, 12.317664492529588, 12.318659203074173, 13.358849294021013, 13.422344037954195, 13.423640058016565, 14.215520284372927]
12 : [3.2202897637586867, 5.88049295777991, 5.880494420042564, 10.6860016121589, 10.693917021615468, 10.693985859211796, 12.317664119923428, 12.318123651099796, 13.35793068603647, 13.422136539806765, 13.422748132979873, 14.215203521850505]
13 : [3.2202897637586823, 5.880492957779911, 5.8804944200425675, 10.686000605455035, 10.69391693107295, 10.693985697513943, 12.317664014897096, 12.317938999074814, 13.357518416006842, 13.422044244223734, 13.422347163371432, 14.215030342088987]
14 : [3.220289763758686, 5.880492957779919, 5.880494420042581, 10.686000374949048, 10.693916910325541, 10.6939856606041, 12.317663977181521, 12.3178753840507, 13.357334419661173, 13.422002992870498, 13.42216577601221, 14.214935325681582]
15 : [3.220289763758685, 5.880492957779887, 5.8804944200425675, 10.686000322135136, 10.693916905560412, 10.693985652164876, 12.31766396299326, 12.31785343752128, 13.35725240457298, 13.421984366948308, 13.42208358275011, 14.214883131428477]
16 : [3.2202897637586854, 5.880492957779936, 5.8804944200424805, 10.686000309998796, 10.693916904465151, 10.693985650228582, 12.317663957904195, 12.317845851967109, 13.357215818399606, 13.421975917978266, 13.422046231164451, 14.214854384320276]
17 : [3.2202897637586845, 5.880492957779782, 5.880494420042707, 10.686000307209094, 10.693916904213078, 10.69398564978409, 12.317663956137588, 12.317843229979236, 13.357199513322408, 13.42197206654501, 13.422029266786867, 14.214838545430915]
18 : [3.2202897637586845, 5.880492957780031, 5.880494420042811, 10.686000306565075, 10.69391690415522, 10.693985649681595, 12.317663955528054, 12.317842318141103, 13.357192198485395, 13.421970312766664, 13.42202152650815, 14.214829753154387]
19 : [3.220289763758685, 5.8804929577804135, 5.8804944200428775, 10.686000306416133, 10.693916904141613, 10.693985649657707, 12.317663955306289, 12.317841994637316, 13.357188648001015, 13.421969483152399, 13.422017977390036, 14.214823277418597]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.220289763758685
5.8804929577804135
5.8804944200428775
10.686000306416133
10.693916904141613
10.693985649657707
12.317663955306289
12.317841994637316
13.357188648001015
13.421969483152399
13.422017977390036
14.214823277418597
[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);
[ ]:

[ ]: