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.229093260825753, 22.636447504214424, 38.75912681473533, 56.956883709480536, 64.68904673231133, 71.61510643781894, 76.89311293495308, 86.36950976414586, 90.92730267997553, 99.90890224292198, 108.07709801699592, 123.53078333441034]
1 : [3.2310126754112583, 6.086234263603961, 6.406899224116552, 12.05453347107921, 12.672421029415977, 13.138599689085096, 14.49003224741752, 17.72368170594841, 18.66258439255411, 19.346120606504748, 23.99916095879927, 25.451852236393055]
2 : [3.220379348186288, 5.886623822600107, 5.905548429297898, 10.85374739057017, 10.901287970041755, 11.033338300620992, 12.49060860982657, 13.713193966591563, 14.184813087232358, 15.348674546715639, 16.52982807743296, 17.964540348690743]
3 : [3.2202534507267284, 5.880768708677462, 5.881908488028285, 10.709133420847996, 10.726052065943048, 10.77865621908992, 12.345960396944236, 12.760676042898405, 13.76028135556215, 14.018326043717467, 14.607718515666967, 16.592214150052033]
4 : [3.2202514669005957, 5.880494279981973, 5.8805569176907975, 10.689894953487622, 10.698983267277733, 10.718281105538138, 12.323621310859162, 12.515134575006034, 13.530915636281819, 13.646352490526864, 14.013638403861949, 15.84404210224441]
5 : [3.2202514327806284, 5.880477114537658, 5.880480023263937, 10.686851036557208, 10.69470089102459, 10.700890407052968, 12.319109793662927, 12.416240501240198, 13.409477889920637, 13.552567887159208, 13.749216311075253, 15.24645246755984]
6 : [3.220251432184982, 5.880473890511275, 5.88047784331886, 10.68616100873123, 10.69403287675778, 10.695820538589528, 12.31800240065238, 12.365108987050705, 13.375722155055717, 13.48849701739742, 13.6139818246655, 14.813746364996375]
7 : [3.220251432174621, 5.880473676484482, 5.880477751399405, 10.685969287196562, 10.693915072917024, 10.69441539784423, 12.317696219797519, 12.338610207912097, 13.365112608788325, 13.452949540229135, 13.53166272879369, 14.548118345897219]
8 : [3.2202514321744444, 5.880473664772772, 5.8804777460152815, 10.685917731763459, 10.693892108137463, 10.694043457650544, 12.317604665085184, 12.326098068789223, 13.360878384967489, 13.435960447017273, 13.480566199374348, 14.398036525517881]
9 : [3.220251432174444, 5.88047366413511, 5.880477745705945, 10.685904872914119, 10.693887275315035, 10.693949433862244, 12.317575537839833, 12.320797130491524, 13.35887798666256, 13.428258106498177, 13.451525424270823, 14.315735096873906]
10 : [3.220251432174439, 5.880473664100326, 5.880477745688514, 10.685901792193484, 10.693886175010858, 10.693926453071182, 12.317565160998472, 12.318715754055436, 13.357892966442634, 13.424793979049443, 13.436298712049943, 14.270543557912326]
11 : [3.2202514321744378, 5.880473664098426, 5.880477745687519, 10.685901066925048, 10.693885916584085, 10.693920952758805, 12.317558831676903, 12.317939833712318, 13.35741456939229, 13.4232280374203, 13.428747518380026, 14.245626495521021]
12 : [3.2202514321744418, 5.880473664098309, 5.880477745687481, 10.685900897346094, 10.693885856844844, 10.693919648675907, 12.317545562849233, 12.317668127026012, 13.357188433830292, 13.422515879330264, 13.42512008783242, 14.231830157982266]
13 : [3.2202514321744404, 5.880473664098394, 5.8804777456870525, 10.685900857789136, 10.693885843152236, 10.693919341040614, 12.317512862637729, 12.317598525676281, 13.35708376945059, 13.422189453955939, 13.42341086172385, 14.224178701165693]
14 : [3.220251432174439, 5.88047366409834, 5.880477745687501, 10.685900848561968, 10.693885840013328, 10.693919268617547, 12.317488493562667, 12.31758642178162, 13.357035953014012, 13.422036713350128, 13.422615978673337, 14.219928315048291]
15 : [3.22025143217444, 5.880473664098331, 5.880477745687434, 10.685900846410991, 10.693885839291505, 10.693919251607431, 12.317478502362691, 12.317583475781976, 13.357014272062736, 13.42195981828595, 13.422254064733059, 14.217569506168001]
16 : [3.220251432174443, 5.880473664098403, 5.880477745687439, 10.68590084590864, 10.693885839125299, 10.693919247608303, 12.317474831053204, 12.317582560292268, 13.357004466076493, 13.421913944057387, 13.422096932460589, 14.216258526640583]
17 : [3.220251432174444, 5.880473664097914, 5.880477745687589, 10.685900845790503, 10.693885839086773, 10.693919246662523, 12.317473500017883, 12.317582246846309, 13.356999996621346, 13.421883255954656, 13.422032062703178, 14.215514665377329]
18 : [3.2202514321744404, 5.880473664098529, 5.880477745687408, 10.685900845763516, 10.693885839077746, 10.69391924644463, 12.31747302795509, 12.317582137194268, 13.356998022700203, 13.421864982792625, 13.42200475695945, 14.215095408295934]
19 : [3.220251432174442, 5.880473664098318, 5.880477745687529, 10.68590084575726, 10.69388583907585, 10.693919246394028, 12.31747287301002, 12.317582102024184, 13.35699712940073, 13.42185608681627, 13.421995141281284, 14.214878078139865]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.220251432174442
5.880473664098318
5.880477745687529
10.68590084575726
10.69388583907585
10.693919246394028
12.31747287301002
12.317582102024184
13.35699712940073
13.42185608681627
13.421995141281284
14.214878078139865
[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);
[ ]:

[ ]: