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.229093261554611, 22.6364474896392, 38.75912680748503, 56.95688371352073, 64.68904674884229, 71.61510645017684, 76.89311294247331, 86.36950975196336, 90.92730267419756, 99.90890223856934, 108.07709801421511, 123.53078333788197]
1 : [3.2310126754204296, 6.086234263672357, 6.406899224207909, 12.0545334738457, 12.672421029004843, 13.138599687690089, 14.490032248169685, 17.72368170692177, 18.662584397409915, 19.346120605257713, 23.999160956743598, 25.451852237887497]
2 : [3.2203793481863903, 5.88662382260831, 5.905548429299967, 10.853747390766365, 10.90128797009822, 11.033338300567749, 12.490608609924944, 13.713193967396064, 14.18481308751501, 15.348674545167993, 16.52982807253465, 17.964540348799492]
3 : [3.2202534507267275, 5.880768708677612, 5.88190848802809, 10.709133420836295, 10.726052065965897, 10.778656219104393, 12.345960396953062, 12.760676042893522, 13.760281355720933, 14.018326042033756, 14.60771851511549, 16.59221415111298]
4 : [3.2202514669005913, 5.880494279981987, 5.880556917690763, 10.689894953486982, 10.6989832672742, 10.718281105552096, 12.323621310859927, 12.515134575001987, 13.530915635865984, 13.646352490735728, 14.013638403716635, 15.84404210365321]
5 : [3.2202514327806284, 5.880477114537624, 5.880480023263904, 10.686851036556439, 10.694700891023524, 10.700890407059525, 12.319109793662513, 12.416240501266042, 13.409477889821686, 13.552567887314392, 13.749216311100092, 15.246452468851881]
6 : [3.220251432184985, 5.88047389051128, 5.880477843318851, 10.686161008730926, 10.69403287675755, 10.695820538591967, 12.318002400651947, 12.36510898708022, 13.375722155024581, 13.488497017477131, 13.61398182478872, 14.813746365907502]
7 : [3.220251432174623, 5.880473676484487, 5.880477751399398, 10.685969287196478, 10.693915072916951, 10.694415397844995, 12.317696219797275, 12.338610207932215, 13.365112608771195, 13.452949540270517, 13.531662728941152, 14.54811834646035]
8 : [3.220251432174444, 5.88047366477277, 5.880477746015271, 10.685917731763416, 10.693892108137472, 10.69404345765087, 12.317604665085062, 12.326098068790804, 13.360878384955393, 13.435960447053013, 13.480566199374087, 14.398036525611557]
9 : [3.22025143217444, 5.880473664135124, 5.880477745705931, 10.685904872914806, 10.693887275313523, 10.693949433865914, 12.317575537838982, 12.320797131229375, 13.358877986619806, 13.428258103823651, 13.451525421427815, 14.315735113603465]
10 : [3.2202514321744435, 5.880473664100323, 5.880477745688504, 10.685901792194187, 10.693886175010086, 10.693926453073283, 12.317565160998456, 12.318715755302266, 13.357892966935895, 13.424793973525574, 13.43629871682694, 14.270543609002766]
11 : [3.2202514321744453, 5.880473664098407, 5.880477745687536, 10.685901066925227, 10.693885916583382, 10.693920952759951, 12.317558831689299, 12.317939834480448, 13.357414569796244, 13.42322803298266, 13.428747519187427, 14.245626538820362]
12 : [3.2202514321744404, 5.880473664098324, 5.880477745687468, 10.685900897346462, 10.693885856844583, 10.69391964867584, 12.317545562945691, 12.317668127557312, 13.357188433599058, 13.422515874379059, 13.42512010463528, 14.231830243083355]
13 : [3.220251432174444, 5.880473664098312, 5.8804777456874735, 10.685900857789315, 10.69388584315269, 10.693919341042074, 12.317512863784449, 12.317598526743975, 13.357083768974688, 13.422189472090365, 13.423411009096112, 14.224178981580206]
14 : [3.220251432174436, 5.880473664098326, 5.880477745687481, 10.685900848564197, 10.693885840013607, 10.693919268635216, 12.317488503265084, 12.317586424834564, 13.357035967530834, 13.422036759537509, 13.422616315999923, 14.219930633550412]
15 : [3.2202514321744418, 5.880473664098288, 5.880477745687387, 10.685900846411414, 10.693885839291628, 10.693919251610897, 12.317478505003328, 12.317583476362064, 13.357014277617727, 13.421959841890756, 13.42225414726219, 14.217570533614884]
16 : [3.2202514321744453, 5.880473664098204, 5.88047774568745, 10.685900845908678, 10.693885839125253, 10.693919247606573, 12.317474825148732, 12.317582558368196, 13.35700446920957, 13.42191393696656, 13.422096706993932, 14.216257442080975]
17 : [3.220251432174441, 5.880473664098545, 5.880477745687504, 10.6859008457914, 10.693885839086875, 10.693919246667209, 12.317473511273889, 12.317582248896846, 13.357000034073389, 13.421883805839924, 13.42203304446062, 14.215528992154134]
18 : [3.220251432174443, 5.880473664098043, 5.880477745687373, 10.68590084576364, 10.693885839077788, 10.693919246445887, 12.317473044100453, 12.317582141440447, 13.356998026173114, 13.421865911793311, 13.42200753044096, 14.215120651704293]
19 : [3.2202514321744418, 5.88047366409831, 5.880477745687449, 10.685900845757374, 10.693885839075902, 10.693919246394127, 12.317472879846449, 12.317582103874495, 13.35699712754991, 13.421856691812918, 13.421996820931644, 14.214895293993145]
[6]:
print ("Eigenvalues")
for lam in evals:
print (lam)
Eigenvalues
3.2202514321744418
5.88047366409831
5.880477745687449
10.685900845757374
10.693885839075902
10.693919246394127
12.317472879846449
12.317582103874495
13.35699712754991
13.421856691812918
13.421996820931644
14.214895293993145
[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);
[ ]:
[ ]: