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.213129797839025, 17.799813245452093, 37.21982565595732, 44.01930032658754, 62.182333374214764, 65.53333261759971, 73.98716308295924, 81.08231275881637, 88.06001748044707, 91.25055693315022, 113.61980595044477, 122.2167339821912]
1 : [3.2380975391312985, 5.99865063057926, 6.161390072266458, 11.485860551139936, 12.037309679498554, 14.413520501380338, 15.234814499724418, 16.218269374465624, 17.346726951925316, 19.89228393076782, 22.638813222742137, 24.80756097207125]
2 : [3.2205426885753754, 5.884876708898107, 5.888226095533873, 10.785471632126656, 10.847156783086653, 11.381986378890433, 13.343238319806527, 13.372706506582368, 14.112393285333393, 14.988639565281712, 15.344080017062911, 17.688329581098987]
3 : [3.2202954145601437, 5.8807134465543145, 5.880839526591737, 10.707466526319951, 10.72006917837889, 10.94477331143944, 12.527973954731271, 13.075036550987846, 13.540697915991807, 14.030695748739928, 14.367494219021696, 15.906741588696338]
4 : [3.220289894102197, 5.880502253291016, 5.880518097321827, 10.69511436307983, 10.69912995835074, 10.787751811974744, 12.354098370699667, 12.871053039250063, 13.44568348780333, 13.649403154551276, 14.21385996148939, 14.607665450743776]
5 : [3.220289766421206, 5.880493375254663, 5.88049596311637, 10.692776367608365, 10.694999509694442, 10.721906317496668, 12.324202578880534, 12.640658946137739, 13.41715795417091, 13.509762214987143, 13.81675927857197, 14.260863114705323]
6 : [3.220289763807026, 5.880492978117791, 5.880494510444949, 10.691097181820993, 10.694161153653463, 10.698643594412285, 12.318883281765697, 12.460746605467866, 13.402394959780207, 13.454795875781764, 13.563293264149312, 14.232709515994964]
7 : [3.2202897637595114, 5.8804929588164905, 5.880494425037932, 10.687861627688866, 10.693985793779873, 10.69455279229991, 12.317899571105619, 12.372364426744639, 13.387450842641297, 13.433939412823628, 13.473271158460577, 14.22343135851116]
8 : [3.2202897637586982, 5.880492957834403, 5.880494420312195, 10.686464368104636, 10.693942498545224, 10.69408526604934, 12.317711055329626, 12.337392167688488, 13.373735427609335, 13.42664006524451, 13.441852510746576, 14.2192678262831]
9 : [3.220289763758683, 5.8804929577828124, 5.880494420057005, 10.686108698962284, 10.693925371312636, 10.694004893999361, 12.31767374021425, 12.324680164955074, 13.365194188502096, 13.423921247793018, 13.43032072976305, 14.217191023840018]
10 : [3.220289763758686, 5.880492957780064, 5.880494420043331, 10.686025216521502, 10.6939190774358, 10.693989757203077, 12.317666113412544, 12.320210206574659, 13.360871916857281, 13.422816823169626, 13.425650479705267, 14.216103594787553]
11 : [3.220289763758685, 5.8804929577799125, 5.8804944200426075, 10.686006011066855, 10.693917414577816, 10.693986570963794, 12.317664492529824, 12.318659205640344, 13.35884929898492, 13.422344037112026, 13.423640057409933, 14.215520285963338]
12 : [3.220289763758687, 5.880492957779897, 5.880494420042563, 10.686001612160362, 10.693917021615881, 10.693985859212793, 12.317664119923794, 12.318123653114077, 13.35793069164358, 13.422136540025255, 13.42274813137757, 14.215203524054079]
13 : [3.2202897637586863, 5.880492957779892, 5.880494420042562, 10.68600060555768, 10.693916931073415, 10.693985697535425, 12.317664014901364, 12.31793902111506, 13.357518493794968, 13.422044246015504, 13.422347234344867, 14.21503039521526]
14 : [3.2202897637586854, 5.88049295777992, 5.880494420042564, 10.686000375023472, 10.69391691032831, 10.693985660620871, 12.317663977207472, 12.317875423188681, 13.35733458833521, 13.422003004818535, 13.422165852303545, 14.214935428083576]
15 : [3.2202897637586827, 5.88049295777999, 5.880494420042588, 10.686000322164292, 10.693916905561629, 10.693985652171548, 12.317663963011585, 12.317853461314138, 13.357252534318727, 13.421984375405671, 13.42208364132284, 14.214883221940733]
16 : [3.220289763758684, 5.88049295777995, 5.880494420042651, 10.686000310008684, 10.693916904465508, 10.693985650230625, 12.317663957913231, 12.317845864239988, 13.357215905851495, 13.421975921987206, 13.422046265523147, 14.214854458051546]
17 : [3.220289763758682, 5.880492957780004, 5.880494420042693, 10.686000307206504, 10.69391690421285, 10.693985649783201, 12.317663956131376, 12.31784323072616, 13.357199510143563, 13.421972045340597, 13.422029200111544, 14.214838451673788]
18 : [3.2202897637586836, 5.880492957780281, 5.880494420042477, 10.68600030656184, 10.693916904155296, 10.69398564968052, 12.317663955520661, 12.317842318907855, 13.357192173880755, 13.421970297827638, 13.422021374487517, 14.214829513815108]
19 : [3.220289763758685, 5.880492957779929, 5.880494420042712, 10.686000306413192, 10.69391690414147, 10.693985649657101, 12.31766395529756, 12.317841982926176, 13.35718856229564, 13.421969448289605, 13.422017866203758, 14.214823991745511]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.220289763758685
5.880492957779929
5.880494420042712
10.686000306413192
10.69391690414147
10.693985649657101
12.31766395529756
12.317841982926176
13.35718856229564
13.421969448289605
13.422017866203758
14.214823991745511
[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);
[ ]:

[ ]: