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.229093260788123, 22.636447511071207, 38.75912681600302, 56.956883725318114, 64.68904673997548, 71.61510643628256, 76.89311294409626, 86.36950976884285, 90.9273026774525, 99.90890224489313, 108.07709802452655, 123.53078334014874]
1 : [3.2310126754147754, 6.0862342638559355, 6.406899224303705, 12.054533473497248, 12.672421029075874, 13.138599689678415, 14.490032248824711, 17.7236817053954, 18.662584390265398, 19.34612060469177, 23.999160967934355, 25.45185223676891]
2 : [3.2203793481862952, 5.886623822608235, 5.905548429306962, 10.853747390756057, 10.90128797035772, 11.03333830068535, 12.490608609970034, 13.713193965032499, 14.184813087312405, 15.348674545954594, 16.529828081487544, 17.96454034866411]
3 : [3.2202534507267297, 5.8807687086776115, 5.881908488028795, 10.709133420861406, 10.72605206601771, 10.778656219106141, 12.34596039696703, 12.760676042574568, 13.760281355645409, 14.018326043753582, 14.60771851594105, 16.5922141507136]
4 : [3.220251466900592, 5.880494279981977, 5.880556917690789, 10.689894953490706, 10.698983267287618, 10.718281105542587, 12.32362131086322, 12.515134574929824, 13.53091563641404, 13.646352490362474, 14.013638403897037, 15.844042103171512]
5 : [3.2202514327806293, 5.880477114537636, 5.8804800232639165, 10.686851036557169, 10.694700891025645, 10.700890407055638, 12.319109793663552, 12.416240501232421, 13.409477889959884, 13.55256788705854, 13.749216311167686, 15.24645246845209]
6 : [3.220251432184983, 5.880473890511292, 5.88047784331888, 10.686161008730954, 10.694032876757888, 10.695820538590763, 12.318002400652421, 12.365108987060129, 13.375722155082535, 13.488497017331632, 13.6139818247799, 14.813746365629271]
7 : [3.2202514321746207, 5.880473676484502, 5.880477751399417, 10.685969287196459, 10.693915072917044, 10.694415397844688, 12.317696219797492, 12.338610207921164, 13.365112608809651, 13.452949540194703, 13.53166272890025, 14.548118346282035]
8 : [3.2202514321744427, 5.880473664772781, 5.88047774601527, 10.685917731763427, 10.693892108137456, 10.694043457650688, 12.317604665085158, 12.326098068791708, 13.360878384983213, 13.43596044699015, 13.480566199393751, 14.398036525710378]
9 : [3.22025143217444, 5.8804736641351365, 5.8804777457059245, 10.685904872914495, 10.693887275312509, 10.693949433865962, 12.317575537840117, 12.320797131411442, 13.35887798664352, 13.428258104433016, 13.451525420810883, 14.31573510927446]
10 : [3.2202514321744404, 5.880473664100308, 5.880477745688507, 10.685901792194025, 10.69388617500938, 10.693926453073349, 12.317565160998464, 12.318715755213493, 13.357892967075102, 13.424793971212322, 13.4362987009656, 14.270543594917024]
11 : [3.2202514321744413, 5.88047366409843, 5.8804777456875295, 10.68590106692528, 10.69388591658333, 10.693920952759937, 12.317558831690768, 12.317939834566186, 13.357414569858395, 13.423228030768632, 13.428747511643218, 14.245626537764528]
12 : [3.2202514321744413, 5.880473664098311, 5.88047774568748, 10.685900897346423, 10.69388585684436, 10.693919648675891, 12.317545562958749, 12.317668127658042, 13.35718843385294, 13.422515873964933, 13.42512010239693, 14.231830240392359]
13 : [3.2202514321744427, 5.880473664098334, 5.880477745687491, 10.685900857788631, 10.693885843151838, 10.693919341041347, 12.317512862863543, 12.317598526100081, 13.357083769772633, 13.42218947176141, 13.4234109050735, 14.22417870622686]
14 : [3.2202514321744404, 5.880473664098455, 5.880477745687593, 10.685900848563357, 10.69388584001291, 10.693919268635312, 12.317488501761694, 12.317586424292143, 13.357035967569066, 13.422036758264873, 13.422616088159732, 14.219930010906149]
15 : [3.2202514321744444, 5.880473664098414, 5.880477745687939, 10.685900846411302, 10.693885839291456, 10.693919251610806, 12.317478504143633, 12.317583476078712, 13.357014278145156, 13.42195984107766, 13.422253979438006, 14.217570048473807]
16 : [3.2202514321744453, 5.880473664098231, 5.8804777456876565, 10.685900845908645, 10.693885839125256, 10.693919247605423, 12.317474825286899, 12.317582558591509, 13.357004471254623, 13.421913902636083, 13.422096827498756, 14.21625706432595]
17 : [3.220251432174443, 5.880473664098395, 5.880477745687471, 10.68590084579092, 10.693885839086697, 10.693919246664057, 12.31747350824187, 12.317582248759201, 13.357000018943838, 13.421883584950073, 13.42203319311098, 14.215524476848207]
18 : [3.22025143217444, 5.880473664098378, 5.880477745687514, 10.685900845763689, 10.693885839077845, 10.693919246444649, 12.317473040644924, 12.317582140689035, 13.35699803013528, 13.421865765399149, 13.422007487639558, 14.21511702322347]
19 : [3.2202514321744404, 5.880473664098393, 5.880477745687515, 10.685900845757352, 10.693885839075888, 10.693919246393927, 12.317472877077174, 12.317582103127172, 13.356997132939199, 13.421856588836272, 13.421996733906285, 14.214892221639927]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.2202514321744404
5.880473664098393
5.880477745687515
10.685900845757352
10.693885839075888
10.693919246393927
12.317472877077174
12.317582103127172
13.356997132939199
13.421856588836272
13.421996733906285
14.214892221639927
[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);
[ ]:

[ ]: