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.229093262876329, 22.636447487565587, 38.75912681001847, 56.95688372853938, 64.68904675565132, 71.61510644608987, 76.89311295325861, 86.36950975405644, 90.92730267438486, 99.90890224369574, 108.07709801262804, 123.53078334327188]
1 : [3.2310126754277753, 6.086234263743855, 6.406899224145533, 12.054533475487139, 12.672421029238283, 13.138599688006257, 14.490032248523168, 17.723681706777707, 18.662584398368857, 19.346120603359985, 23.999160956893114, 25.45185223960607]
2 : [3.220379348186461, 5.886623822611986, 5.905548429292463, 10.85374739088961, 10.90128797017727, 11.03333830053192, 12.490608609972018, 13.71319396687997, 14.184813087430836, 15.348674545243378, 16.529828073123493, 17.964540349234802]
3 : [3.2202534507267293, 5.880768708677722, 5.881908488027432, 10.709133420839269, 10.72605206599649, 10.778656219076398, 12.345960396968824, 12.760676042672692, 13.760281355614573, 14.018326042188463, 14.607718515278563, 16.59221415155412]
4 : [3.2202514669005913, 5.880494279982001, 5.8805569176907095, 10.689894953487666, 10.69898326728038, 10.718281105540633, 12.323621310865752, 12.515134574919047, 13.530915635881428, 13.646352490696568, 14.013638403796959, 15.844042104281364]
5 : [3.2202514327806324, 5.88047711453763, 5.880480023263911, 10.68685103655679, 10.694700891024665, 10.700890407056102, 12.319109793664833, 12.416240501239523, 13.409477889828345, 13.552567887279137, 13.749216311159621, 15.246452469476559]
6 : [3.2202514321849853, 5.880473890511284, 5.880477843318856, 10.686161008731, 10.694032876757753, 10.695820538591141, 12.318002400652855, 12.365108987075766, 13.375722155025212, 13.488497017461171, 13.613981824842254, 14.813746366380316]
7 : [3.2202514321746216, 5.880473676484499, 5.880477751399413, 10.685969287196503, 10.69391507291704, 10.69441539784488, 12.317696219797602, 12.33861020793381, 13.365112608770934, 13.452949540265518, 13.531662728975995, 14.548118346741209]
8 : [3.220251432174447, 5.8804736647727776, 5.880477746015272, 10.68591773176343, 10.693892108137453, 10.694043457650816, 12.317604665085199, 12.3260980687952, 13.360878384955647, 13.4359604470533, 13.480566199412669, 14.398036525787905]
9 : [3.220251432174438, 5.88047366413509, 5.880477745705945, 10.68590487291503, 10.693887275313731, 10.693949433865612, 12.317575537839401, 12.32079713122775, 13.358877986624323, 13.42825810396055, 13.451525423671008, 14.31573511533798]
10 : [3.2202514321744418, 5.880473664100352, 5.88047774568851, 10.685901792193844, 10.693886175010048, 10.693926453073033, 12.317565160998443, 12.31871575399306, 13.357892966475195, 13.42479397589448, 13.436298695091162, 14.270543553012201]
11 : [3.2202514321744427, 5.880473664098444, 5.88047774568755, 10.685901066924874, 10.693885916583701, 10.693920952758898, 12.317558831674791, 12.317939833571321, 13.357414569280134, 13.423228035087485, 13.428747502775414, 14.24562647805286]
12 : [3.220251432174442, 5.8804736640983295, 5.8804777456874815, 10.685900897346285, 10.69388585684459, 10.693919648675418, 12.31754556286531, 12.317668127085224, 13.35718843326108, 13.422515875012309, 13.425120084398063, 14.231830179433517]
13 : [3.2202514321744395, 5.880473664098294, 5.880477745687464, 10.685900857786184, 10.693885843152392, 10.693919341026042, 12.317512856847157, 12.317598521824603, 13.357083765130387, 13.422189470576487, 13.423410669317917, 14.224177483318888]
14 : [3.220251432174437, 5.880473664098302, 5.880477745687488, 10.68590084856304, 10.693885840013989, 10.693919268628326, 12.3174884980931, 12.317586423197216, 13.35703596261163, 13.42203675431574, 13.422616157985994, 14.219929393971334]
15 : [3.2202514321744387, 5.880473664098296, 5.880477745687494, 10.685900846411165, 10.693885839291848, 10.693919251608818, 12.317478502868317, 12.31758347582379, 13.35701427458835, 13.421959827976062, 13.422254111603811, 14.217569862698568]
16 : [3.2202514321744404, 5.880473664098484, 5.880477745687631, 10.68590084590883, 10.693885839125228, 10.693919247604926, 12.317474823117934, 12.317582558005139, 13.35700446907883, 13.421913905058688, 13.422096747780834, 14.216256592816936]
17 : [3.220251432174443, 5.880473664098284, 5.8804777456875, 10.68590084579131, 10.69388583908604, 10.69391924666551, 12.317473510677312, 12.317582248669993, 13.357000024060463, 13.421883788920391, 13.422031957240703, 14.215527874207035]
18 : [3.2202514321744413, 5.880473664098225, 5.8804777456874255, 10.685900845763852, 10.693885839077629, 10.693919246445631, 12.31747304530691, 12.317582141451174, 13.35699802250174, 13.4218660295448, 13.42200646188931, 14.215123159531064]
19 : [3.2202514321744395, 5.880473664098352, 5.880477745687378, 10.685900845757168, 10.693885839075799, 10.693919246393296, 12.31747287824808, 12.31758210351184, 13.35699708625949, 13.421856580658305, 13.421996175615922, 14.214892177134315]
[6]:
print ("Eigenvalues")
for lam in evals:
print (lam)
Eigenvalues
3.2202514321744395
5.880473664098352
5.880477745687378
10.685900845757168
10.693885839075799
10.693919246393296
12.31747287824808
12.31758210351184
13.35699708625949
13.421856580658305
13.421996175615922
14.214892177134315
[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);
[ ]:
[ ]: