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.229093265431429, 22.63644748624693, 38.75912681277979, 56.956883745494416, 64.68904676538875, 71.61510643119607, 76.89311297596701, 86.36950974989206, 90.9273026770357, 99.90890224466601, 108.07709801676015, 123.5307833513511]
1 : [3.231012675445081, 6.086234263661832, 6.406899224416665, 12.054533476679154, 12.672421029633274, 13.138599687420221, 14.490032249945607, 17.723681705532353, 18.662584401443475, 19.346120603107337, 23.99916095815965, 25.451852241235475]
2 : [3.220379348186626, 5.886623822612569, 5.905548429306398, 10.853747390884129, 10.901287970308402, 11.033338300547909, 12.490608610109934, 13.713193967334522, 14.184813087325994, 15.348674545895319, 16.529828072175835, 17.96454034900095]
3 : [3.2202534507267364, 5.8807687086777385, 5.881908488028335, 10.709133420833043, 10.726052066023941, 10.778656219072031, 12.345960397006216, 12.760676042676552, 13.760281355477884, 14.018326042073053, 14.607718515886157, 16.5922141504765]
4 : [3.2202514669005926, 5.880494279981991, 5.8805569176907655, 10.689894953489647, 10.698983267283085, 10.718281105536521, 12.323621310879329, 12.515134574884039, 13.530915635860442, 13.646352490594621, 14.01363840424707, 15.844042102896635]
5 : [3.2202514327806293, 5.880477114537658, 5.880480023263933, 10.68685103655748, 10.694700891025052, 10.700890407053741, 12.319109793670023, 12.416240501203152, 13.409477889843542, 13.552567887214414, 13.749216311354122, 15.246452468186863]
6 : [3.2202514321849827, 5.880473890511284, 5.880477843318871, 10.686161008731123, 10.69403287675782, 10.695820538590128, 12.318002400654754, 12.365108987046254, 13.375722155038533, 13.488497017437435, 13.613981824869418, 14.81374636545481]
7 : [3.2202514321746234, 5.88047367648447, 5.880477751399407, 10.685969287196505, 10.693915072917019, 10.69441539784448, 12.317696219798247, 12.338610207915378, 13.365112608777826, 13.452949540256368, 13.53166272893751, 14.548118346172409]
8 : [3.220251432174443, 5.880473664772747, 5.88047774601528, 10.685917731763432, 10.69389210813747, 10.694043457650679, 12.31760466508539, 12.326098068784827, 13.360878384959607, 13.435960447019513, 13.480566199302642, 14.398036525501098]
9 : [3.22025143217444, 5.8804736641351045, 5.88047774570594, 10.685904872913536, 10.693887275314228, 10.69394943386198, 12.31757553784004, 12.32079712998159, 13.358877986661598, 13.428258105641804, 13.451525409891975, 14.315735079880326]
10 : [3.2202514321744418, 5.880473664100323, 5.880477745688497, 10.685901792193222, 10.69388617501039, 10.69392645307094, 12.31756516099854, 12.318715753593949, 13.357892966507858, 13.424793976987855, 13.436298692054242, 14.270543535803954]
11 : [3.220251432174441, 5.880473664098408, 5.880477745687533, 10.6859010669246, 10.69388591658373, 10.693920952758132, 12.317558831671604, 12.31793983347186, 13.357414569279625, 13.423228035817806, 13.42874750363161, 14.245626469825607]
12 : [3.2202514321744435, 5.880473664098331, 5.880477745687474, 10.685900897346423, 10.693885856844636, 10.693919648674685, 12.317545562884426, 12.317668127233814, 13.357188432850016, 13.422515873553676, 13.425120100750757, 14.23183021812065]
13 : [3.2202514321744413, 5.880473664098321, 5.880477745687473, 10.68590085778774, 10.69388584315191, 10.693919341037192, 12.3175128607869, 12.31759852459238, 13.35708376819319, 13.422189470814196, 13.423410794487584, 14.224178202833404]
14 : [3.2202514321744395, 5.880473664098315, 5.88047774568749, 10.685900848563557, 10.693885840014058, 10.693919268633131, 12.317488500605442, 12.31758642397424, 13.357035965832342, 13.422036765011915, 13.42261622114451, 14.21992990523892]
15 : [3.2202514321744427, 5.8804736640983615, 5.880477745687601, 10.685900846411348, 10.693885839291866, 10.693919251609763, 12.317478502582299, 12.31758347560982, 13.357014277772109, 13.421959839617069, 13.42225406842197, 14.217569855132277]
16 : [3.2202514321744413, 5.88047366409845, 5.8804777456875765, 10.685900845908684, 10.693885839124874, 10.693919247605255, 12.317474822406098, 12.317582557342611, 13.357004467573322, 13.42191393351635, 13.4220967524678, 14.216257231647642]
17 : [3.2202514321744427, 5.88047366409852, 5.880477745687003, 10.685900845791268, 10.693885839086409, 10.693919246667255, 12.317473511529839, 12.317582248960647, 13.357000038520875, 13.42188380831262, 13.42203294715666, 14.215528784093507]
18 : [3.220251432174444, 5.880473664098393, 5.880477745687501, 10.685900845763461, 10.693885839077812, 10.693919246443745, 12.317473038821388, 12.317582140398233, 13.356997996798407, 13.421865604271819, 13.422007364110339, 14.215113136174029]
19 : [3.2202514321744435, 5.880473664098375, 5.8804777456874975, 10.685900845757148, 10.693885839075774, 10.693919246392788, 12.317472877130145, 12.31758210333777, 13.356997103836651, 13.42185647247887, 13.421996705054951, 14.214889050034586]
[6]:
print ("Eigenvalues")
for lam in evals:
print (lam)
Eigenvalues
3.2202514321744435
5.880473664098375
5.8804777456874975
10.685900845757148
10.693885839075774
10.693919246392788
12.317472877130145
12.31758210333777
13.356997103836651
13.42185647247887
13.421996705054951
14.214889050034586
[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);
[ ]:
[ ]: