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.229093261317333, 22.636447506193196, 38.759126805245636, 56.956883726320754, 64.68904674413668, 71.61510644246093, 76.89311292105116, 86.36950975284923, 90.92730267051611, 99.90890223761907, 108.07709801127926, 123.53078334171252]
1 : [3.2310126754184894, 6.086234263720918, 6.4068992240254445, 12.05453347262448, 12.672421028772245, 13.138599686575384, 14.490032246648457, 17.723681705551144, 18.66258439708563, 19.346120606680945, 23.999160962880506, 25.45185223882124]
2 : [3.2203793481863836, 5.886623822604202, 5.905548429300085, 10.853747390698919, 10.901287970003706, 11.033338300454293, 12.49060860983734, 13.713193968344335, 14.184813087400945, 15.348674547999487, 16.5298280783659, 17.964540348629555]
3 : [3.2202534507267324, 5.880768708677544, 5.881908488028703, 10.709133420830913, 10.726052065970897, 10.778656219140117, 12.345960396952954, 12.76067604348762, 13.760281355676272, 14.018326044055202, 14.607718516090214, 16.592214149670962]
4 : [3.2202514669005904, 5.880494279982008, 5.880556917690812, 10.689894953484401, 10.6989832672793, 10.718281105567723, 12.323621310861073, 12.515134575259319, 13.530915636420394, 13.646352490527272, 14.013638404207095, 15.844042101844824]
5 : [3.220251432780627, 5.880477114537634, 5.880480023263906, 10.68685103655619, 10.694700891024514, 10.700890407063284, 12.319109793663129, 12.416240501359443, 13.409477889977705, 13.552567887150385, 13.749216311311354, 15.246452467187106]
6 : [3.2202514321849844, 5.880473890511277, 5.8804778433188645, 10.686161008730997, 10.694032876757726, 10.69582053859247, 12.31800240065231, 12.36510898710488, 13.375722155084864, 13.488497017387017, 13.613981824805826, 14.813746364703292]
7 : [3.2202514321746207, 5.8804736764844865, 5.880477751399398, 10.685969287196537, 10.693915072917003, 10.694415397845024, 12.317696219797451, 12.3386102079349, 13.365112608803718, 13.452949540222912, 13.531662728865571, 14.548118345696844]
8 : [3.2202514321744458, 5.880473664772777, 5.880477746015279, 10.685917731763428, 10.693892108137456, 10.694043457650787, 12.317604665085163, 12.326098068793819, 13.360878384975402, 13.435960447010515, 13.480566199329216, 14.398036525292671]
9 : [3.2202514321744413, 5.880473664135123, 5.8804777457059485, 10.685904872914517, 10.69388727531334, 10.693949433864669, 12.317575537839137, 12.320797130516812, 13.358877986650185, 13.428258103489428, 13.451525408573445, 14.315735092559665]
10 : [3.2202514321744378, 5.880473664100317, 5.8804777456885216, 10.685901792193718, 10.693886175010158, 10.693926453072473, 12.317565160998594, 12.318715754177237, 13.357892966506197, 13.424793976864283, 13.436298702020714, 14.27054355842269]
11 : [3.220251432174442, 5.880473664098419, 5.880477745687536, 10.68590106692482, 10.693885916583547, 10.693920952759182, 12.317558831676244, 12.317939833747362, 13.357414569373615, 13.423228036020825, 13.428747508176308, 14.245626485671481]
12 : [3.2202514321744435, 5.880473664098315, 5.880477745687481, 10.685900897346347, 10.69388585684455, 10.69391964867538, 12.317545562889753, 12.317668127265433, 13.357188433171215, 13.422515875585256, 13.425120099177756, 14.23183020956734]
13 : [3.2202514321744413, 5.8804736640983, 5.880477745687484, 10.685900857785512, 10.69388584315246, 10.693919341025433, 12.317512860226712, 12.317598524731777, 13.357083756762826, 13.42218944831303, 13.423410894727745, 14.224177959610156]
14 : [3.220251432174438, 5.8804736640983055, 5.880477745687475, 10.685900848562582, 10.693885840013957, 10.693919268628418, 12.317488500708437, 12.317586424305912, 13.35703595725314, 13.422036744025098, 13.42261635722612, 14.219929766432028]
15 : [3.2202514321744404, 5.880473664098848, 5.880477745687826, 10.685900846411261, 10.693885839291553, 10.693919251608822, 12.317478503376465, 12.317583476031189, 13.357014272993576, 13.421959824067816, 13.422254231242153, 14.217569943480939]
16 : [3.220251432174441, 5.880473664098762, 5.880477745687895, 10.685900845908627, 10.69388583912425, 10.693919247595357, 12.317474801633939, 12.317582551279274, 13.357004460052496, 13.421913808249926, 13.422096321406416, 14.216252271587598]
17 : [3.220251432174441, 5.880473664098262, 5.880477745687491, 10.685900845791188, 10.693885839086445, 10.693919246664613, 12.317473501324532, 12.317582245933842, 13.357000039445706, 13.421883708356743, 13.42203282651158, 14.215525870005782]
18 : [3.2202514321744413, 5.88047366409767, 5.880477745687326, 10.685900845763717, 10.693885839077625, 10.693919246445631, 12.317473041897204, 12.31758214054047, 13.3569980394484, 13.421865999813802, 13.422007504536893, 14.215122574341423]
19 : [3.2202514321744418, 5.880473664098049, 5.880477745687264, 10.685900845757164, 10.693885839075898, 10.693919246393557, 12.31747287395527, 12.31758210231319, 13.356997125328832, 13.421856374476718, 13.421996513630242, 14.21488566324688]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.2202514321744418
5.880473664098049
5.880477745687264
10.685900845757164
10.693885839075898
10.693919246393557
12.31747287395527
12.31758210231319
13.356997125328832
13.421856374476718
13.421996513630242
14.21488566324688
[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);
[ ]:

[ ]: