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.213129814810902, 17.799813243902125, 37.219825643854705, 44.01930033100344, 62.18233337164603, 65.53333261764435, 73.98716308482844, 81.08231277458596, 88.06001748067199, 91.2505569586067, 113.61980591470005, 122.21673392687063]
1 : [3.2380975391803895, 5.998650630552989, 6.161390072060389, 11.485860551539993, 12.037309680038101, 14.413520502950268, 15.23481450053956, 16.218269373059275, 17.34672695377123, 19.892283927649178, 22.638813214614746, 24.807560957625814]
2 : [3.220542688575918, 5.8848767088998395, 5.888226095524639, 10.785471632156153, 10.847156783175318, 11.381986378749364, 13.343238319642175, 13.372706506894474, 14.112393286410274, 14.98863956319568, 15.344080016830755, 17.688329576942074]
3 : [3.2202954145601557, 5.880713446554535, 5.880839526591103, 10.7074665263108, 10.720069178399498, 10.944773311455418, 12.52797395481206, 13.075036551246054, 13.54069791591705, 14.030695748797296, 14.367494218159667, 15.906741589800728]
4 : [3.220289894102197, 5.8805022532910565, 5.880518097321826, 10.69511436307729, 10.699129958355796, 10.787751812017861, 12.35409837070715, 12.871053039767272, 13.445683487795225, 13.649403154566546, 14.213859960544362, 14.60766545301286]
5 : [3.220289766421211, 5.880493375254687, 5.880495963116369, 10.692776367609241, 10.69499950969574, 10.721906317529887, 12.324202578881533, 12.640658946735144, 13.417157954203887, 13.509762215045255, 13.816759278689219, 14.260863115192743]
6 : [3.2202897638070196, 5.880492978117846, 5.880494510444939, 10.691097181825258, 10.69416115365378, 10.698643594421696, 12.31888328176604, 12.46074660581755, 13.402394959854009, 13.454795875839899, 13.56329326412352, 14.232709516141071]
7 : [3.220289763759514, 5.880492958816452, 5.880494425037931, 10.687861627692293, 10.693985793780012, 10.69455279230062, 12.317899571105842, 12.372364426895516, 13.387450842731804, 13.433939412852363, 13.47327115838431, 14.223431358575272]
8 : [3.2202897637586982, 5.880492957834425, 5.880494420312181, 10.686464368105478, 10.693942498545255, 10.694085266049356, 12.317711055329656, 12.337392167734635, 13.373735427664625, 13.426640065252624, 13.441852510692199, 14.219267826315182]
9 : [3.2202897637586867, 5.880492957782823, 5.880494420056991, 10.686108698959108, 10.693925371312508, 10.694004893998786, 12.317673740213866, 12.324680164846331, 13.365194188477124, 13.423921247804426, 13.430320729553117, 14.21719102389813]
10 : [3.2202897637586863, 5.880492957780064, 5.880494420043325, 10.68602521651859, 10.69391907743549, 10.693989757202536, 12.317666113412294, 12.320210206208825, 13.360871916429053, 13.422816823081954, 13.425650479488674, 14.216103594756133]
11 : [3.2202897637586863, 5.880492957779908, 5.880494420042602, 10.686006011068569, 10.693917414577362, 10.693986570962709, 12.317664492530282, 12.318659204151569, 13.358849295633044, 13.422344037669054, 13.423640059898736, 14.215520284772131]
12 : [3.2202897637586867, 5.880492957779926, 5.880494420042568, 10.686001612168873, 10.69391702161663, 10.693985859213564, 12.317664119922918, 12.31812365389094, 13.35793069195704, 13.42213654088213, 13.422748135476905, 14.21520352309928]
13 : [3.2202897637586836, 5.8804929577798895, 5.8804944200425595, 10.686000605582155, 10.693916931076263, 10.693985697537817, 12.31766401489975, 12.317939028394433, 13.357518505607366, 13.422044243213497, 13.422347248667066, 14.215030384982496]
14 : [3.220289763758684, 5.88049295777991, 5.88049442004257, 10.686000374992538, 10.693916910320464, 10.69398566061332, 12.317663977184607, 12.31787538934751, 13.357334463879312, 13.422002958126363, 13.422165845597721, 14.214935392774631]
15 : [3.220289763758688, 5.880492957779812, 5.88049442004256, 10.686000322147441, 10.693916905558172, 10.693985652167449, 12.317663962996681, 12.31785343638188, 13.357252416916277, 13.421984334517326, 13.422083627089636, 14.214883183997506]
16 : [3.220289763758685, 5.880492957779705, 5.880494420042822, 10.686000309999375, 10.69391690446416, 10.693985650228866, 12.317663957903463, 12.317845846645765, 13.357215796563702, 13.421975896393997, 13.422046249826701, 14.214854408976791]
17 : [3.220289763758684, 5.880492957779747, 5.880494420042659, 10.686000307203559, 10.693916904212555, 10.693985649783151, 12.317663956129909, 12.317843217685617, 13.357199424652908, 13.421972051407526, 13.42202923633643, 14.214838534027349]
18 : [3.220289763758687, 5.88049295777996, 5.880494420042344, 10.686000306553487, 10.693916904154458, 10.693985649679183, 12.31766395550558, 12.317842299904912, 13.357192014881209, 13.421970310479683, 13.422021287000995, 14.214829581167146]
19 : [3.220289763758684, 5.880492957779892, 5.880494420042291, 10.686000306411504, 10.693916904141007, 10.693985649657087, 12.31766395530728, 12.317841990362073, 13.357188810990765, 13.421969484221943, 13.42201782652132, 14.214824741956553]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.220289763758684
5.880492957779892
5.880494420042291
10.686000306411504
10.693916904141007
10.693985649657087
12.31766395530728
12.317841990362073
13.357188810990765
13.421969484221943
13.42201782652132
14.214824741956553
[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);
[ ]:

[ ]: