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.229093262742357, 22.636447487082457, 38.759126811260714, 56.95688372610865, 64.68904675678684, 71.61510644515525, 76.89311294580897, 86.36950975189993, 90.92730267385816, 99.90890224403493, 108.07709800775478, 123.53078333996818]
1 : [3.2310126754252364, 6.086234263798401, 6.406899224282229, 12.054533475983463, 12.672421028987106, 13.138599687230903, 14.490032246475936, 17.72368170633454, 18.66258439795801, 19.346120603426726, 23.999160958038928, 25.45185223911863]
2 : [3.220379348186447, 5.8866238226129415, 5.905548429302214, 10.853747390986122, 10.90128797002589, 11.033338300459315, 12.490608609869618, 13.71319396675308, 14.184813087347859, 15.348674545293902, 16.529828073428245, 17.964540349229754]
3 : [3.2202534507267324, 5.880768708677698, 5.881908488028058, 10.709133420827273, 10.726052066009512, 10.778656219066178, 12.34596039696001, 12.760676042531856, 13.760281355535582, 14.018326042226303, 14.60771851543339, 16.592214151423132]
4 : [3.2202514669005917, 5.880494279982013, 5.88055691769075, 10.689894953485737, 10.698983267282863, 10.718281105538072, 12.3236213108656, 12.515134574823952, 13.53091563582474, 13.646352490746482, 14.013638403847624, 15.84404210412403]
5 : [3.220251432780628, 5.880477114537649, 5.88048002326392, 10.686851036556689, 10.69470089102504, 10.700890407054851, 12.319109793665321, 12.416240501181647, 13.40947788980739, 13.552567887309012, 13.74921631114937, 15.246452469356662]
6 : [3.2202514321849796, 5.8804738905112774, 5.880477843318853, 10.686161008731055, 10.69403287675782, 10.695820538590626, 12.318002400653155, 12.36510898704469, 13.375722155011434, 13.48849701748386, 13.61398182481423, 14.813746366314236]
7 : [3.2202514321746234, 5.880473676484506, 5.8804777513994075, 10.68596928719651, 10.693915072917031, 10.69441539784467, 12.317696219797716, 12.338610207919041, 13.365112608761436, 13.452949540276474, 13.531662728950614, 14.54811834671489]
8 : [3.2202514321744453, 5.880473664772772, 5.880477746015252, 10.685917731763428, 10.693892108137437, 10.69404345765072, 12.317604665085241, 12.326098068792652, 13.360878384949988, 13.435960447048478, 13.480566199426455, 14.398036525874012]
9 : [3.2202514321744444, 5.880473664135113, 5.880477745705941, 10.685904872914211, 10.693887275313946, 10.693949433863834, 12.317575537840002, 12.320797130880626, 13.358877986670121, 13.428258105497022, 13.451525422111038, 14.31573510151862]
10 : [3.2202514321744404, 5.880473664100319, 5.880477745688501, 10.685901792193658, 10.6938861750102, 10.69392645307208, 12.317565160998594, 12.31871575466271, 13.357892966512107, 13.424793976298137, 13.436298707968017, 14.270543569804094]
11 : [3.220251432174442, 5.880473664098429, 5.8804777456875135, 10.685901066924986, 10.693885916583577, 10.693920952759138, 12.31755883167907, 12.317939834071112, 13.357414569407176, 13.423228035164515, 13.428747512516054, 14.245626499871918]
12 : [3.2202514321744378, 5.880473664098329, 5.880477745687477, 10.685900897346382, 10.69388585684447, 10.693919648675339, 12.317545562914855, 12.317668127464893, 13.35718843332411, 13.422515875405576, 13.42512010605169, 14.231830225665949]
13 : [3.2202514321744404, 5.8804736640983055, 5.880477745687496, 10.685900857787665, 10.69388584315273, 10.693919341029341, 12.317512860207339, 12.31759852436993, 13.357083765796764, 13.422189463790536, 13.423410925377706, 14.224178325801285]
14 : [3.2202514321744395, 5.880473664098285, 5.880477745687475, 10.68590084856365, 10.693885840014259, 10.693919268630022, 12.317488500557563, 12.317586423977382, 13.35703596353972, 13.422036746794133, 13.422616344936062, 14.219930248889204]
15 : [3.2202514321744413, 5.880473664098303, 5.880477745687415, 10.685900846411444, 10.693885839291845, 10.693919251608952, 12.31747850343528, 12.317583475906806, 13.357014275762648, 13.42195982861365, 13.422254199944314, 14.217570386331738]
16 : [3.220251432174444, 5.880473664098352, 5.880477745687438, 10.685900845908778, 10.6938858391253, 10.693919247608402, 12.317474830415037, 12.317582560034838, 13.357004470523439, 13.42191395876631, 13.42209690719736, 14.216258707583462]
17 : [3.220251432174442, 5.880473664098118, 5.88047774568736, 10.685900845791087, 10.69388583908702, 10.693919246666004, 12.317473504969113, 12.317582247468916, 13.357000041256029, 13.421883599556823, 13.42203236355034, 14.215522731520002]
18 : [3.2202514321744427, 5.880473664098519, 5.880477745687389, 10.685900845763722, 10.693885839077602, 10.693919246445438, 12.31747304186677, 12.317582140693956, 13.356998032333534, 13.421865804472779, 13.422006172738055, 14.215117200123622]
19 : [3.220251432174444, 5.880473664098277, 5.880477745687555, 10.685900845757219, 10.693885839075755, 10.693919246394133, 12.317472878578654, 12.317582103358811, 13.356997123352352, 13.421856616664742, 13.421995777421708, 14.214893092129445]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.220251432174444
5.880473664098277
5.880477745687555
10.685900845757219
10.693885839075755
10.693919246394133
12.317472878578654
12.317582103358811
13.356997123352352
13.421856616664742
13.421995777421708
14.214893092129445
[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);
[ ]:

[ ]: