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.2290932600982085, 22.63644750064766, 38.759126806003586, 56.956883729826316, 64.68904673929455, 71.61510645149218, 76.8931129370363, 86.36950975400062, 90.92730267092381, 99.90890223731255, 108.07709800742664, 123.53078333985856]
1 : [3.231012675414463, 6.086234263779472, 6.406899224125987, 12.054533473374166, 12.672421028766514, 13.138599688224003, 14.490032245492651, 17.723681705572243, 18.66258439550712, 19.346120605408707, 23.9991609574633, 25.45185223934599]
2 : [3.2203793481863583, 5.886623822607993, 5.905548429301577, 10.853747390888794, 10.901287969810129, 11.033338300499802, 12.490608609756432, 13.713193967326182, 14.184813087236222, 15.348674547601398, 16.529828076899612, 17.964540350183988]
3 : [3.220253450726735, 5.880768708677667, 5.881908488028644, 10.709133420835416, 10.726052065965542, 10.778656219096169, 12.345960396947037, 12.760676042941226, 13.760281355627269, 14.018326044268433, 14.607718515444464, 16.592214151352163]
4 : [3.220251466900594, 5.880494279982013, 5.880556917690815, 10.68989495348153, 10.698983267284625, 10.718281105545595, 12.323621310863677, 12.515134575011572, 13.530915636425691, 13.646352490652179, 14.013638403879586, 15.844042103976731]
5 : [3.2202514327806324, 5.8804771145376495, 5.880480023263936, 10.686851036555902, 10.694700891026182, 10.700890407056217, 12.319109793665222, 12.416240501262315, 13.409477889981362, 13.552567887201173, 13.749216311234289, 15.24645246920046]
6 : [3.2202514321849844, 5.88047389051128, 5.880477843318861, 10.686161008731071, 10.694032876758136, 10.695820538590818, 12.31800240065333, 12.365108987078608, 13.375722155077574, 13.488497017403533, 13.613981824882925, 14.813746366173842]
7 : [3.220251432174622, 5.8804736764844865, 5.880477751399404, 10.68596928719655, 10.693915072917102, 10.694415397844669, 12.31769621979785, 12.338610207932499, 13.365112608794512, 13.452949540228428, 13.531662729002376, 14.548118346629236]
8 : [3.220251432174446, 5.880473664772776, 5.880477746015271, 10.685917731763457, 10.693892108137465, 10.694043457650729, 12.31760466508529, 12.326098068794337, 13.360878384968268, 13.435960447013805, 13.480566199422057, 14.398036525810834]
9 : [3.220251432174443, 5.880473664135115, 5.8804777457059485, 10.685904872914655, 10.693887275314943, 10.693949433863288, 12.31757553783886, 12.320797130193794, 13.358877986665796, 13.428258104309563, 13.451525413158912, 14.31573509295545]
10 : [3.2202514321744387, 5.88047366410035, 5.8804777456885144, 10.68590179219398, 10.693886175010691, 10.693926453071706, 12.317565160998472, 12.318715754671292, 13.357892966788738, 13.424793974842684, 13.436298710383335, 14.270543584914366]
11 : [3.220251432174441, 5.880473664098428, 5.880477745687523, 10.685901066924849, 10.693885916583938, 10.693920952757653, 12.317558831677923, 12.317939833911655, 13.357414569517541, 13.42322803413977, 13.428747515910759, 14.245626506398507]
12 : [3.220251432174442, 5.8804736640983135, 5.880477745687473, 10.685900897346352, 10.693885856844643, 10.693919648674717, 12.31754556289751, 12.317668127325572, 13.357188433300282, 13.422515874390236, 13.425120104468698, 14.231830221943852]
13 : [3.22025143217444, 5.880473664098331, 5.880477745687477, 10.6859008577883, 10.693885843152465, 10.693919341037084, 12.317512860698969, 12.317598524346206, 13.35708376988285, 13.422189472078832, 13.423410836482482, 14.224178385277835]
14 : [3.2202514321744413, 5.880473664098318, 5.88047774568747, 10.685900848562275, 10.693885840013444, 10.693919268621753, 12.317488491941603, 12.317586420984828, 13.357035964254052, 13.422036752322292, 13.422615843646886, 14.219928026214834]
15 : [3.2202514321744404, 5.880473664098296, 5.880477745687489, 10.685900846410851, 10.693885839291614, 10.693919251606548, 12.31747849914713, 12.317583474825758, 13.357014277401708, 13.421959827521462, 13.422253922026647, 14.217568607014107]
16 : [3.220251432174442, 5.880473664098395, 5.880477745687444, 10.685900845907693, 10.693885839124693, 10.693919247607232, 12.317474828857291, 12.317582559524578, 13.357004434880315, 13.421913869421301, 13.42209688283527, 14.216256934361406]
17 : [3.2202514321744444, 5.880473664098548, 5.880477745687858, 10.685900845788284, 10.693885839086569, 10.693919246655192, 12.317473499388399, 12.317582248191327, 13.356999779715625, 13.421882538079709, 13.422033075782478, 14.215502659702361]
18 : [3.2202514321744418, 5.880473664098304, 5.880477745687449, 10.685900845762445, 10.693885839077671, 10.693919246441288, 12.317473036283097, 12.31758214088091, 13.356997808114476, 13.421864710279193, 13.42200735007939, 14.215091282405693]
19 : [3.22025143217444, 5.880473664098828, 5.8804777456877995, 10.685900845756853, 10.693885839075703, 10.693919246392799, 12.317472875717511, 12.317582103355843, 13.356997003781094, 13.421855980796893, 13.421996672745061, 14.21487478713152]
[6]:
print ("Eigenvalues")
for lam in evals:
print (lam)
Eigenvalues
3.22025143217444
5.880473664098828
5.8804777456877995
10.685900845756853
10.693885839075703
10.693919246392799
12.317472875717511
12.317582103355843
13.356997003781094
13.421855980796893
13.421996672745061
14.21487478713152
[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);
[ ]:
[ ]: