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.229093261108188, 22.636447501186336, 38.75912680385797, 56.956883716875716, 64.68904674720346, 71.61510644501536, 76.89311292946938, 86.36950975843092, 90.92730267212008, 99.90890224404565, 108.0770980128544, 123.53078334310587]
1 : [3.231012675415519, 6.086234263883559, 6.40689922403135, 12.054533474056596, 12.67242102840479, 13.138599689594514, 14.490032245655255, 17.723681705952995, 18.662584394432145, 19.34612060356715, 23.999160962897008, 25.451852235962964]
2 : [3.2203793481863316, 5.8866238226111225, 5.905548429297897, 10.853747390873782, 10.901287970151879, 11.033338300632131, 12.490608609873444, 13.713193965881622, 14.184813087556183, 15.348674546748924, 16.529828079161593, 17.964540348789157]
3 : [3.2202534507267333, 5.880768708677749, 5.88190848802824, 10.70913342085047, 10.726052066012093, 10.778656219146773, 12.345960396958263, 12.760676042763343, 13.760281355812134, 14.018326044031317, 14.607718515333216, 16.592214150435414]
4 : [3.2202514669005926, 5.8804942799820115, 5.880556917690768, 10.689894953487256, 10.698983267288853, 10.718281105564618, 12.323621310863624, 12.515134575004648, 13.530915636424629, 13.646352490630678, 14.013638403709352, 15.844042102796697]
5 : [3.2202514327806293, 5.880477114537647, 5.8804800232639325, 10.68685103655663, 10.694700891026436, 10.700890407063051, 12.319109793664438, 12.416240501267305, 13.409477889964103, 13.552567887194364, 13.749216311140856, 15.246452468056656]
6 : [3.2202514321849796, 5.880473890511256, 5.880477843318864, 10.686161008730982, 10.694032876758136, 10.695820538592827, 12.318002400652869, 12.365108987075779, 13.375722155074552, 13.488497017389504, 13.613981824812155, 14.813746365313278]
7 : [3.220251432174621, 5.880473676484469, 5.880477751399396, 10.685969287196489, 10.69391507291711, 10.694415397845216, 12.317696219797677, 12.338610207927474, 13.3651126087979, 13.452949540219949, 13.53166272893567, 14.548118346073622]
8 : [3.220251432174443, 5.8804736647727935, 5.88047774601527, 10.685917731763448, 10.693892108137481, 10.694043457650842, 12.317604665085245, 12.32609806879582, 13.360878384973265, 13.435960447014288, 13.480566199481084, 14.398036525624832]
9 : [3.2202514321744395, 5.880473664135128, 5.880477745705958, 10.685904872915037, 10.693887275313203, 10.69394943386606, 12.317575537838524, 12.320797131217752, 13.358877986617687, 13.428258103369288, 13.451525420361559, 14.315735114758777]
10 : [3.22025143217444, 5.8804736641003, 5.880477745688498, 10.685901792194061, 10.693886175009885, 10.693926453073338, 12.317565160998505, 12.318715755005398, 13.35789296653311, 13.424793975017652, 13.436298711063287, 14.27054358949142]
11 : [3.220251432174442, 5.880473664098408, 5.880477745687526, 10.68590106692526, 10.693885916583488, 10.693920952760175, 12.317558831686782, 12.31793983449287, 13.357414569443575, 13.423228034388766, 13.428747521148582, 14.245626531348888]
12 : [3.220251432174442, 5.880473664098334, 5.880477745687479, 10.685900897346519, 10.693885856844503, 10.693919648675768, 12.317545562956864, 12.317668127682033, 13.357188433374326, 13.422515874924517, 13.425120112779986, 14.231830253684384]
13 : [3.220251432174443, 5.880473664098308, 5.880477745687497, 10.685900857789354, 10.69388584315275, 10.69391934104015, 12.317512863215237, 12.317598526267526, 13.357083769980973, 13.422189472508741, 13.423411006993891, 14.22417895061426]
14 : [3.2202514321744387, 5.880473664098297, 5.880477745687495, 10.68590084856374, 10.69388584001332, 10.693919268633385, 12.317488500544021, 12.317586423906894, 13.357035967588292, 13.42203676407834, 13.422616126088853, 14.219929872035609]
15 : [3.220251432174439, 5.880473664098381, 5.88047774568763, 10.685900846411261, 10.693885839290823, 10.693919251604408, 12.317478497452267, 12.317583474107774, 13.357014278772441, 13.421959821941563, 13.422253963162737, 14.21756905338916]
16 : [3.2202514321744418, 5.880473664098284, 5.880477745687424, 10.685900845908716, 10.693885839125034, 10.693919247607452, 12.317474828893932, 12.317582559566905, 13.357004471507743, 13.421913951137102, 13.422096883618254, 14.216258416125703]
17 : [3.2202514321744444, 5.880473664098354, 5.880477745687523, 10.685900845790648, 10.693885839086734, 10.693919246664224, 12.317473505592911, 12.317582248231663, 13.35700001789587, 13.42188346531371, 13.42203282581387, 14.215520795568146]
18 : [3.220251432174443, 5.880473664098403, 5.880477745687771, 10.685900845763504, 10.69388583907808, 10.693919246443025, 12.317473039729757, 12.31758214075261, 13.356998003494276, 13.421865673242104, 13.4220074250728, 14.215115046688888]
19 : [3.2202514321744427, 5.880473664098323, 5.880477745687494, 10.685900845757356, 10.69388583907595, 10.693919246393179, 12.317472876663103, 12.31758210310319, 13.356997117280802, 13.421856563228983, 13.421996665466686, 14.214891451190615]
[6]:
print ("Eigenvalues")
for lam in evals:
print (lam)
Eigenvalues
3.2202514321744427
5.880473664098323
5.880477745687494
10.685900845757356
10.69388583907595
10.693919246393179
12.317472876663103
12.31758210310319
13.356997117280802
13.421856563228983
13.421996665466686
14.214891451190615
[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);
[ ]:
[ ]: