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.229093263839168, 22.63644748734098, 38.75912681224863, 56.95688374137755, 64.68904675291994, 71.61510643787639, 76.89311293503295, 86.36950976023813, 90.92730267735544, 99.90890224530513, 108.07709801512493, 123.53078334590454]
1 : [3.2310126754251636, 6.086234263835054, 6.4068992240998766, 12.054533476035267, 12.672421029155068, 13.138599689456745, 14.490032245669784, 17.723681707537096, 18.662584396162483, 19.346120603096477, 23.999160957391624, 25.451852237988717]
2 : [3.2203793481864307, 5.886623822613393, 5.905548429292531, 10.85374739097596, 10.901287970152401, 11.03333830047825, 12.49060860986689, 13.713193966350621, 14.184813087365864, 15.348674545295893, 16.52982807367577, 17.96454034913878]
3 : [3.2202534507267297, 5.880768708677716, 5.881908488027533, 10.709133420836551, 10.726052066022191, 10.778656219029383, 12.345960396956688, 12.760676042495307, 13.760281355531383, 14.018326042547244, 14.60771851526271, 16.592214150535586]
4 : [3.22025146690059, 5.880494279982008, 5.88055691769071, 10.68989495348563, 10.698983267286664, 10.718281105519486, 12.32362131086342, 12.515134574821568, 13.530915635949876, 13.646352490633857, 14.013638403727809, 15.844042102996415]
5 : [3.220251432780632, 5.880477114537658, 5.8804800232639245, 10.686851036555485, 10.694700891025787, 10.700890407048323, 12.319109793664323, 12.416240501169606, 13.409477889837655, 13.552567887216725, 13.749216311070654, 15.246452468316745]
6 : [3.2202514321849853, 5.880473890511272, 5.88047784331886, 10.68616100873034, 10.694032876757984, 10.695820538588658, 12.318002400652771, 12.365108987029187, 13.375722155022899, 13.488497017425269, 13.613981824731987, 14.813746365562965]
7 : [3.2202514321746234, 5.880473676484505, 5.88047775139941, 10.685969287196269, 10.693915072917063, 10.694415397844141, 12.317696219797568, 12.338610207907692, 13.365112608768202, 13.452949540247504, 13.531662728879008, 14.548118346245296]
8 : [3.2202514321744453, 5.880473664772779, 5.880477746015295, 10.685917731763348, 10.693892108137481, 10.694043457650649, 12.317604665085154, 12.326098068777435, 13.360878384954756, 13.435960447033882, 13.48056619926108, 14.398036525432218]
9 : [3.2202514321744387, 5.880473664135134, 5.880477745705939, 10.685904872914435, 10.69388727531388, 10.693949433864265, 12.317575537840076, 12.320797131251927, 13.35887798665495, 13.428258105681124, 13.451525428935355, 14.315735113454076]
10 : [3.2202514321744387, 5.880473664100308, 5.880477745688505, 10.685901792193652, 10.693886175010334, 10.693926453072285, 12.31756516099848, 12.318715754257394, 13.357892966665661, 13.424793978712652, 13.436298712520644, 14.270543566597665]
11 : [3.2202514321744395, 5.880473664098412, 5.880477745687531, 10.68590106692489, 10.69388591658381, 10.6939209527587, 12.317558831676097, 12.317939833667664, 13.357414569440259, 13.423228037227132, 13.428747514879385, 14.245626488326304]
12 : [3.2202514321744413, 5.880473664098318, 5.880477745687488, 10.685900897346409, 10.693885856844972, 10.693919648675177, 12.31754556288799, 12.317668127252396, 13.357188433243488, 13.422515877049543, 13.425120112331003, 14.231830221040372]
13 : [3.22025143217444, 5.880473664098287, 5.880477745687478, 10.685900857789434, 10.69388584315232, 10.693919341041969, 12.31751286382646, 12.317598526763536, 13.357083769189531, 13.422189467062681, 13.423410995411551, 14.22417900304734]
14 : [3.2202514321744427, 5.880473664098283, 5.880477745687463, 10.685900848564332, 10.693885840013285, 10.693919268635273, 12.317488504552697, 12.317586425314897, 13.357035968538563, 13.422036759182065, 13.422616416793721, 14.219931017016473]
15 : [3.220251432174442, 5.8804736640982345, 5.880477745687334, 10.685900846411718, 10.693885839291523, 10.693919251609692, 12.317478504462308, 12.31758347614738, 13.357014279956338, 13.42195983994086, 13.422254204621659, 14.217570733321136]
16 : [3.2202514321744387, 5.880473664098343, 5.880477745687579, 10.685900845908757, 10.693885839124453, 10.693919247606779, 12.317474824563325, 12.317582557464545, 13.357004434577604, 13.421913967180851, 13.422096835300861, 14.216258885530202]
17 : [3.220251432174442, 5.8804736640981545, 5.880477745687476, 10.685900845791197, 10.69388583908663, 10.693919246667583, 12.317473512411278, 12.317582248908558, 13.3570000075939, 13.421883849980667, 13.422033034211175, 14.215530255713718]
18 : [3.2202514321744418, 5.88047366409812, 5.8804777456874575, 10.685900845763504, 10.69388583907758, 10.693919246446312, 12.317473045419717, 12.31758214153424, 13.356998018332423, 13.421866000560307, 13.422007367513322, 14.215122771367394]
19 : [3.22025143217444, 5.880473664098318, 5.880477745687349, 10.685900845757136, 10.693885839075863, 10.693919246393712, 12.317472876223109, 12.317582103008208, 13.356997105435882, 13.421856387150116, 13.421996445997337, 14.21488654749208]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.22025143217444
5.880473664098318
5.880477745687349
10.685900845757136
10.693885839075863
10.693919246393712
12.317472876223109
12.317582103008208
13.356997105435882
13.421856387150116
13.421996445997337
14.21488654749208
[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);
[ ]:

[ ]: