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.229093260305384, 22.63644749635922, 38.75912681396046, 56.95688371672044, 64.68904673674517, 71.61510644762376, 76.8931129371339, 86.36950975729417, 90.92730267500315, 99.90890224072902, 108.07709801365425, 123.53078333581357]
1 : [3.2310126754103723, 6.0862342637059985, 6.406899224427608, 12.054533473484705, 12.672421028702022, 13.138599689003435, 14.490032245598229, 17.72368170592679, 18.66258439425658, 19.346120605233015, 23.999160962054304, 25.45185223698652]
2 : [3.220379348186308, 5.8866238226053715, 5.905548429314656, 10.853747390868104, 10.901287969904597, 11.033338300547996, 12.490608609758333, 13.713193966565187, 14.184813087365871, 15.348674546152195, 16.529828076221374, 17.964540349587633]
3 : [3.220253450726728, 5.880768708677558, 5.881908488028968, 10.709133420829904, 10.726052065980793, 10.778656219076055, 12.345960396925788, 12.76067604270894, 13.760281355610768, 14.01832604309339, 14.607718515556686, 16.59221415167846]
4 : [3.2202514669005953, 5.8804942799819795, 5.8805569176908215, 10.68989495348226, 10.698983267282767, 10.718281105538821, 12.323621310851681, 12.515134574905602, 13.530915636067942, 13.646352490683286, 14.013638403874413, 15.844042104361721]
5 : [3.2202514327806293, 5.880477114537628, 5.880480023263912, 10.686851036555282, 10.694700891025441, 10.700890407055272, 12.319109793659777, 12.41624050121284, 13.40947788985575, 13.552567887271534, 13.749216311178595, 15.246452469560348]
6 : [3.220251432184984, 5.880473890511292, 5.880477843318868, 10.68616100873064, 10.694032876757943, 10.695820538590858, 12.318002400651105, 12.36510898705756, 13.375722155026644, 13.488497017459016, 13.61398182483936, 14.813746366447656]
7 : [3.2202514321746216, 5.880473676484483, 5.880477751399416, 10.685969287196412, 10.693915072917068, 10.694415397844779, 12.317696219797005, 12.33861020792461, 13.365112608771454, 13.452949540261022, 13.531662728963138, 14.548118346782127]
8 : [3.220251432174442, 5.880473664772753, 5.880477746015289, 10.685917731763421, 10.693892108137485, 10.694043457650716, 12.317604665084946, 12.32609806879304, 13.360878384958086, 13.435960447020808, 13.480566199385025, 14.398036525920471]
9 : [3.2202514321744427, 5.880473664135114, 5.880477745705945, 10.68590487291481, 10.69388727531386, 10.693949433865491, 12.317575537839925, 12.320797131385484, 13.358877986677411, 13.428258104544035, 13.451525426714483, 14.315735116049805]
10 : [3.2202514321744404, 5.880473664100326, 5.88047774568849, 10.68590179219396, 10.693886175010316, 10.693926453073015, 12.317565160998564, 12.318715755032123, 13.357892966599684, 13.424793976650633, 13.436298719583588, 14.270543590988797]
11 : [3.2202514321744427, 5.880473664098406, 5.880477745687546, 10.685901066924973, 10.693885916583774, 10.693920952758557, 12.317558831677628, 12.317939834007737, 13.357414569390365, 13.423228035661285, 13.428747518162169, 14.245626504020048]
12 : [3.2202514321744395, 5.880473664098283, 5.880477745687461, 10.685900897346402, 10.693885856844906, 10.693919648675138, 12.317545562890814, 12.317668127295601, 13.357188433244811, 13.422515875622347, 13.425120105470167, 14.2318302163783]
13 : [3.2202514321744418, 5.880473664098324, 5.880477745687484, 10.685900857789427, 10.69388584315201, 10.693919341041646, 12.317512863853956, 12.317598526802549, 13.357083769936821, 13.422189468127094, 13.423410995589553, 14.224178992482965]
14 : [3.2202514321744373, 5.880473664098292, 5.880477745687452, 10.685900848563557, 10.69388584001285, 10.693919268634241, 12.317488502366167, 12.31758642469448, 13.35703596365735, 13.422036761672008, 13.422616218868201, 14.219930145383374]
15 : [3.2202514321744387, 5.8804736640983375, 5.880477745687548, 10.68590084641123, 10.69388583929139, 10.693919251610453, 12.3174785051684, 12.317583476560685, 13.35701427563945, 13.421959840064481, 13.422254161361467, 14.217570301077572]
16 : [3.22025143217444, 5.880473664098275, 5.8804777456873545, 10.685900845908066, 10.693885839120426, 10.693919247565587, 12.31747472708531, 12.317582527697715, 13.357004453999672, 13.421913517053433, 13.42209430990009, 14.216235605226204]
17 : [3.220251432174439, 5.880473664098304, 5.880477745687099, 10.685900845790409, 10.693885839085208, 10.693919246645194, 12.317473444537082, 12.317582229993869, 13.35699999670657, 13.421882719554704, 13.422031501512217, 14.215498380667068]
18 : [3.220251432174443, 5.880473664098921, 5.880477745687719, 10.685900845763433, 10.693885839077469, 10.693919246438634, 12.317473014571293, 12.317582133057293, 13.356997979148902, 13.421865162519405, 13.422006343479222, 14.215099719414894]
19 : [3.220251432174439, 5.88047366409921, 5.880477745687379, 10.68590084575696, 10.693885839075506, 10.69391924639091, 12.31747285576686, 12.317582097497924, 13.356997080655775, 13.421855289196541, 13.42199494462754, 14.214852989878487]
[6]:
print ("Eigenvalues")
for lam in evals:
print (lam)
Eigenvalues
3.220251432174439
5.88047366409921
5.880477745687379
10.68590084575696
10.693885839075506
10.69391924639091
12.31747285576686
12.317582097497924
13.356997080655775
13.421855289196541
13.42199494462754
14.214852989878487
[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);
[ ]:
[ ]: