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 = 42562
[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 : [6.229093264786797, 22.636447499331588, 38.75912680943771, 56.95688373026947, 64.68904675320641, 71.61510643467192, 76.89311294238219, 86.36950975567923, 90.92730267905193, 99.90890224196235, 108.07709801259793, 123.53078334428623]
1 : [3.231012675434386, 6.086234263756916, 6.4068992240264135, 12.054533474399893, 12.672421029188746, 13.138599688057516, 14.490032248052795, 17.723681706145108, 18.66258439847152, 19.34612060373716, 23.999160959618045, 25.451852238804417]
2 : [3.220379348186532, 5.886623822610176, 5.9055484292926375, 10.853747390778707, 10.901287970138593, 11.033338300582972, 12.49060860995537, 13.713193967572463, 14.184813087423747, 15.348674545838609, 16.529828075792313, 17.964540349012445]
3 : [3.2202534507267346, 5.880768708677711, 5.8819084880277615, 10.709133420835416, 10.726052065985265, 10.778656219112557, 12.345960396969497, 12.76067604305868, 13.76028135563162, 14.018326042873763, 14.607718515744535, 16.592214150835765]
4 : [3.2202514669005917, 5.880494279982, 5.880556917690733, 10.689894953485688, 10.698983267280639, 10.718281105552185, 12.323621310865793, 12.515134575083456, 13.530915636120715, 13.646352490591044, 14.013638403995914, 15.844042103319989]
5 : [3.2202514327806298, 5.8804771145376336, 5.880480023263912, 10.68685103655616, 10.694700891024908, 10.700890407059001, 12.319109793664776, 12.416240501302422, 13.409477889898788, 13.552567887202976, 13.74921631124324, 15.246452468550515]
6 : [3.2202514321849858, 5.880473890511299, 5.880477843318854, 10.686161008730863, 10.694032876757815, 10.695820538591711, 12.318002400652832, 12.365108987094823, 13.375722155054579, 13.48849701741845, 13.613981824841636, 14.81374636569212]
7 : [3.220251432174624, 5.880473676484472, 5.880477751399403, 10.685969287196476, 10.693915072917044, 10.694415397844935, 12.317696219797625, 12.338610207937437, 13.365112608787404, 13.45294954024236, 13.531662728947383, 14.548118346320106]
8 : [3.2202514321744435, 5.880473664772776, 5.880477746015275, 10.685917731763425, 10.693892108137462, 10.694043457650851, 12.317604665085188, 12.326098068794307, 13.360878384965286, 13.435960447046861, 13.480566199399348, 14.398036525551326]
9 : [3.220251432174438, 5.880473664135101, 5.88047774570594, 10.685904872914914, 10.693887275314236, 10.693949433864502, 12.317575537838119, 12.320797130649877, 13.358877986622415, 13.428258104112944, 13.45152541697752, 14.315735102773868]
10 : [3.2202514321744427, 5.8804736641003235, 5.880477745688486, 10.685901792194208, 10.693886175010338, 10.693926453072484, 12.317565160998337, 12.318715754951024, 13.357892966935353, 13.424793973345572, 13.43629871091776, 14.270543597884595]
11 : [3.2202514321744435, 5.880473664098424, 5.8804777456875215, 10.685901066925268, 10.693885916583676, 10.693920952759527, 12.317558831688624, 12.31793983446391, 13.357414569658118, 13.423228032755313, 13.42874752200064, 14.245626540942153]
12 : [3.2202514321744404, 5.880473664098306, 5.880477745687478, 10.685900897346439, 10.693885856844595, 10.69391964867569, 12.317545562951777, 12.31766812760482, 13.357188433549366, 13.422515874489237, 13.42512010873811, 14.231830246772033]
13 : [3.2202514321744418, 5.880473664098295, 5.880477745687438, 10.685900857783231, 10.693885843152593, 10.693919341011567, 12.317512857731625, 12.31759852340779, 13.357083745405962, 13.422189421102173, 13.423410850672042, 14.224177354948992]
14 : [3.2202514321744364, 5.880473664098297, 5.880477745687466, 10.685900848561475, 10.693885840014234, 10.693919268621174, 12.3174884939662, 12.317586421941511, 13.357035950186237, 13.422036722124416, 13.422615989066411, 14.219928043222946]
15 : [3.220251432174443, 5.880473664097927, 5.88047774568732, 10.685900846410446, 10.693885839292212, 10.693919251606733, 12.317478500002291, 12.317583475139696, 13.357014268837533, 13.421959802552376, 13.422253968945626, 14.217568510968347]
16 : [3.22025143217444, 5.880473664098405, 5.880477745687624, 10.685900845908563, 10.693885839125452, 10.69391924760757, 12.317474829001524, 12.31758255975832, 13.357004464788968, 13.421913917911771, 13.422096883928033, 14.216257562147337]
17 : [3.220251432174439, 5.880473664098572, 5.88047774568745, 10.685900845788792, 10.693885839086304, 10.693919246652875, 12.317473505472828, 12.317582249192798, 13.356999723023694, 13.421882834441845, 13.422032549615096, 14.21551257604324]
18 : [3.2202514321744404, 5.8804736640983135, 5.880477745687569, 10.685900845762703, 10.69388583907813, 10.693919246440261, 12.317473036954762, 12.31758214028254, 13.356997790964977, 13.421865175158825, 13.422006571841175, 14.215103407638045]
19 : [3.220251432174444, 5.880473664098494, 5.880477745687685, 10.685900845757118, 10.693885839075776, 10.693919246392888, 12.317472876526603, 12.317582103160868, 13.35699701366591, 13.421856349482027, 13.421996197410696, 14.21488542402547]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.220251432174444
5.880473664098494
5.880477745687685
10.685900845757118
10.693885839075776
10.693919246392888
12.317472876526603
12.317582103160868
13.35699701366591
13.421856349482027
13.421996197410696
14.21488542402547
[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);
[ ]:

[ ]: