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.229093260396661, 22.636447507829757, 38.75912680877755, 56.956883722363564, 64.68904674247491, 71.61510643846283, 76.89311293212295, 86.36950975772999, 90.92730267698535, 99.90890224050986, 108.07709802318475, 123.53078333904877]
1 : [3.2310126754123716, 6.086234263761918, 6.4068992243008775, 12.054533472330963, 12.672421029168293, 13.138599689966297, 14.490032246651685, 17.723681705696244, 18.66258439541556, 19.346120605077232, 23.99916096114943, 25.45185223988894]
2 : [3.220379348186298, 5.886623822606213, 5.905548429301719, 10.85374739065348, 10.901287970013126, 11.033338300811026, 12.490608609873144, 13.713193967647314, 14.184813087325727, 15.348674545490137, 16.529828076212894, 17.964540348826297]
3 : [3.2202534507267297, 5.880768708677535, 5.881908488028286, 10.7091334208513, 10.72605206591037, 10.778656219180665, 12.345960396956587, 12.760676043206018, 13.760281355569163, 14.018326042664066, 14.607718515944775, 16.59221415099389]
4 : [3.2202514669005904, 5.880494279982025, 5.880556917690788, 10.68989495348639, 10.698983267268, 10.718281105576885, 12.323621310862313, 12.515134575165895, 13.530915636017323, 13.646352490601803, 14.013638404060288, 15.844042103663288]
5 : [3.2202514327806284, 5.880477114537629, 5.880480023263916, 10.686851036556765, 10.694700891022675, 10.700890407067263, 12.319109793663646, 12.416240501351044, 13.40947788986183, 13.552567887227433, 13.74921631124587, 15.246452468949428]
6 : [3.220251432184984, 5.880473890511287, 5.88047784331886, 10.686161008731336, 10.694032876757369, 10.69582053859416, 12.318002400652501, 12.365108987122092, 13.375722155040718, 13.488497017441288, 13.613981824830917, 14.813746366000283]
7 : [3.2202514321746243, 5.880473676484502, 5.880477751399428, 10.685969287196642, 10.693915072916942, 10.694415397845649, 12.317696219797513, 12.338610207951282, 13.365112608782576, 13.452949540258194, 13.531662728936784, 14.548118346505388]
8 : [3.2202514321744404, 5.8804736647727704, 5.880477746015278, 10.685917731763475, 10.693892108137446, 10.694043457650974, 12.317604665085183, 12.326098068802494, 13.360878384965604, 13.43596044702812, 13.480566199370621, 14.398036525743851]
9 : [3.220251432174439, 5.8804736641351365, 5.880477745705927, 10.68590487291415, 10.693887275312528, 10.693949433865214, 12.317575537839534, 12.320797130674139, 13.358877986620138, 13.428258103285627, 13.451525409156785, 14.315735095823491]
10 : [3.220251432174441, 5.88047366410034, 5.880477745688504, 10.68590179219364, 10.693886175009638, 10.693926453072773, 12.317565160998493, 12.31871575437014, 13.357892966537275, 13.424793974747583, 13.436298696502254, 14.270543562966761]
11 : [3.2202514321744404, 5.880473664098436, 5.880477745687524, 10.685901066924874, 10.693885916583334, 10.693920952759218, 12.317558831677777, 12.317939833844072, 13.357414569343355, 13.423228034155242, 13.428747502771833, 14.245626489429707]
12 : [3.2202514321744395, 5.880473664098323, 5.8804777456874655, 10.6859008973463, 10.693885856844453, 10.693919648675573, 12.317545562883161, 12.317668127203381, 13.357188433314096, 13.422515874396618, 13.425120084244432, 14.23183018902249]
13 : [3.22025143217444, 5.880473664098307, 5.88047774568749, 10.68590085778901, 10.693885843152517, 10.693919341041742, 12.317512863703719, 12.31759852667209, 13.357083769198743, 13.422189468736697, 13.423410996316557, 14.2241789429622]
14 : [3.2202514321744378, 5.880473664098325, 5.8804777456874975, 10.685900848563966, 10.693885840014135, 10.693919268633564, 12.317488502636113, 12.317586424666743, 13.357035966487002, 13.422036759344733, 13.422616386742618, 14.219930569002681]
15 : [3.22025143217444, 5.880473664098371, 5.880477745687513, 10.685900846411391, 10.693885839291891, 10.693919251610216, 12.31747850432522, 12.317583476150716, 13.357014277922937, 13.421959840109471, 13.422254221378731, 14.217570511002606]
16 : [3.2202514321744413, 5.880473664098278, 5.880477745687487, 10.685900845908343, 10.69388583912399, 10.69391924759071, 12.317474796203365, 12.317582549965826, 13.357004471489372, 13.421913752764338, 13.422096288196327, 14.216250738908732]
17 : [3.2202514321744395, 5.880473664098203, 5.88047774568757, 10.685900845788382, 10.693885839085182, 10.693919246614222, 12.317473357820846, 12.317582209257523, 13.356999975559484, 13.421880176271921, 13.422028662275297, 14.215430637473878]
18 : [3.2202514321744373, 5.880473664098218, 5.880477745687664, 10.68590084576264, 10.69388583907748, 10.693919246427866, 12.317472962721535, 12.317582120587106, 13.356997985426148, 13.421862754509181, 13.422004839694768, 14.215033923383757]
19 : [3.2202514321744378, 5.880473664098296, 5.880477745687385, 10.685900845757113, 10.69388583907562, 10.693919246389479, 12.317472847814312, 12.317582095739601, 13.356997105258737, 13.421854904738588, 13.421995515602898, 14.214840499150508]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.2202514321744378
5.880473664098296
5.880477745687385
10.685900845757113
10.69388583907562
10.693919246389479
12.317472847814312
12.317582095739601
13.356997105258737
13.421854904738588
13.421995515602898
14.214840499150508
[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);
[ ]:

[ ]: