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 = 44384
[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 : [5.921595127433707, 16.354970677863818, 19.60560289013183, 42.584719746277685, 47.467538034935494, 53.51420727552569, 67.47212897355165, 73.88926358371477, 82.35531731758853, 91.53567639160576, 101.12543372502425, 108.83008164864181]
1 : [3.2324570677895967, 6.017167187364033, 6.219423864761252, 11.719693734564197, 12.177834962044797, 12.681820735354565, 14.477850968266194, 15.596821198947303, 16.85069621997704, 17.290963314218008, 20.775093723940735, 24.103812492296395]
2 : [3.220424486541087, 5.884453035252117, 5.895380849667886, 10.805685741565204, 10.948315771148252, 11.003449900607443, 12.687851537996634, 13.149348631050373, 13.947130545437773, 14.165023180049966, 15.461066387944742, 18.936725209485676]
3 : [3.220285796313135, 5.880622718639504, 5.881215482937746, 10.703348608250549, 10.751617496013164, 10.770399530229296, 12.441884160883141, 12.578846769862833, 13.575857444962361, 13.652071371338257, 14.353411534065934, 16.7012254546565]
4 : [3.2202839009329463, 5.880495053921094, 5.880522638349874, 10.688469331968724, 10.708276465930716, 10.714592150445931, 12.36355262071544, 12.40931798315846, 13.44283323906035, 13.50927096491628, 13.962225966058028, 15.353134427727682]
5 : [3.2202838706905514, 5.880486677009537, 5.880490627128294, 10.68637370616118, 10.69736443288755, 10.699377121645709, 12.334461399135103, 12.351675241307893, 13.3985866228698, 13.452338702911295, 13.73057491961914, 14.724791874474311]
6 : [3.220283870199275, 5.88048483270961, 5.880490431251677, 10.686060068934589, 10.694716372781372, 10.69528956659008, 12.323561396443731, 12.330081934910233, 13.377256158849724, 13.432622030889625, 13.579314947434296, 14.45353970737025]
7 : [3.2202838701913827, 5.880484735957958, 5.880490423365247, 10.686008883815202, 10.694102643839017, 10.694266634968495, 12.319694706542975, 12.322036880510716, 13.366371151462532, 13.426058619779393, 13.496332828806223, 14.334088332300103]
8 : [3.220283870191257, 5.880484730924003, 5.880490423013741, 10.685999809324988, 10.69395395603713, 10.694032711789937, 12.31834761559844, 12.319158839115449, 13.361210678067012, 13.423642238758552, 13.45579809379924, 14.276849045975545]
9 : [3.220283870191253, 5.88048473066162, 5.880490422997119, 10.68599808484814, 10.693915270697072, 10.693983669156404, 12.3178697548861, 12.318171918643861, 13.358856073446686, 13.422668881524144, 13.437079418676516, 14.247753582893575]
10 : [3.2202838701912557, 5.880484730647931, 5.880490422996332, 10.685997738746785, 10.693905985112476, 10.6939730012205, 12.317678925977235, 12.317860021098326, 13.357800345637251, 13.42225454345881, 13.428645258731164, 14.232476437393547]
11 : [3.220283870191251, 5.880484730647212, 5.8804904229963055, 10.685997666374515, 10.693903874825349, 10.693970583354016, 12.317599091168594, 12.317767809236793, 13.357330167561553, 13.422072584833519, 13.424886367364797, 14.22431694921677]
12 : [3.2202838701912544, 5.880484730647127, 5.880490422996282, 10.68599765080436, 10.693903402355513, 10.693970028456224, 12.317569385704447, 12.317738635955068, 13.357121483378618, 13.421991026968305, 13.423219930206388, 14.219922058401139]
13 : [3.220283870191256, 5.880484730647169, 5.88049042299629, 10.685997647387051, 10.693903296789484, 10.693969900680258, 12.317559010797355, 12.317728834339029, 13.35702893246204, 13.421953300457925, 13.422483205990845, 14.217543264246371]
14 : [3.2202838701912553, 5.880484730647147, 5.88049042299627, 10.685997646627062, 10.693903273194879, 10.693969871223935, 12.317555460211334, 12.31772546225146, 13.356987909605069, 13.421933653303986, 13.422160048480064, 14.21625287831094]
15 : [3.2202838701912553, 5.880484730646912, 5.880490422996122, 10.68599764645643, 10.693903267909304, 10.693969864428594, 12.317554251041942, 12.317724292298514, 13.356969718577622, 13.421918431346006, 13.422023438694598, 14.215551621704071]
16 : [3.2202838701912557, 5.880484730647208, 5.880490422996352, 10.685997646417958, 10.693903266721925, 10.693969862859094, 12.317553839174483, 12.317723884598182, 13.356961646682613, 13.421901166103321, 13.421973206133405, 14.215169884285647]
17 : [3.2202838701912557, 5.880484730647443, 5.88049042299651, 10.685997646409504, 10.693903266455553, 10.693969862497186, 12.317553698737848, 12.31772374256367, 13.356958056817465, 13.421886952616042, 13.421957469223917, 14.214961455344698]
18 : [3.2202838701912553, 5.880484730647171, 5.88049042299624, 10.685997646407358, 10.693903266394857, 10.69396986241371, 12.317553649870781, 12.317723692840161, 13.356956470232978, 13.421878653766155, 13.42195188082979, 14.214845036461229]
19 : [3.2202838701912517, 5.88048473064725, 5.880490422996389, 10.685997646406843, 10.693903266381488, 10.693969862394358, 12.317553633631668, 12.317723675421995, 13.356955759559845, 13.421874678931134, 13.421949746919065, 14.214781293310333]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.2202838701912517
5.88048473064725
5.880490422996389
10.685997646406843
10.693903266381488
10.693969862394358
12.317553633631668
12.317723675421995
13.356955759559845
13.421874678931134
13.421949746919065
14.214781293310333
[7]:
gfu = GridFunction(fes, multidim=len(evecs))
for i in range(len(evecs)):
    gfu.vecs[i].data = evecs[i]
Draw (Norm(gfu), mesh, order=2, min=0, max=2);
[ ]:

[ ]: