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.22909326277295, 22.63644750133328, 38.75912681486428, 56.95688372426084, 64.68904673845074, 71.61510643501143, 76.89311293799405, 86.36950975847232, 90.92730268043374, 99.90890224262897, 108.07709801930287, 123.53078333949155]
1 : [3.231012675419242, 6.086234263860723, 6.406899224299488, 12.054533473324685, 12.67242102950673, 13.138599689146584, 14.49003224695817, 17.723681707037805, 18.662584394441726, 19.346120605246114, 23.999160957994743, 25.451852240204566]
2 : [3.2203793481863747, 5.886623822612182, 5.905548429299338, 10.85374739067019, 10.901287970044736, 11.033338300818537, 12.49060860988545, 13.71319396688432, 14.184813087402851, 15.348674545145036, 16.529828074305065, 17.96454034897031]
3 : [3.2202534507267337, 5.880768708677791, 5.881908488028093, 10.709133420837478, 10.726052065943975, 10.77865621913533, 12.34596039695859, 12.76067604283116, 13.760281355679535, 14.018326042449164, 14.607718515212229, 16.592214151006782]
4 : [3.220251466900596, 5.880494279982014, 5.880556917690771, 10.689894953487464, 10.698983267272087, 10.718281105553489, 12.32362131086424, 12.515134575011151, 13.530915636059797, 13.646352490567041, 14.01363840365514, 15.844042103684131]
5 : [3.220251432780627, 5.880477114537636, 5.880480023263922, 10.686851036557078, 10.694700891023244, 10.700890407058774, 12.319109793664557, 12.41624050128127, 13.409477889886984, 13.552567887203788, 13.749216311043782, 15.246452468972288]
6 : [3.220251432184982, 5.880473890511281, 5.8804778433188725, 10.686161008731277, 10.694032876757507, 10.695820538591585, 12.31800240065281, 12.365108987090558, 13.375722155058876, 13.488497017412326, 13.613981824745187, 14.813746366027013]
7 : [3.2202514321746216, 5.8804736764845, 5.880477751399393, 10.685969287196617, 10.69391507291697, 10.69441539784487, 12.317696219797622, 12.338610207937627, 13.365112608793881, 13.452949540235789, 13.53166272891588, 14.548118346549218]
8 : [3.220251432174441, 5.8804736647728095, 5.880477746015279, 10.685917731763475, 10.69389210813746, 10.694043457650809, 12.31760466508517, 12.326098068794662, 13.360878384970704, 13.435960447023033, 13.480566199328521, 14.39803652567749]
9 : [3.2202514321744413, 5.880473664135113, 5.88047774570593, 10.685904872914527, 10.693887275314736, 10.69394943386357, 12.317575537839659, 12.3207971307262, 13.358877986673944, 13.42825810549477, 13.451525421387064, 14.315735100087968]
10 : [3.2202514321744427, 5.880473664100312, 5.880477745688503, 10.685901792193793, 10.693886175010588, 10.69392645307193, 12.317565160998551, 12.318715754477678, 13.357892966624219, 13.424793977124276, 13.436298710205435, 14.270543569089048]
11 : [3.2202514321744435, 5.880473664098434, 5.880477745687538, 10.685901066925027, 10.693885916583715, 10.693920952759237, 12.31755883168034, 12.317939834017203, 13.357414569481238, 13.423228035833215, 13.428747515755472, 14.245626503056767]
12 : [3.2202514321744418, 5.880473664098332, 5.88047774568747, 10.685900897346398, 10.693885856844753, 10.693919648675678, 12.3175455628975, 12.317668127279468, 13.35718843333748, 13.42251587531427, 13.425120093285877, 14.231830203136765]
13 : [3.2202514321744387, 5.8804736640983055, 5.880477745687512, 10.685900857787777, 10.693885843152792, 10.69391934102894, 12.317512858674727, 12.317598522817882, 13.357083768264408, 13.422189469152414, 13.423410829848535, 14.224178080823085]
14 : [3.220251432174442, 5.880473664098312, 5.880477745687467, 10.685900848563758, 10.693885840012678, 10.693919268627582, 12.317488500755067, 12.31758642403266, 13.357035966415577, 13.422036755907506, 13.42261646595886, 14.219930463562509]
15 : [3.220251432174442, 5.8804736640982735, 5.880477745687422, 10.685900846411387, 10.693885839291438, 10.693919251608108, 12.31747850373043, 12.31758347602676, 13.35701427772075, 13.42195983130787, 13.42225431995139, 14.217570593654886]
16 : [3.2202514321744435, 5.880473664098408, 5.880477745687412, 10.685900845908627, 10.693885839123041, 10.693919247591666, 12.317474792210232, 12.317582547632663, 13.357004463053881, 13.421913829543168, 13.42209619701162, 14.216252118441451]
17 : [3.220251432174446, 5.880473664098097, 5.8804777456873625, 10.685900845790187, 10.693885839083602, 10.693919246662434, 12.317473496520437, 12.317582243944509, 13.356999926955767, 13.421883467745221, 13.422032513251542, 14.215520085591166]
18 : [3.220251432174446, 5.880473664098437, 5.880477745687466, 10.6859008457633, 10.69388583907719, 10.69391924644506, 12.317473039438843, 12.317582139634872, 13.356997968956499, 13.421865756343248, 13.42200715502037, 14.215116603756288]
19 : [3.2202514321744444, 5.880473664098404, 5.880477745687549, 10.68590084575732, 10.693885839075671, 10.69391924639401, 12.317472877739032, 12.31758210306491, 13.356997101270382, 13.421856575487604, 13.421996434035723, 14.214891989385505]
[6]:
print ("Eigenvalues")
for lam in evals:
print (lam)
Eigenvalues
3.2202514321744444
5.880473664098404
5.880477745687549
10.68590084575732
10.693885839075671
10.69391924639401
12.317472877739032
12.31758210306491
13.356997101270382
13.421856575487604
13.421996434035723
14.214891989385505
[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);
[ ]:
[ ]: