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.213129806727155, 17.799813246171507, 37.21982563181586, 44.019300328118916, 62.18233337208704, 65.5333326153055, 73.98716309443921, 81.08231276453367, 88.06001748423351, 91.25055695535691, 113.61980592261085, 122.21673391862669]
1 : [3.2380975391374576, 5.998650630630013, 6.16139007189183, 11.48586055144063, 12.037309680533786, 14.413520502306, 15.234814500358011, 16.218269371904988, 17.346726951675386, 19.892283931578973, 22.638813208067248, 24.807560965834046]
2 : [3.220542688575506, 5.884876708903952, 5.888226095522143, 10.785471632207017, 10.847156783313123, 11.381986378280175, 13.343238318428432, 13.372706506974556, 14.112393286000641, 14.988639564949683, 15.344080017850985, 17.68832957934073]
3 : [3.22029541456015, 5.880713446554668, 5.88083952659149, 10.707466526328178, 10.720069178438738, 10.94477331125635, 12.527973954553072, 13.075036551605512, 13.540697916112041, 14.030695748899324, 14.367494218678663, 15.906741594788988]
4 : [3.220289894102199, 5.880502253291036, 5.880518097321885, 10.695114363080434, 10.699129958366255, 10.787751812017818, 12.354098370663175, 12.871053040490537, 13.44568348781062, 13.649403154701782, 14.213859960809515, 14.607665457626746]
5 : [3.220289766421207, 5.880493375254656, 5.880495963116404, 10.69277636760827, 10.694999509697581, 10.72190631757275, 12.324202578872818, 12.640658947660075, 13.417157954157135, 13.50976221512091, 13.816759280630588, 14.260863115684613]
6 : [3.2202897638070223, 5.880492978117817, 5.880494510444956, 10.691097181825924, 10.694161153654136, 10.698643594442832, 12.318883281764153, 12.460746606407637, 13.402394959795418, 13.454795875894115, 13.563293264978826, 14.232709516292475]
7 : [3.2202897637595145, 5.880492958816466, 5.880494425037958, 10.687861627695765, 10.693985793780174, 10.694552792303863, 12.317899571105366, 12.372364427153565, 13.38745084271235, 13.43393941288889, 13.473271158726254, 14.223431358646096]
8 : [3.220289763758698, 5.8804929578344, 5.88049442031217, 10.686464368106531, 10.693942498545386, 10.694085266049875, 12.317711055329557, 12.33739216783273, 13.373735427674877, 13.42664006527272, 13.441852510818235, 14.21926782635009]
9 : [3.2202897637586836, 5.880492957782835, 5.880494420057004, 10.68610869895862, 10.693925371312664, 10.694004893998287, 12.317673740213948, 12.324680164776533, 13.365194188327177, 13.423921247809812, 13.430320729577739, 14.217191023796335]
10 : [3.220289763758683, 5.8804929577800635, 5.880494420043331, 10.686025216521022, 10.693919077435737, 10.693989757202777, 12.317666113412441, 12.320210206430138, 13.360871916655551, 13.422816823121673, 13.425650479700197, 14.216103594721151]
11 : [3.220289763758682, 5.880492957779913, 5.8804944200425995, 10.686006011052623, 10.693917414577653, 10.693986570960782, 12.31766449252758, 12.318659203909224, 13.358849295046701, 13.422344037202956, 13.423640055465357, 14.215520283974962]
12 : [3.2202897637586836, 5.880492957779899, 5.880494420042562, 10.686001612157838, 10.69391702161518, 10.693985859212356, 12.317664119923972, 12.318123652329817, 13.357930690157687, 13.422136538865997, 13.422748130719544, 14.215203523793761]
13 : [3.2202897637586863, 5.880492957779903, 5.880494420042577, 10.686000605476783, 10.693916931073963, 10.693985697516844, 12.317664014895511, 12.317939002975033, 13.357518423445503, 13.422044243309292, 13.422347180847847, 14.2150303376657]
14 : [3.2202897637586863, 5.880492957779899, 5.880494420042547, 10.686000374954808, 10.693916910325335, 10.693985660605478, 12.317663977169047, 12.31787538877574, 13.357334444515795, 13.422002999211218, 13.422165786145719, 14.214935341364308]
15 : [3.220289763758689, 5.880492957779872, 5.880494420042906, 10.686000322134246, 10.693916905560181, 10.693985652165074, 12.317663962989554, 12.317853437575721, 13.357252411328673, 13.421984367410966, 13.42208358509874, 14.214883144179222]
16 : [3.220289763758684, 5.880492957779866, 5.880494420042446, 10.686000309998184, 10.693916904465329, 10.693985650228392, 12.317663957900072, 12.317845852519826, 13.3572158238935, 13.421975924836888, 13.42204622821272, 14.214854391875019]
17 : [3.220289763758682, 5.880492957779864, 5.880494420042542, 10.68600030720412, 10.693916904213056, 10.693985649782793, 12.31766395612568, 12.31784321999844, 13.357199419923088, 13.42197207107731, 13.422029236995844, 14.214838446695854]
18 : [3.220289763758684, 5.8804929577798, 5.88049442004275, 10.686000306543866, 10.693916904153703, 10.693985649677682, 12.317663955496478, 12.317842272949713, 13.357191755050422, 13.421970288578986, 13.422021164028722, 14.214829459786174]
19 : [3.2202897637586863, 5.880492957779875, 5.880494420042554, 10.686000306407118, 10.693916904141, 10.693985649656122, 12.317663955296958, 12.317841973795131, 13.357188566594925, 13.421969499917257, 13.42201767051331, 14.214824581298927]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.2202897637586863
5.880492957779875
5.880494420042554
10.686000306407118
10.693916904141
10.693985649656122
12.317663955296958
12.317841973795131
13.357188566594925
13.421969499917257
13.42201767051331
14.214824581298927
[7]:
gfu = GridFunction(fes, multidim=len(evecs))
for i in range(len(evecs)):
    gfu.vecs[i].data = evecs[i]
Draw (Norm(gfu), mesh, order=4, min=0, max=2);
[ ]:

[ ]: