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.229093262639734, 22.63644749127682, 38.75912680969035, 56.95688373190037, 64.68904676258606, 71.61510644722989, 76.89311293509657, 86.36950974947118, 90.92730267498584, 99.90890224396665, 108.07709801705813, 123.53078334498154]
1 : [3.2310126754253767, 6.086234263854522, 6.406899224130492, 12.054533475711166, 12.672421029017055, 13.138599687822394, 14.490032246951575, 17.72368170638431, 18.662584398567805, 19.346120603004405, 23.99916095969146, 25.45185223910228]
2 : [3.2203793481864165, 5.886623822614714, 5.905548429291377, 10.853747390925164, 10.901287970076371, 11.033338300578984, 12.490608609971908, 13.713193966933975, 14.184813087324057, 15.348674545161035, 16.529828074775843, 17.964540348766697]
3 : [3.220253450726734, 5.880768708677871, 5.881908488027514, 10.7091334208336, 10.726052065994445, 10.778656219098288, 12.345960396979617, 12.760676042735712, 13.760281355535879, 14.018326042394563, 14.607718515427896, 16.59221415122042]
4 : [3.220251466900596, 5.880494279981991, 5.880556917690728, 10.689894953485389, 10.698983267280838, 10.718281105547245, 12.323621310870394, 12.515134574954242, 13.530915635973292, 13.646352490575143, 14.013638403835918, 15.844042103984052]
5 : [3.2202514327806275, 5.880477114537649, 5.880480023263933, 10.686851036556238, 10.694700891024691, 10.700890407058054, 12.319109793666696, 12.416240501254112, 13.409477889863263, 13.552567887187044, 13.749216311188183, 15.246452469239573]
6 : [3.220251432184986, 5.880473890511279, 5.880477843318857, 10.68616100873091, 10.694032876757761, 10.695820538591677, 12.31800240065361, 12.365108987080175, 13.375722155043858, 13.488497017408182, 13.613981824846599, 14.813746366217728]
7 : [3.220251432174625, 5.8804736764844785, 5.880477751399398, 10.685969287196492, 10.69391507291702, 10.694415397844978, 12.317696219797892, 12.338610207934481, 13.365112608782235, 13.45294954023641, 13.531662728976242, 14.54811834666346]
8 : [3.2202514321744413, 5.880473664772785, 5.880477746015268, 10.68591773176342, 10.69389210813749, 10.694043457650775, 12.317604665085282, 12.326098068801846, 13.36087838496241, 13.435960447023636, 13.480566199535543, 14.39803652598774]
9 : [3.2202514321744404, 5.8804736641351125, 5.8804777457059405, 10.685904872915948, 10.69388727531451, 10.693949433867068, 12.317575537837463, 12.320797131622296, 13.358877986626178, 13.42825810327373, 13.451525429180494, 14.31573512972308]
10 : [3.2202514321744418, 5.88047366410033, 5.880477745688505, 10.685901792194446, 10.693886175010501, 10.693926453073903, 12.31756516099863, 12.318715754553233, 13.357892966533829, 13.424793975184102, 13.436298708721539, 14.270543582725102]
11 : [3.2202514321744435, 5.880473664098412, 5.880477745687519, 10.685901066924883, 10.693885916583664, 10.6939209527589, 12.317558831676221, 12.317939833676936, 13.357414569316814, 13.423228034605524, 13.428747506553261, 14.245626488552775]
12 : [3.220251432174442, 5.880473664098311, 5.8804777456874735, 10.685900897346302, 10.693885856844723, 10.693919648675342, 12.317545562872342, 12.317668127135626, 13.357188433207686, 13.422515875633044, 13.425120095365767, 14.231830196272416]
13 : [3.2202514321744404, 5.880473664098323, 5.880477745687439, 10.685900857787303, 10.693885843152662, 10.693919341028389, 12.317512857456222, 12.31759852186882, 13.357083768513062, 13.422189472204007, 13.42341071875816, 14.224177742889045]
14 : [3.220251432174441, 5.8804736640982975, 5.88047774568747, 10.685900848562522, 10.693885840013513, 10.69391926862793, 12.317488496368826, 12.317586422655875, 13.357035964368777, 13.422036769589749, 13.422616027759755, 14.219928666143364]
15 : [3.220251432174443, 5.88047366409832, 5.880477745687346, 10.68590084641085, 10.69388583929154, 10.693919251608465, 12.317478501257403, 12.317583475538015, 13.35701427530452, 13.421959830455794, 13.422254018494982, 14.217568817128404]
16 : [3.220251432174445, 5.880473664098342, 5.880477745687527, 10.685900845908465, 10.693885839122162, 10.693919247589925, 12.317474778078992, 12.317582542356366, 13.357004416266545, 13.421913820210936, 13.42209559044494, 14.216249506572243]
17 : [3.220251432174444, 5.880473664098397, 5.880477745687628, 10.685900845791137, 10.693885839086168, 10.693919246662853, 12.317473488719916, 12.31758224127598, 13.357000007340833, 13.421883660865403, 13.422032400404287, 14.215523273103647]
18 : [3.220251432174438, 5.880473664098328, 5.8804777456874575, 10.685900845762262, 10.693885839077346, 10.693919246434236, 12.317472990876134, 12.317582128309319, 13.356997856023877, 13.421863063727823, 13.422003262223054, 14.215043833089307]
19 : [3.220251432174443, 5.880473664098522, 5.880477745687441, 10.68590084575685, 10.693885839075666, 10.693919246388853, 12.317472850014914, 12.317582096857757, 13.356996983100466, 13.421854258033463, 13.421994109273728, 14.214822480127104]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.220251432174443
5.880473664098522
5.880477745687441
10.68590084575685
10.693885839075666
10.693919246388853
12.317472850014914
12.317582096857757
13.356996983100466
13.421854258033463
13.421994109273728
14.214822480127104
[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);
[ ]:

[ ]: