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 = 43192
[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)
/usr/lib/python3/dist-packages/scipy/__init__.py:146: UserWarning: A NumPy version >=1.17.3 and <1.25.0 is required for this version of SciPy (detected version 1.26.4
  warnings.warn(f"A NumPy version >={np_minversion} and <{np_maxversion}"
0 : [8.213129806662726, 17.799813247320422, 37.21982566943738, 44.0193003176738, 62.18233338007569, 65.53333261138647, 73.98716308730715, 81.0823127560318, 88.06001747633435, 91.25055694121237, 113.61980592134275, 122.21673396514657]
1 : [3.2380975391532516, 5.998650630676201, 6.161390072434293, 11.48586055068807, 12.037309679802714, 14.413520500297881, 15.234814499483306, 16.21826937226508, 17.346726951532553, 19.892283927168613, 22.638813223334143, 24.807560952554116]
2 : [3.220542688575675, 5.884876708903982, 5.888226095537274, 10.785471632103887, 10.847156783123758, 11.381986377827792, 13.343238319893644, 13.372706506237517, 14.11239328473971, 14.988639564364197, 15.34408001554822, 17.688329578805053]
3 : [3.2202954145601543, 5.880713446554603, 5.880839526591969, 10.707466526321701, 10.720069178386925, 10.944773310950346, 12.527973954720935, 13.075036550913648, 13.540697915881182, 14.030695748002074, 14.36749421859367, 15.906741588783488]
4 : [3.220289894102194, 5.880502253291043, 5.880518097321879, 10.695114363079908, 10.699129958352293, 10.78775181180235, 12.354098370697058, 12.871053039297564, 13.445683487798288, 13.649403154146059, 14.213859961104433, 14.607665451656969]
5 : [3.2202897664212093, 5.88049337525466, 5.8804959631163785, 10.692776367607363, 10.694999509694908, 10.721906317449033, 12.32420257888011, 12.640658946261182, 13.417157954192813, 13.5097622148142, 13.816759278701044, 14.260863114904971]
6 : [3.2202897638070236, 5.8804929781178235, 5.88049451044494, 10.691097181818124, 10.694161153653535, 10.698643594403919, 12.318883281765709, 12.460746605561056, 13.402394959817533, 13.45479587571817, 13.563293264194959, 14.232709516054495]
7 : [3.2202897637595154, 5.880492958816454, 5.880494425037921, 10.687861627687209, 10.69398579377986, 10.694552792299424, 12.317899571105752, 12.372364426791044, 13.387450842682792, 13.43393941280189, 13.473271158454446, 14.223431358536795]
8 : [3.2202897637586982, 5.880492957834384, 5.880494420312183, 10.686464368104195, 10.693942498545185, 10.694085266049292, 12.317711055329653, 12.337392167704833, 13.373735427638152, 13.426640065236331, 13.441852510732128, 14.219267826295189]
9 : [3.2202897637586827, 5.8804929577828275, 5.880494420056976, 10.686108698957005, 10.693925371312494, 10.694004893998423, 12.317673740214232, 12.324680164793959, 13.365194188415282, 13.423921247789114, 13.430320729442093, 14.217191023860313]
10 : [3.2202897637586845, 5.880492957780065, 5.8804944200433535, 10.686025216519425, 10.69391907743575, 10.693989757202822, 12.317666113412542, 12.320210206543557, 13.360871916890284, 13.422816823195596, 13.42565047944317, 14.216103594821192]
11 : [3.220289763758685, 5.880492957779892, 5.880494420042592, 10.686006011027187, 10.69391741457624, 10.693986570955898, 12.317664492530746, 12.318659200540916, 13.358849290077444, 13.422344037999789, 13.423640051553738, 14.21552028324461]
12 : [3.2202897637586854, 5.880492957779903, 5.88049442004258, 10.686001612142082, 10.693917021614238, 10.69398585920967, 12.317664119923707, 12.318123649766568, 13.35793068593903, 13.4221365401517, 13.422748126002961, 14.215203523693638]
13 : [3.220289763758684, 5.880492957779901, 5.88049442004258, 10.686000605540046, 10.693916931075526, 10.693985697528339, 12.317664014898458, 12.317939017877379, 13.35751846821583, 13.422044242233461, 13.422347224164689, 14.215030357743023]
14 : [3.2202897637586836, 5.88049295777993, 5.880494420042571, 10.686000374957876, 10.69391691031485, 10.693985660606163, 12.317663977188824, 12.317875366088522, 13.357334382372857, 13.422002931161972, 13.422165817662501, 14.21493536325286]
15 : [3.2202897637586827, 5.880492957779908, 5.880494420042718, 10.686000322134197, 10.693916905555612, 10.693985652164555, 12.317663962997816, 12.317853422209481, 13.357252355363546, 13.421984311165387, 13.422083607176628, 14.21488316409479]
16 : [3.220289763758686, 5.8804929577799445, 5.8804944200422105, 10.686000309999313, 10.693916904463954, 10.693985650228601, 12.317663957904443, 12.317845847600688, 13.35721580739596, 13.421975885139316, 13.422046240578887, 14.21485440631397]
17 : [3.2202897637586814, 5.880492957780259, 5.880494420042958, 10.686000307201386, 10.69391690421217, 10.693985649782805, 12.317663956131295, 12.317843215146942, 13.35719940334385, 13.42197203894705, 13.422029197812318, 14.214838473942988]
18 : [3.2202897637586827, 5.880492957780043, 5.880494420042497, 10.686000306526424, 10.693916904147438, 10.693985649675579, 12.317663955498062, 12.317842216396107, 13.35719144281174, 13.421969767888832, 13.422020792755646, 14.214829494070095]
19 : [3.220289763758687, 5.8804929577800396, 5.880494420042625, 10.68600030639704, 10.693916904138748, 10.693985649654921, 12.31766395529176, 12.317841937109678, 13.357188265808071, 13.421969196469503, 13.422017174642297, 14.214824609575201]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.220289763758687
5.8804929577800396
5.880494420042625
10.68600030639704
10.693916904138748
10.693985649654921
12.31766395529176
12.317841937109678
13.357188265808071
13.421969196469503
13.422017174642297
14.214824609575201
[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);
[ ]:

[ ]: