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.229093262517495, 22.63644750632995, 38.75912681430485, 56.95688372933606, 64.68904675132764, 71.61510644428388, 76.89311295249415, 86.36950976165458, 90.9273026714243, 99.90890224146798, 108.07709801531121, 123.53078333910108]
1 : [3.231012675426837, 6.086234263817299, 6.406899224206793, 12.054533474166774, 12.672421029216602, 13.13859968714201, 14.490032250010394, 17.723681707619278, 18.662584394722224, 19.346120604234173, 23.99916096142672, 25.451852238521344]
2 : [3.220379348186442, 5.886623822610013, 5.905548429302344, 10.85374739076303, 10.901287970247983, 11.033338300443676, 12.490608609966753, 13.71319396711491, 14.18481308739666, 15.348674545569267, 16.529828077771203, 17.96454034926258]
3 : [3.2202534507267337, 5.880768708677767, 5.8819084880283645, 10.709133420841953, 10.72605206601449, 10.778656219069564, 12.345960396954961, 12.760676043057934, 13.76028135565477, 14.018326043262483, 14.607718515565123, 16.592214150763912]
4 : [3.2202514669005926, 5.88049427998201, 5.880556917690791, 10.689894953488857, 10.69898326728709, 10.718281105540221, 12.32362131085862, 12.515134575100829, 13.530915636187254, 13.646352490620588, 14.01363840372658, 15.844042103146585]
5 : [3.220251432780628, 5.880477114537637, 5.880480023263897, 10.686851036557366, 10.694700891026075, 10.70089040705614, 12.319109793661749, 12.416240501308316, 13.409477889893983, 13.552567887210971, 13.749216311063835, 15.246452468394242]
6 : [3.2202514321849853, 5.8804738905112695, 5.880477843318861, 10.686161008731233, 10.694032876758072, 10.695820538591027, 12.318002400651679, 12.36510898709511, 13.37572215504247, 13.488497017417007, 13.613981824733164, 14.813746365590685]
7 : [3.2202514321746256, 5.880473676484487, 5.88047775139941, 10.685969287196565, 10.693915072917099, 10.694415397844738, 12.317696219797185, 12.338610207936414, 13.365112608779517, 13.452949540237624, 13.531662728883324, 14.5481183462721]
8 : [3.220251432174443, 5.880473664772799, 5.880477746015269, 10.685917731763459, 10.693892108137478, 10.694043457650755, 12.317604665085051, 12.32609806879665, 13.360878384961556, 13.435960447028737, 13.480566199376645, 14.398036525606516]
9 : [3.22025143217444, 5.8804736641351125, 5.880477745705938, 10.685904872913605, 10.693887275314228, 10.693949433861311, 12.317575537839353, 12.320797129990506, 13.358877986679097, 13.42825810641051, 13.451525413330248, 14.315735079603325]
10 : [3.220251432174441, 5.880473664100313, 5.8804777456885, 10.685901792193395, 10.693886175010544, 10.693926453070588, 12.317565160998477, 12.31871575293806, 13.357892966813361, 13.424793975445633, 13.436298683754973, 14.27054352772191]
11 : [3.2202514321744418, 5.880473664098407, 5.880477745687548, 10.6859010669249, 10.693885916583845, 10.693920952758392, 12.31755883167803, 12.317939832988229, 13.357414569625197, 13.423228034257354, 13.428747493194535, 14.245626463249023]
12 : [3.220251432174442, 5.880473664098297, 5.880477745687462, 10.685900897346375, 10.693885856844426, 10.693919648675092, 12.3175455628595, 12.317668126916569, 13.357188433348298, 13.42251587359887, 13.425120082924108, 14.231830187836255]
13 : [3.220251432174443, 5.880473664098313, 5.880477745687493, 10.685900857787267, 10.693885843152573, 10.693919341021765, 12.31751285861197, 12.317598523120234, 13.35708376529633, 13.422189453469977, 13.423410911314376, 14.22417814124833]
14 : [3.2202514321744395, 5.880473664098298, 5.88047774568746, 10.685900848563467, 10.693885840013811, 10.693919268626493, 12.317488500521115, 12.317586424092632, 13.357035963055875, 13.422036737962623, 13.422616421966596, 14.219930270704435]
15 : [3.2202514321744427, 5.880473664098311, 5.880477745687424, 10.685900846411288, 10.693885839291744, 10.693919251608214, 12.317478504032687, 12.31758347616797, 13.357014275173075, 13.421959819853926, 13.422254262938146, 14.217570459186538]
16 : [3.220251432174442, 5.880473664098284, 5.880477745687473, 10.685900845908693, 10.693885839123201, 10.693919247592403, 12.317474798075212, 12.317582549679791, 13.35700445876812, 13.421913827120406, 13.422096388971607, 14.21625316173742]
17 : [3.2202514321744427, 5.880473664098329, 5.880477745687498, 10.685900845791089, 10.693885839085766, 10.693919246663498, 12.317473499090676, 12.317582244536204, 13.35699998781013, 13.421883729977935, 13.422032889585273, 14.215527046182116]
18 : [3.220251432174441, 5.880473664098336, 5.880477745687504, 10.685900845763745, 10.693885839077643, 10.693919246445526, 12.317473040618616, 12.317582139844776, 13.356998015851262, 13.42186600099224, 13.422007423000245, 14.215122621556766]
19 : [3.220251432174444, 5.880473664098343, 5.880477745687457, 10.685900845757393, 10.693885839075802, 10.693919246394143, 12.317472878408607, 12.31758210318456, 13.356997126620561, 13.421856778181393, 13.421996761827367, 14.214897515518203]
[6]:
print ("Eigenvalues")
for lam in evals:
print (lam)
Eigenvalues
3.220251432174444
5.880473664098343
5.880477745687457
10.685900845757393
10.693885839075802
10.693919246394143
12.317472878408607
12.31758210318456
13.356997126620561
13.421856778181393
13.421996761827367
14.214897515518203
[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);
[ ]:
[ ]: