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.229093261725584, 22.636447502925527, 38.75912680847177, 56.95688372411647, 64.68904674012049, 71.61510643928402, 76.89311293004317, 86.3695097610718, 90.92730267780436, 99.90890224461792, 108.07709801776006, 123.5307833450522]
1 : [3.231012675420079, 6.086234263775162, 6.406899223929082, 12.05453347271104, 12.67242102913964, 13.138599688225048, 14.490032248047317, 17.723681706037084, 18.662584398394397, 19.346120604880635, 23.99916096028202, 25.451852239845703]
2 : [3.2203793481863365, 5.886623822608356, 5.90554842928832, 10.853747390668602, 10.901287970050628, 11.03333830063201, 12.49060860998015, 13.713193967571486, 14.184813087398146, 15.348674546398728, 16.529828076134436, 17.964540348913886]
3 : [3.2202534507267315, 5.8807687086776514, 5.8819084880277215, 10.709133420842376, 10.726052065944689, 10.778656219121268, 12.345960396977086, 12.760676043064276, 13.760281355632227, 14.018326043316856, 14.607718515424581, 16.592214150182784]
4 : [3.2202514669005966, 5.880494279982014, 5.880556917690755, 10.689894953485542, 10.698983267276013, 10.718281105552569, 12.32362131086821, 12.515134575082364, 13.530915636250663, 13.646352490534861, 14.01363840382269, 15.844042102425819]
5 : [3.2202514327806298, 5.8804771145376495, 5.880480023263927, 10.686851036556519, 10.694700891024258, 10.700890407058456, 12.319109793665637, 12.416240501290643, 13.409477889939762, 13.552567887149255, 13.749216311133909, 15.246452467706835]
6 : [3.2202514321849836, 5.88047389051128, 5.880477843318865, 10.686161008731133, 10.694032876757673, 10.69582053859132, 12.318002400653185, 12.365108987080763, 13.375722155074284, 13.48849701737984, 13.613981824746965, 14.813746365082354]
7 : [3.220251432174623, 5.880473676484496, 5.880477751399412, 10.68596928719655, 10.693915072917022, 10.69441539784476, 12.317696219797758, 12.338610207927344, 13.36511260879855, 13.452949540217714, 13.531662728870074, 14.548118345953416]
8 : [3.220251432174445, 5.880473664772785, 5.880477746015278, 10.685917731763464, 10.693892108137456, 10.694043457650766, 12.317604665085256, 12.326098068791534, 13.360878384972052, 13.43596044701871, 13.480566199349814, 14.398036525408747]
9 : [3.2202514321744418, 5.8804736641350965, 5.880477745705944, 10.68590487291474, 10.69388727531493, 10.693949433864342, 12.317575537839279, 12.320797131162985, 13.358877986675687, 13.428258105314796, 13.451525429580274, 14.315735116336713]
10 : [3.2202514321744378, 5.880473664100327, 5.88047774568851, 10.685901792193908, 10.693886175010707, 10.693926453072349, 12.317565160998639, 12.318715754737125, 13.357892966577941, 13.424793976616753, 13.436298715828007, 14.270543584112927]
11 : [3.2202514321744418, 5.880473664098412, 5.880477745687527, 10.685901066925169, 10.693885916583886, 10.69392095275922, 12.317558831681591, 12.317939834050062, 13.357414569513582, 13.423228035489313, 13.428747518451, 14.245626513127972]
12 : [3.220251432174443, 5.880473664098317, 5.880477745687471, 10.685900897346617, 10.693885856844695, 10.693919648675141, 12.317545562923646, 12.317668127465119, 13.357188433164879, 13.422515874263855, 13.425120111554126, 14.231830248353468]
13 : [3.22025143217444, 5.8804736640983055, 5.880477745687484, 10.685900857788216, 10.69388584315234, 10.693919341037757, 12.317512861498601, 12.31759852505846, 13.357083768438594, 13.422189471509352, 13.423410866967071, 14.224178450198309]
14 : [3.22025143217444, 5.880473664098295, 5.880477745687465, 10.68590084856367, 10.693885840013557, 10.693919268633502, 12.31748850182339, 12.317586424415154, 13.357035966475221, 13.42203676769861, 13.422616314446309, 14.219930233229169]
15 : [3.220251432174443, 5.880473664098323, 5.880477745687433, 10.685900846411375, 10.693885839291747, 10.693919251610435, 12.317478504355183, 12.317583476211048, 13.357014278271668, 13.421959842462993, 13.422254170312812, 14.21757028538443]
16 : [3.220251432174442, 5.88047366409832, 5.880477745687447, 10.68590084590865, 10.693885839124986, 10.693919247605269, 12.317474825012548, 12.317582558428088, 13.357004470011526, 13.421913924833758, 13.422096859240346, 14.216257589555047]
17 : [3.2202514321744444, 5.880473664098573, 5.880477745687517, 10.68590084579134, 10.693885839086757, 10.69391924666673, 12.317473511823845, 12.317582249295718, 13.357000037907751, 13.421883771989068, 13.422033154593608, 14.215528503518573]
18 : [3.2202514321744444, 5.880473664098204, 5.880477745687227, 10.68590084576382, 10.693885839077803, 10.693919246446006, 12.317473045042885, 12.317582141564964, 13.356998039773659, 13.421866005849681, 13.422007505621186, 14.215122886044039]
19 : [3.2202514321744413, 5.880473664098503, 5.880477745687307, 10.685900845757331, 10.693885839075907, 10.693919246394394, 12.31747288023969, 12.317582103918982, 13.35699713520967, 13.421856754633277, 13.421996816625574, 14.214896993482704]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.2202514321744413
5.880473664098503
5.880477745687307
10.685900845757331
10.693885839075907
10.693919246394394
12.31747288023969
12.317582103918982
13.35699713520967
13.421856754633277
13.421996816625574
14.214896993482704
[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);
[ ]:

[ ]: