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)
/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.3
  warnings.warn(f"A NumPy version >={np_minversion} and <{np_maxversion}"
0 : [6.503661638503196, 18.28000559369274, 32.400955190973164, 44.74211013521184, 55.66625104305, 67.56144913171728, 71.85226456144288, 82.45261447404701, 97.5414714716846, 104.20791440265421, 119.25461398998388, 129.00962192061905]
1 : [3.2326201032780078, 6.0151252723548865, 6.324162915050844, 11.656799045091251, 12.517366573955506, 13.889682738826073, 14.294353803531983, 15.617357216861217, 17.564164297124115, 21.911459497937773, 26.481037175286403, 27.483341390005872]
2 : [3.220503347683658, 5.885411818036616, 5.90262303649664, 10.803247858828053, 11.111948042316046, 11.209812137741086, 12.517607876335523, 13.229455378616068, 14.3696165795129, 14.948471531707883, 15.80014401098003, 20.315857789415883]
3 : [3.2202888617094514, 5.880790648651129, 5.881426598481317, 10.714002290455909, 10.766498296193681, 10.878783445822542, 12.366527319686988, 12.688898772128908, 13.690746873885258, 13.721687127679086, 14.029886469489933, 18.18108094157683]
4 : [3.2202839873351463, 5.880508093129893, 5.880532268530896, 10.698241349675158, 10.702746554419676, 10.76281694137739, 12.33145777559133, 12.485118352638246, 13.453049432807363, 13.509949313861357, 13.706220372864562, 16.898411414843977]
5 : [3.2202838727289125, 5.8804867249811945, 5.880492065885029, 10.690534467985483, 10.694911834211219, 10.71729137156061, 12.321829566483077, 12.398794181809311, 13.389546622496297, 13.45213108385587, 13.585960290581392, 15.939362641672425]
6 : [3.2202838702452636, 5.880484860739922, 5.8804905020718525, 10.687741310963828, 10.694166748335096, 10.700719606524807, 12.318960198358939, 12.356656803997575, 13.369680205755012, 13.433174459106155, 13.52297715037836, 15.257581550910478]
7 : [3.2202838701954106, 5.880484738540208, 5.880490426990247, 10.686622972664448, 10.693998088206504, 10.695684550617573, 12.318056494315819, 12.335059518181838, 13.362274276301813, 13.426410200837351, 13.482628077127172, 14.815662119756935]
8 : [3.2202838701944856, 5.880484731105718, 5.880490423203137, 10.686175995680907, 10.69395775088781, 10.694366167658064, 12.317759245830757, 12.32482398017813, 13.359253513508094, 13.423793347725464, 13.456168925442133, 14.55207130891542]
9 : [3.220283870194469, 5.880484730673671, 5.88049042300778, 10.686042393372992, 10.693940614711014, 10.694041984405645, 12.317658269508922, 12.320419793704009, 13.35796370602923, 13.422713305838716, 13.439932743369006, 14.401613901907893]
10 : [3.2202838701944656, 5.880484730649438, 5.8804904229975765, 10.686008320345087, 10.693918504699578, 10.693980628723882, 12.31762067509824, 12.31868545889796, 13.357400735717823, 13.422248076286305, 13.430896559382983, 14.317764875453177]
11 : [3.2202838701944674, 5.880484730648127, 5.8804904229970605, 10.686000138343053, 10.693907238097516, 10.69397197082405, 12.317600284138411, 12.318044956006553, 13.357152467731876, 13.42204281651889, 13.426233947965704, 14.271387955798021]
12 : [3.220283870194466, 5.880484730648042, 5.8804904229970205, 10.685998223458785, 10.693904207606149, 10.693970330803573, 12.317581475885602, 12.31782555964091, 13.357042691169129, 13.421951041274585, 13.423953958796032, 14.245855706253863]
13 : [3.220283870194468, 5.880484730648021, 5.880490422997009, 10.685997779518887, 10.693903484688981, 10.693969969389412, 12.31756640165118, 12.31775627541096, 13.356994055614338, 13.421909642390752, 13.422873481289885, 14.231814840479199]
14 : [3.2202838701944665, 5.880484730648055, 5.880490422996988, 10.685997677086638, 10.693903316722658, 10.693969886972988, 12.317558577885235, 12.317734530844211, 13.356972498397626, 13.421890753905123, 13.42237134116572, 14.224102609309758]
15 : [3.22028387019447, 5.88048473064805, 5.880490422996971, 10.685997653478006, 10.693903277974753, 10.693969868049033, 12.317555416705574, 12.317727377662642, 13.356962912442953, 13.421881839237013, 13.422140496732185, 14.219868262886944]
16 : [3.220283870194471, 5.8804847306481625, 5.880490422997122, 10.685997648040832, 10.693903269057177, 10.693969863700076, 12.31755425837452, 12.317724948995185, 13.356958645294505, 13.421877330315679, 13.422035333332488, 14.21754569739388]
17 : [3.220283870194472, 5.880484730648496, 5.880490422996954, 10.685997646785971, 10.693903266999657, 10.69396986269786, 12.31755384732725, 12.317724110841205, 13.35695674128071, 13.421874862444403, 13.42198764591848, 14.216267058883215]
18 : [3.2202838701944687, 5.880484730648013, 5.880490422996572, 10.685997646494634, 10.693903266522065, 10.693969862465265, 12.317553702341911, 12.317723817872698, 13.356955889714076, 13.421873484639043, 13.42196584538761, 14.21555272584488]
19 : [3.220283870194471, 5.880484730648379, 5.880490422996613, 10.685997646429774, 10.693903266416008, 10.69396986241319, 12.317553652778441, 12.31772371867834, 13.356955462831765, 13.421872736130078, 13.421956336327534, 14.215169375547026]
[6]:
print ("Eigenvalues")
for lam in evals:
    print (lam)
Eigenvalues
3.220283870194471
5.880484730648379
5.880490422996613
10.685997646429774
10.693903266416008
10.69396986241319
12.317553652778441
12.31772371867834
13.356955462831765
13.421872736130078
13.421956336327534
14.215169375547026
[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);
[ ]:

[ ]: