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 = 43192
[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 : [8.213129807433281, 17.799813243476642, 37.219825630917704, 44.01930033575626, 62.182333369610014, 65.5333326176935, 73.98716308827396, 81.08231276789282, 88.06001748202364, 91.25055695394083, 113.61980592240572, 122.21673393013504]
1 : [3.2380975391467786, 5.99865063057431, 6.16139007181632, 11.485860551737803, 12.037309680452932, 14.413520502387778, 15.234814501024623, 16.218269371414898, 17.3467269521971, 19.892283929703467, 22.638813212885495, 24.807560965930858]
2 : [3.2205426885754815, 5.884876708902131, 5.888226095519312, 10.785471632227031, 10.847156783249767, 11.38198637861422, 13.343238318587616, 13.372706507212673, 14.112393286221048, 14.988639564518211, 15.344080017587014, 17.688329578353013]
3 : [3.2202954145601446, 5.880713446554476, 5.880839526591397, 10.707466526329547, 10.720069178416418, 10.944773311351312, 12.52797395467153, 13.075036551480183, 13.540697916007712, 14.030695748945433, 14.367494218474043, 15.906741590634507]
4 : [3.2202898941021916, 5.880502253291045, 5.880518097321864, 10.695114363079703, 10.699129958360583, 10.787751811973864, 12.35409837068683, 12.87105304001649, 13.445683487770072, 13.649403154679991, 14.213859960611483, 14.607665453521825]
5 : [3.220289766421209, 5.880493375254633, 5.880495963116381, 10.692776367607264, 10.694999509696416, 10.72190631751917, 12.32420257887725, 12.640658946929793, 13.417157954135064, 13.509762215105345, 13.816759278855978, 14.260863115275317]
6 : [3.2202897638070254, 5.880492978117805, 5.880494510444929, 10.691097181820915, 10.69416115365393, 10.69864359442297, 12.31888328176491, 12.460746605913451, 13.402394959764473, 13.454795875880757, 13.563293264230754, 14.232709516172516]
7 : [3.220289763759515, 5.880492958816478, 5.880494425037907, 10.687861627690577, 10.693985793780072, 10.694552792301552, 12.317899571105425, 12.372364426933725, 13.387450842658794, 13.433939412880058, 13.473271158453052, 14.223431358591203]
8 : [3.220289763758697, 5.8804929578344085, 5.880494420312166, 10.686464368105307, 10.693942498545319, 10.694085266049619, 12.317711055329555, 12.33739216776379, 13.373735427632298, 13.4266400652706, 13.441852510733863, 14.219267826322259]
9 : [3.2202897637586854, 5.880492957782831, 5.88049442005701, 10.686108698959009, 10.693925371312535, 10.694004893998665, 12.31767374021372, 12.324680164828361, 13.365194188430557, 13.423921247812626, 13.430320729564615, 14.217191023874257]
10 : [3.220289763758686, 5.880492957780075, 5.880494420043351, 10.686025216519429, 10.693919077435595, 10.693989757202642, 12.317666113412312, 12.32021020633851, 13.360871916602646, 13.42281682312806, 13.425650479541517, 14.21610359477233]
11 : [3.220289763758685, 5.880492957779919, 5.880494420042608, 10.686006011076943, 10.693917414577859, 10.69398657096486, 12.317664492529525, 12.318659205572484, 13.358849298224087, 13.422344037382345, 13.423640060505733, 14.215520285379744]
12 : [3.2202897637586836, 5.880492957779927, 5.88049442004257, 10.68600161216282, 10.6939170216158, 10.693985859212951, 12.317664119923569, 12.318123652906175, 13.357930691154058, 13.422136540650106, 13.422748133512812, 14.215203523697973]
13 : [3.2202897637586854, 5.88049295777992, 5.880494420042562, 10.686000605618666, 10.693916931077856, 10.6939856975428, 12.317664014899247, 12.317939034884924, 13.35751851674705, 13.422044241634865, 13.422347279527322, 14.215030378880067]
14 : [3.2202897637586863, 5.88049295777991, 5.880494420042553, 10.686000375053458, 10.69391691032861, 10.6939856606241, 12.317663977200707, 12.317875424508623, 13.357334571222188, 13.42200296728519, 13.42216589187846, 14.214935395967437]
15 : [3.2202897637586827, 5.8804929577800795, 5.880494420042441, 10.686000322174694, 10.693916905561904, 10.693985652172644, 12.317663963009258, 12.317853460444992, 13.357252513331781, 13.421984344896025, 13.422083665420265, 14.214883191195433]
16 : [3.220289763758683, 5.880492957780017, 5.880494420042982, 10.686000310011858, 10.693916904465882, 10.693985650231035, 12.317663957912664, 12.317845866527689, 13.357215906216771, 13.421975914350574, 13.422046272907036, 14.21485443952646]
17 : [3.220289763758679, 5.880492957780194, 5.88049442004258, 10.68600030721254, 10.69391690421292, 10.693985649784034, 12.317663956143507, 12.317843236667555, 13.357199582231589, 13.421972066676357, 13.422029278480078, 14.214838619553099]
18 : [3.2202897637586863, 5.8804929577798095, 5.880494420042483, 10.68600030651636, 10.693916904151129, 10.693985649673108, 12.317663955462683, 12.317842205368072, 13.35719113278918, 13.421970184896898, 13.4220206308138, 14.214829215771095]
19 : [3.220289763758685, 5.880492957779563, 5.880494420042597, 10.686000306398654, 10.693916904140057, 10.693985649655463, 12.31766395528269, 12.317841939006831, 13.357188135991336, 13.421969414213303, 13.422017376192283, 14.214824380618872]
[6]:
print ("Eigenvalues")
for lam in evals:
print (lam)
Eigenvalues
3.220289763758685
5.880492957779563
5.880494420042597
10.686000306398654
10.693916904140057
10.693985649655463
12.31766395528269
12.317841939006831
13.357188135991336
13.421969414213303
13.422017376192283
14.214824380618872
[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);
[ ]:
[ ]: