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.229093263030939, 22.63644749845294, 38.759126815050365, 56.95688371923037, 64.68904675582496, 71.6151064432098, 76.89311294819163, 86.36950975289437, 90.9273026799541, 99.90890223866596, 108.0770980132416, 123.53078334034487]
1 : [3.2310126754305757, 6.086234263735, 6.406899224352031, 12.054533474487707, 12.67242102885243, 13.138599687551444, 14.490032249298535, 17.72368170680729, 18.66258439894642, 19.346120604297926, 23.999160961390903, 25.45185223956675]
2 : [3.220379348186483, 5.8866238226103, 5.90554842930596, 10.853747390754, 10.901287970279357, 11.033338300578407, 12.490608610003358, 13.713193967667966, 14.184813087541915, 15.348674544385062, 16.52982807703296, 17.964540349879922]
3 : [3.220253450726734, 5.88076870867771, 5.881908488028316, 10.709133420846037, 10.72605206598382, 10.778656219125608, 12.345960396963092, 12.760676043043883, 13.760281355763947, 14.0183260424985, 14.607718515773605, 16.59221415230251]
4 : [3.220251466900596, 5.880494279982016, 5.880556917690788, 10.689894953486759, 10.698983267278889, 10.718281105561807, 12.323621310860354, 12.515134575107359, 13.530915635996354, 13.64635249085355, 14.013638403795598, 15.844042105010148]
5 : [3.220251432780629, 5.880477114537635, 5.88048002326391, 10.686851036555717, 10.694700891024496, 10.700890407064655, 12.31910979366205, 12.416240501347415, 13.409477889851164, 13.552567887384276, 13.749216311153418, 15.246452470061332]
6 : [3.2202514321849836, 5.880473890511282, 5.8804778433188565, 10.686161008730622, 10.694032876757744, 10.695820538594074, 12.318002400651658, 12.365108987134413, 13.375722155026091, 13.488497017515277, 13.613981824870164, 14.81374636675156]
7 : [3.2202514321746247, 5.880473676484506, 5.880477751399411, 10.685969287196391, 10.69391507291702, 10.694415397845718, 12.317696219797133, 12.338610207962233, 13.365112608769744, 13.452949540289367, 13.531662729017814, 14.54811834696442]
8 : [3.2202514321744466, 5.880473664772797, 5.880477746015278, 10.685917731763409, 10.693892108137469, 10.694043457650968, 12.317604665085003, 12.326098068813366, 13.360878384956358, 13.435960447039632, 13.480566199523068, 14.398036526128198]
9 : [3.22025143217444, 5.880473664135133, 5.880477745705952, 10.68590487291429, 10.693887275313566, 10.693949433864619, 12.317575537840787, 12.32079713122898, 13.358877986663096, 13.42825810540399, 13.45152542611742, 14.315735109619283]
10 : [3.2202514321744435, 5.88047366410032, 5.880477745688505, 10.685901792193738, 10.69388617501005, 10.693926453072502, 12.317565160998509, 12.318715755177262, 13.357892966798289, 13.424793976545763, 13.436298718818202, 14.270543591547746]
11 : [3.2202514321744418, 5.880473664098408, 5.880477745687529, 10.685901066924655, 10.693885916583636, 10.693920952758111, 12.317558831675962, 12.317939834213142, 13.35741456944422, 13.423228035365167, 13.428747521079597, 14.245626509055091]
12 : [3.22025143217444, 5.880473664098318, 5.880477745687476, 10.685900897346155, 10.693885856844647, 10.693919648675092, 12.317545562899507, 12.317668127421765, 13.357188433396633, 13.422515876524006, 13.425120106915218, 14.231830211266848]
13 : [3.220251432174441, 5.880473664098311, 5.880477745687476, 10.685900857789028, 10.693885843152678, 10.69391934104149, 12.317512863593556, 12.317598526607053, 13.357083768265296, 13.422189471672299, 13.42341100106627, 14.22417892540382]
14 : [3.220251432174439, 5.880473664098303, 5.880477745687477, 10.685900848564097, 10.693885840013888, 10.693919268635298, 12.317488503472934, 12.317586424917593, 13.357035966260534, 13.422036764411985, 13.422616400627325, 14.219930733442473]
15 : [3.2202514321744387, 5.8804736640983055, 5.880477745687442, 10.685900846411457, 10.693885839291761, 10.693919251611023, 12.317478505020366, 12.317583476309355, 13.35701427777467, 13.421959844151383, 13.422254226356847, 14.217570713545529]
16 : [3.2202514321744458, 5.880473664098419, 5.880477745688002, 10.68590084590902, 10.693885839124665, 10.69391924760548, 12.317474818733748, 12.31758255587327, 13.35700445525507, 13.421913950609072, 13.422096676645124, 14.216257434695079]
17 : [3.220251432174441, 5.8804736640982584, 5.88047774568752, 10.685900845791354, 10.693885839086414, 10.693919246666695, 12.317473507525222, 12.31758224760535, 13.35700002370025, 13.421883774081689, 13.422032752017804, 14.215527678975553]
18 : [3.220251432174442, 5.880473664098246, 5.880477745687463, 10.685900845763898, 10.69388583907781, 10.693919246446244, 12.31747304417359, 12.317582141127962, 13.356998029616255, 13.421866018872809, 13.422007361730227, 14.215123082243357]
19 : [3.220251432174442, 5.8804736640983375, 5.880477745687503, 10.68590084575733, 10.693885839075572, 10.693919246393957, 12.31747288028328, 12.317582103811933, 13.356997132711847, 13.421856807774267, 13.421996783361767, 14.214898434710602]
[6]:
print ("Eigenvalues")
for lam in evals:
print (lam)
Eigenvalues
3.220251432174442
5.8804736640983375
5.880477745687503
10.68590084575733
10.693885839075572
10.693919246393957
12.31747288028328
12.317582103811933
13.356997132711847
13.421856807774267
13.421996783361767
14.214898434710602
[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);
[ ]:
[ ]: