This page was generated from unit-2.4-Maxwell/Maxwell.ipynb.

2.4 Maxwell’s Equations

[Peter Monk: "Finite Elements for Maxwell’s Equations"]

Magnetostatic field generated by a permanent magnet

magnetic flux \(B\), magnetic field \(H\), given magnetization \(M\):

\[\DeclareMathOperator{\Grad}{grad} \DeclareMathOperator{\Curl}{curl} \DeclareMathOperator{\Div}{div} B = \mu (H + M), \quad \Div B = 0, \quad \Curl H = 0\]

Introducing a vector-potential \(A\) such that \(B = \Curl A\), and putting equations together we get

\[\Curl \mu^{-1} \Curl A = \Curl M\]

In weak form: Find \(A \in H(\Curl)\) such that

\[\int \mu^{-1} \Curl A \Curl v = \int M \Curl v \qquad \forall \, v \in H(\Curl)\]

Usually, the permeability \(\mu\) is given as \(\mu = \mu_r \mu_0\), with \(\mu_0 = 4 \pi 10^{-7}\) the permeability of vacuum.

[1]:
from ngsolve import *
from netgen.csg import *
import netgen.gui

Geometric model and meshing of a bar magnet:

[2]:
box = OrthoBrick(Pnt(-3,-3,-3),Pnt(3,3,3)).bc("outer")
magnet = Cylinder(Pnt(-1,0,0),Pnt(1,0,0), 0.3) * OrthoBrick(Pnt(-1,-3,-3),Pnt(1,3,3))
air = box - magnet

geo = CSGeometry()
geo.Add (air.mat("air"), transparent=True)
geo.Add (magnet.mat("magnet").maxh(1), col=(0.3,0.3,0.1))
geo.Draw()
mesh = Mesh(geo.GenerateMesh(maxh=2, curvaturesafety=1))
mesh.Curve(3)
[3]:
mesh.GetMaterials(), mesh.GetBoundaries()
[3]:
(('air', 'magnet'),
 ('outer',
  'outer',
  'outer',
  'outer',
  'outer',
  'outer',
  'default',
  'default',
  'default'))

Define space, forms and preconditioner.

  • To obtain a regular system matrix, we regularize by adding a very small \(L_2\) term.

  • We solve magnetostatics, so we can gauge by adding and arbitrary gradient field. A cheap possibility is to delete all basis-functions which are gradients (flag 'nograds')

[4]:
fes = HCurl(mesh, order=3, dirichlet="outer", nograds=True)
print ("ndof =", fes.ndof)
u,v = fes.TnT()

from math import pi
mu0 = 4*pi*1e-7
mur = CoefficientFunction( [1000 if mat== "magnet" else 1
                            for mat in mesh.GetMaterials()])

a = BilinearForm(fes)
a += 1/(mu0*mur)*curl(u)*curl(v)*dx + 1e-8/(mu0*mur)*u*v*dx
c = Preconditioner(a, "bddc")

f = LinearForm(fes)
mag = CoefficientFunction((1,0,0)) * \
    CoefficientFunction( [1 if mat == "magnet" else 0 for mat in mesh.GetMaterials()])
f += SymbolicLFI(mag*curl(v), definedon=mesh.Materials("magnet"))
ndof = 23791

Assemble system and setup preconditioner using task-parallelization:

[5]:
with TaskManager():
    a.Assemble()
    f.Assemble()

Finally, declare GridFunction and solve by preconditioned CG iteration:

[6]:
gfu = GridFunction(fes)
with TaskManager():
    solvers.CG(sol=gfu.vec, rhs=f.vec, mat=a.mat, pre=c.mat)
iteration 0 error = 0.004855037018856693
iteration 1 error = 0.0027293684352409264
iteration 2 error = 0.0019178169752064082
iteration 3 error = 0.0011868200578276729
iteration 4 error = 0.0008381444433524454
iteration 5 error = 0.0005352958317220561
iteration 6 error = 0.00041878588073819854
iteration 7 error = 0.00024966729220306623
iteration 8 error = 0.00014794441567466277
iteration 9 error = 9.829281774119915e-05
iteration 10 error = 6.796526225372915e-05
iteration 11 error = 3.914655734049878e-05
iteration 12 error = 2.2824616831835196e-05
iteration 13 error = 1.4024099850047262e-05
iteration 14 error = 9.1099147199686e-06
iteration 15 error = 5.472491313755574e-06
iteration 16 error = 3.3322895627949107e-06
iteration 17 error = 1.9728342120026605e-06
iteration 18 error = 1.1871550342888096e-06
iteration 19 error = 7.339626800917806e-07
iteration 20 error = 4.3623176671614723e-07
iteration 21 error = 2.566836359188603e-07
iteration 22 error = 1.607474278874305e-07
iteration 23 error = 1.0206447162818923e-07
iteration 24 error = 6.217499592128519e-08
iteration 25 error = 3.72653787722966e-08
iteration 26 error = 2.181501210285165e-08
iteration 27 error = 1.2445629160114031e-08
iteration 28 error = 7.3237410812179445e-09
iteration 29 error = 4.85304570680509e-09
iteration 30 error = 3.1476230725430698e-09
iteration 31 error = 1.7144938912260811e-09
iteration 32 error = 1.0204511653609662e-09
iteration 33 error = 6.285157101665242e-10
iteration 34 error = 3.8201081097035006e-10
iteration 35 error = 2.2615096781894397e-10
iteration 36 error = 1.3334272929782343e-10
iteration 37 error = 7.592566297542565e-11
iteration 38 error = 4.462527879173254e-11
iteration 39 error = 2.5689939202764936e-11
iteration 40 error = 1.538049632156891e-11
iteration 41 error = 9.117313517763468e-12
iteration 42 error = 5.339750141809667e-12
iteration 43 error = 3.270507144264021e-12
iteration 44 error = 2.086813020223854e-12
iteration 45 error = 1.5020857314956551e-12
iteration 46 error = 9.294642081130611e-13
iteration 47 error = 5.434729221858305e-13
iteration 48 error = 3.062070108432163e-13
iteration 49 error = 2.1357473876334865e-13
iteration 50 error = 1.420593753317279e-13
iteration 51 error = 7.98333396290258e-14
iteration 52 error = 4.425430557291156e-14
iteration 53 error = 2.595628779325309e-14
iteration 54 error = 1.532069846262541e-14
iteration 55 error = 8.925511330456275e-15
iteration 56 error = 6.555378743393437e-15
iteration 57 error = 3.9495062638444604e-15
[7]:
# the vector potential is not supposed to look nice
Draw (gfu, mesh, "vector-potential", draw_surf=False)

Draw (curl(gfu), mesh, "B-field", draw_surf=False)
Draw (1/(mu0*mur)*curl(gfu)-mag, mesh, "H-field", draw_surf=False)