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.229093264277982, 22.636447499802912, 38.75912681158575, 56.956883735260334, 64.68904674711972, 71.61510642967278, 76.89311293512647, 86.36950975433056, 90.92730268043083, 99.90890224554285, 108.07709802023624, 123.53078334425769]
1 : [3.231012675429531, 6.086234263778436, 6.406899223992413, 12.05453347408023, 12.672421029304479, 13.138599688429236, 14.490032247483782, 17.723681706679734, 18.662584397579266, 19.346120605774185, 23.9991609554252, 25.45185223957955]
2 : [3.2203793481864733, 5.886623822609621, 5.905548429289152, 10.853747390803784, 10.901287970057092, 11.033338300641457, 12.490608609974915, 13.713193967690085, 14.184813087415044, 15.348674546287787, 16.529828073929174, 17.964540348841762]
3 : [3.2202534507267355, 5.88076870867766, 5.881908488027745, 10.70913342083713, 10.726052065970118, 10.778656219141713, 12.345960396976684, 12.760676043090939, 13.760281355656884, 14.018326042938515, 14.60771851542716, 16.592214149734218]
4 : [3.220251466900595, 5.880494279981976, 5.880556917690763, 10.689894953484808, 10.698983267278567, 10.718281105562331, 12.323621310867868, 12.515134575094493, 13.530915636167515, 13.646352490530296, 14.013638403886539, 15.844042101722918]
5 : [3.2202514327806324, 5.880477114537657, 5.880480023263946, 10.686851036555982, 10.694700891024514, 10.700890407061122, 12.31910979366545, 12.416240501289057, 13.409477889914681, 13.552567887170524, 13.749216311116836, 15.246452467003156]
6 : [3.2202514321849827, 5.880473890511275, 5.880477843318859, 10.686161008730823, 10.694032876757708, 10.695820538591848, 12.318002400653109, 12.365108987072768, 13.37572215506604, 13.488497017400233, 13.613981824686737, 14.813746364568955]
7 : [3.220251432174622, 5.880473676484473, 5.8804777513994155, 10.685969287196446, 10.693915072917001, 10.694415397844834, 12.317696219797714, 12.338610207920315, 13.36511260879578, 13.452949540232874, 13.531662728814423, 14.548118345633787]
8 : [3.220251432174444, 5.8804736647727704, 5.880477746015265, 10.685917731763409, 10.693892108137465, 10.694043457650732, 12.31760466508523, 12.326098068790687, 13.360878384971578, 13.435960447022438, 13.480566199374312, 14.39803652533662]
9 : [3.22025143217444, 5.880473664135106, 5.880477745705927, 10.68590487291513, 10.69388727531327, 10.693949433866038, 12.317575537838156, 12.320797130965559, 13.358877986636571, 13.428258102720054, 13.451525414021564, 14.315735107549965]
10 : [3.2202514321744435, 5.880473664100354, 5.880477745688493, 10.685901792194315, 10.69388617500982, 10.693926453073361, 12.317565160998566, 12.318715755553415, 13.357892966964478, 13.424793972346164, 13.4362987136111, 14.270543611463273]
11 : [3.2202514321744418, 5.880473664098424, 5.880477745687531, 10.685901066925096, 10.693885916583477, 10.693920952759415, 12.317558831686492, 12.317939834804788, 13.357414569617045, 13.423228032156642, 13.428747524677883, 14.245626548623898]
12 : [3.220251432174445, 5.880473664098313, 5.8804777456874735, 10.68590089734641, 10.693885856844544, 10.693919648675603, 12.317545562962225, 12.317668127737036, 13.35718843347244, 13.422515873282645, 13.425120107087324, 14.231830249653564]
13 : [3.220251432174443, 5.880473664098291, 5.8804777456874655, 10.685900857789308, 10.693885843151593, 10.693919341040987, 12.317512863770894, 12.317598526717084, 13.35708376978225, 13.42218946844301, 13.423410993750789, 14.224178986102617]
14 : [3.220251432174439, 5.880473664098318, 5.880477745687464, 10.68590084856409, 10.693885840012808, 10.69391926863515, 12.317488503877529, 12.317586425101025, 13.35703596775081, 13.422036762213263, 13.422616319376095, 14.219930648164189]
15 : [3.220251432174444, 5.880473664098365, 5.880477745687594, 10.68590084641106, 10.693885839291417, 10.693919251608833, 12.317478502628653, 12.317583475761403, 13.357014277601195, 13.421959834223605, 13.422254042329437, 14.217569780666647]
16 : [3.2202514321744418, 5.880473664098143, 5.8804777456875295, 10.685900845908657, 10.693885839124954, 10.693919247606422, 12.317474826992704, 12.317582559069916, 13.357004470681831, 13.421913926246695, 13.422096803452252, 14.21625740703359]
17 : [3.220251432174443, 5.880473664098409, 5.880477745687611, 10.685900845791183, 10.693885839086823, 10.693919246667397, 12.317473512026504, 12.317582249355176, 13.357000038821333, 13.421883763082478, 13.422032803234428, 14.215527738131698]
18 : [3.220251432174443, 5.880473664098331, 5.880477745687434, 10.685900845763813, 10.693885839077879, 10.693919246446185, 12.317473045477639, 12.317582141696441, 13.356998033912353, 13.421865993170648, 13.422007220272743, 14.215122448191794]
19 : [3.220251432174442, 5.880473664098471, 5.8804777456874415, 10.685900845757295, 10.693885839075817, 10.693919246394232, 12.317472880594622, 12.317582104006283, 13.356997128915816, 13.421856769706874, 13.421996658119243, 14.214897366478427]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.220251432174442
5.880473664098471
5.8804777456874415
10.685900845757295
10.693885839075817
10.693919246394232
12.317472880594622
12.317582104006283
13.356997128915816
13.421856769706874
13.421996658119243
14.214897366478427
[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);
[ ]:

[ ]: