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.2290932640676715, 22.636447498078578, 38.75912680734666, 56.95688372486332, 64.6890467454977, 71.61510644320737, 76.89311293050173, 86.36950975139492, 90.92730267366093, 99.90890224081902, 108.07709801039714, 123.53078334313035]
1 : [3.2310126754320967, 6.086234263756573, 6.4068992239104, 12.054533474383147, 12.672421028710486, 13.138599685848519, 14.490032247820023, 17.723681707636434, 18.662584400145292, 19.346120604991633, 23.999160961128645, 25.451852239368158]
2 : [3.2203793481865093, 5.886623822609898, 5.905548429291832, 10.853747390847467, 10.901287970102448, 11.033338300442013, 12.490608609970772, 13.713193968068527, 14.184813087521942, 15.348674546132722, 16.52982807771912, 17.964540348493433]
3 : [3.2202534507267346, 5.880768708677604, 5.8819084880279435, 10.709133420836263, 10.726052065998536, 10.778656219105256, 12.345960396974762, 12.760676043232321, 13.760281355795362, 14.018326043413813, 14.607718515287978, 16.59221415049482]
4 : [3.220251466900589, 5.880494279982001, 5.880556917690774, 10.68989495348608, 10.698983267283632, 10.718281105552407, 12.323621310867818, 12.515134575186735, 13.530915636339655, 13.646352490592061, 14.013638403669344, 15.844042102923868]
5 : [3.2202514327806298, 5.88047711453763, 5.880480023263914, 10.686851036556737, 10.694700891025324, 10.700890407059173, 12.31910979366547, 12.416240501356379, 13.409477889966896, 13.552567887179052, 13.749216311103094, 15.246452468154086]
6 : [3.2202514321849813, 5.880473890511276, 5.880477843318869, 10.686161008731137, 10.694032876757896, 10.695820538591713, 12.318002400653105, 12.365108987118434, 13.375722155084715, 13.48849701738329, 13.613981824783764, 14.813746365384908]
7 : [3.2202514321746234, 5.880473676484461, 5.880477751399416, 10.685969287196578, 10.693915072917068, 10.694415397844924, 12.317696219797673, 12.338610207946488, 13.365112608803019, 13.452949540217983, 13.531662728912307, 14.54811834611296]
8 : [3.2202514321744475, 5.88047366477278, 5.8804777460152655, 10.685917731763451, 10.693892108137478, 10.694043457650848, 12.317604665085243, 12.326098068795398, 13.360878384974317, 13.435960447019184, 13.480566199321775, 14.398036525405868]
9 : [3.2202514321744413, 5.880473664135103, 5.880477745705951, 10.685904872914751, 10.693887275313685, 10.693949433865388, 12.317575537839284, 12.320797131334485, 13.358877986635425, 13.428258104801749, 13.451525426650582, 14.315735114849959]
10 : [3.2202514321744418, 5.88047366410032, 5.880477745688513, 10.68590179219398, 10.693886175010373, 10.693926453072976, 12.317565160998557, 12.318715754647506, 13.357892966673496, 13.42479397678198, 13.43629871672295, 14.270543585464363]
11 : [3.220251432174443, 5.880473664098409, 5.880477745687522, 10.685901066925057, 10.693885916583662, 10.6939209527596, 12.31755883168096, 12.317939833965738, 13.357414569528464, 13.423228035681403, 13.428747516943364, 14.245626508849927]
12 : [3.2202514321744373, 5.880473664098311, 5.8804777456874575, 10.685900897346421, 10.693885856844929, 10.693919648675772, 12.317545562900166, 12.317668127297976, 13.35718843341429, 13.422515876436583, 13.425120104787608, 14.231830218560765]
13 : [3.220251432174442, 5.880473664098299, 5.88047774568747, 10.685900857789196, 10.693885843152872, 10.693919341041417, 12.317512863379847, 12.317598526477818, 13.357083769440584, 13.422189473618312, 13.423410995540578, 14.224178899618686]
14 : [3.2202514321744404, 5.880473664098306, 5.880477745687475, 10.685900848563929, 10.693885840013806, 10.69391926863341, 12.317488502368565, 12.317586424627294, 13.357035965441455, 13.42203676027171, 13.422616341726592, 14.219930431459408]
15 : [3.220251432174441, 5.880473664098103, 5.880477745687305, 10.68590084641107, 10.693885839291658, 10.693919251608584, 12.317478501739599, 12.317583475392096, 13.357014277660108, 13.421959835357232, 13.422254058432934, 14.217569698808399]
16 : [3.2202514321744444, 5.880473664098259, 5.88047774568718, 10.685900845908385, 10.693885839123526, 10.693919247588138, 12.317474780343856, 12.31758254432825, 13.357004463296505, 13.421913737980107, 13.42209550674888, 14.216247202503867]
17 : [3.2202514321744435, 5.880473664098062, 5.880477745687318, 10.685900845791055, 10.693885839086164, 10.69391924666083, 12.317473489692066, 12.317582242333822, 13.35700003182617, 13.421883566425674, 13.42203248400186, 14.215521671695729]
18 : [3.22025143217444, 5.880473664098331, 5.88047774568727, 10.685900845763435, 10.693885839077549, 10.693919246442888, 12.317473032339883, 12.317582138021306, 13.356998029315406, 13.421865705548868, 13.422007266696488, 14.215114701782642]
19 : [3.2202514321744413, 5.880473664097927, 5.880477745687271, 10.6859008457572, 10.693885839075747, 10.693919246393076, 12.317472874214646, 12.317582102276639, 13.35699712930589, 13.421856546368021, 13.421996637888197, 14.214890601219734]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.2202514321744413
5.880473664097927
5.880477745687271
10.6859008457572
10.693885839075747
10.693919246393076
12.317472874214646
12.317582102276639
13.35699712930589
13.421856546368021
13.421996637888197
14.214890601219734
[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);
[ ]:

[ ]: