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 = 25073

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.004801904878959067
iteration 1 error = 0.0025206764324378892
iteration 2 error = 0.002055186406944168
iteration 3 error = 0.0021126409485256524
iteration 4 error = 0.0011634336768156861
iteration 5 error = 0.000825747848035533
iteration 6 error = 0.000648942215942733
iteration 7 error = 0.00044904321604994144
iteration 8 error = 0.00036444673646794694
iteration 9 error = 0.00022036140002118737
iteration 10 error = 0.00017327061351427595
iteration 11 error = 0.00010853931819067265
iteration 12 error = 7.085992463910584e-05
iteration 13 error = 5.17033271184761e-05
iteration 14 error = 3.274197469731269e-05
iteration 15 error = 2.349897182060458e-05
iteration 16 error = 1.4876027757752212e-05
iteration 17 error = 1.0121351561742535e-05
iteration 18 error = 7.011616090259097e-06
iteration 19 error = 4.674823018972441e-06
iteration 20 error = 3.096490835605277e-06
iteration 21 error = 2.148575211781618e-06
iteration 22 error = 1.4412744203830105e-06
iteration 23 error = 9.734237027499186e-07
iteration 24 error = 6.309688388522009e-07
iteration 25 error = 6.82305942602368e-07
iteration 26 error = 3.9813693169625365e-07
iteration 27 error = 2.353174492006475e-07
iteration 28 error = 1.6521415875422527e-07
iteration 29 error = 1.0719783880228633e-07
iteration 30 error = 7.126597145494043e-08
iteration 31 error = 5.1028363142973635e-08
iteration 32 error = 3.177702387272794e-08
iteration 33 error = 2.222551259876551e-08
iteration 34 error = 1.4397290742256472e-08
iteration 35 error = 9.498340037887657e-09
iteration 36 error = 6.209554706501953e-09
iteration 37 error = 4.196986658147114e-09
iteration 38 error = 2.8069609484027617e-09
iteration 39 error = 1.8294020978695653e-09
iteration 40 error = 1.197384953521398e-09
iteration 41 error = 7.552282277299711e-10
iteration 42 error = 5.348179674542248e-10
iteration 43 error = 3.5622103139960756e-10
iteration 44 error = 2.4743153321277744e-10
iteration 45 error = 2.541568636088183e-10
iteration 46 error = 1.4688766080773762e-10
iteration 47 error = 1.1397502832424828e-10
iteration 48 error = 7.799716393733422e-11
iteration 49 error = 4.953872780853822e-11
iteration 50 error = 2.9518968348924264e-11
iteration 51 error = 2.2554389814794794e-11
iteration 52 error = 1.3711820877416213e-11
iteration 53 error = 9.064507377621764e-12
iteration 54 error = 5.948248283103065e-12
iteration 55 error = 4.098994615961808e-12
iteration 56 error = 2.50507167706437e-12
iteration 57 error = 1.634937068690429e-12
iteration 58 error = 1.1707784881029839e-12
iteration 59 error = 7.566339518917745e-13
iteration 60 error = 4.924144391564219e-13
iteration 61 error = 3.22442869377653e-13
iteration 62 error = 2.2278482147429225e-13
iteration 63 error = 1.3499433235037164e-13
iteration 64 error = 8.609388815911895e-14
iteration 65 error = 7.439089038193713e-14
iteration 66 error = 5.982447215560622e-14
iteration 67 error = 3.45035616187308e-14
iteration 68 error = 2.254589247383458e-14
iteration 69 error = 1.4855186472458122e-14
iteration 70 error = 9.322190242660997e-15
iteration 71 error = 5.983640880361689e-15
iteration 72 error = 3.767939076851137e-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)