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.229093261626065, 22.63644749338552, 38.75912681375783, 56.956883728123564, 64.68904675576682, 71.61510645059002, 76.89311293857456, 86.36950975583576, 90.9273026766007, 99.90890224161701, 108.0770980124297, 123.53078334086594]
1 : [3.231012675423183, 6.086234263795656, 6.406899224154681, 12.054533475318026, 12.672421029004822, 13.138599687726135, 14.490032249028651, 17.723681706627644, 18.66258439727946, 19.346120603377074, 23.999160958049387, 25.451852238023942]
2 : [3.2203793481863734, 5.886623822611853, 5.905548429296161, 10.853747390859365, 10.90128797022997, 11.03333830042761, 12.490608609958674, 13.713193966979425, 14.184813087330264, 15.348674543986524, 16.5298280738234, 17.964540348625516]
3 : [3.2202534507267275, 5.880768708677661, 5.881908488027613, 10.709133420844385, 10.726052065981984, 10.778656219061084, 12.345960396955327, 12.760676042694097, 13.76028135549549, 14.018326041670313, 14.60771851555034, 16.592214151523518]
4 : [3.220251466900592, 5.880494279981982, 5.880556917690745, 10.689894953487357, 10.698983267276608, 10.718281105541074, 12.323621310858588, 12.515134574909926, 13.530915635669155, 13.646352490774788, 14.013638403773562, 15.844042104325233]
5 : [3.2202514327806293, 5.880477114537646, 5.880480023263944, 10.68685103655631, 10.69470089102384, 10.700890407057795, 12.319109793661559, 12.416240501230623, 13.409477889754987, 13.552567887351042, 13.749216311055427, 15.24645246955614]
6 : [3.2202514321849853, 5.880473890511292, 5.8804778433188645, 10.68616100873086, 10.694032876757584, 10.695820538591938, 12.318002400651507, 12.365108987070226, 13.37572215498713, 13.488497017519162, 13.613981824743423, 14.813746366466985]
7 : [3.2202514321746234, 5.880473676484483, 5.880477751399428, 10.685969287196475, 10.693915072916983, 10.694415397845098, 12.317696219797071, 12.338610207930854, 13.365112608748332, 13.452949540298915, 13.53166272892073, 14.548118346833268]
8 : [3.2202514321744435, 5.880473664772775, 5.880477746015279, 10.68591773176341, 10.693892108137447, 10.694043457650864, 12.317604665084954, 12.326098068792684, 13.360878384943904, 13.43596044705922, 13.480566199356353, 14.398036525861524]
9 : [3.2202514321744413, 5.880473664135104, 5.880477745705942, 10.685904872914922, 10.693887275313799, 10.693949433866349, 12.317575537840352, 12.320797131749098, 13.358877986626432, 13.428258104684414, 13.451525433116297, 14.315735126454044]
10 : [3.2202514321744427, 5.880473664100304, 5.8804777456884905, 10.685901792194066, 10.693886175010226, 10.693926453073486, 12.31756516099841, 12.31871575522362, 13.357892966596037, 13.424793975517378, 13.436298720850566, 14.270543601183784]
11 : [3.2202514321744418, 5.880473664098407, 5.8804777456875374, 10.685901066925085, 10.693885916583774, 10.693920952759292, 12.31755883168186, 12.317939834440978, 13.357414569367428, 13.423228034658907, 13.428747526542784, 14.245626530643408]
12 : [3.220251432174441, 5.8804736640983375, 5.880477745687479, 10.685900897346311, 10.693885856844904, 10.693919648675703, 12.317545562922447, 12.317668127503731, 13.357188433418125, 13.422515876219691, 13.425120109048084, 14.231830226165545]
13 : [3.22025143217444, 5.880473664098289, 5.880477745687487, 10.685900857787997, 10.693885843152815, 10.693919341027518, 12.31751285780382, 12.317598522119853, 13.357083767128994, 13.422189469842703, 13.42341080079679, 14.224177892712358]
14 : [3.220251432174439, 5.880473664098315, 5.880477745687479, 10.685900848563882, 10.693885840012998, 10.693919268630141, 12.317488500577282, 12.317586423926784, 13.357035966734996, 13.422036750343985, 13.422616243223562, 14.219930125125831]
15 : [3.220251432174442, 5.880473664098296, 5.880477745687449, 10.68590084641146, 10.693885839291502, 10.693919251609502, 12.317478504145496, 12.317583476113894, 13.35701427866216, 13.421959834728918, 13.422254123068234, 14.217570392314398]
16 : [3.220251432174441, 5.880473664098252, 5.880477745687451, 10.685900845908828, 10.693885839124798, 10.693919247606278, 12.317474825430661, 12.317582558243652, 13.357004463513023, 13.421913949886536, 13.422096849221253, 14.216258242544551]
17 : [3.220251432174442, 5.8804736640983775, 5.880477745687469, 10.68590084579095, 10.69388583908635, 10.693919246666312, 12.31747350776662, 12.31758224832303, 13.357000021528743, 13.421883597515663, 13.422032534934887, 14.215522998388288]
18 : [3.2202514321744413, 5.880473664098203, 5.880477745687358, 10.68590084576359, 10.693885839077486, 10.693919246445132, 12.317473042696466, 12.31758214085251, 13.35699800920354, 13.421865815057622, 13.422005678611717, 14.215117407489116]
19 : [3.220251432174445, 5.880473664098269, 5.880477745687479, 10.685900845757395, 10.693885839075726, 10.693919246393579, 12.31747287633408, 12.317582102944929, 13.356997110334094, 13.421856408276646, 13.421995508266118, 14.214887099701665]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.220251432174445
5.880473664098269
5.880477745687479
10.685900845757395
10.693885839075726
10.693919246393579
12.31747287633408
12.317582102944929
13.356997110334094
13.421856408276646
13.421995508266118
14.214887099701665
[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);
[ ]:

[ ]: