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

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)
it =  0  err =  0.004840564319042956
it =  1  err =  0.002877746735453422
it =  2  err =  0.0018540826001404458
it =  3  err =  0.001201924591628051
it =  4  err =  0.0008069147490216403
it =  5  err =  0.0005450296846235297
it =  6  err =  0.0003675266044562467
it =  7  err =  0.0002434589451927929
it =  8  err =  0.0001435746661404234
it =  9  err =  9.620650815911942e-05
it =  10  err =  5.87331405801642e-05
it =  11  err =  3.586599439923349e-05
it =  12  err =  2.236083879867455e-05
it =  13  err =  1.3049751993777967e-05
it =  14  err =  8.43490296134633e-06
it =  15  err =  5.399805488388646e-06
it =  16  err =  3.0315968419612513e-06
it =  17  err =  1.826457133757794e-06
it =  18  err =  1.14341432597726e-06
it =  19  err =  7.016189041097574e-07
it =  20  err =  4.542643111994101e-07
it =  21  err =  2.6888162144515167e-07
it =  22  err =  1.5855127448291939e-07
it =  23  err =  9.730336619012949e-08
it =  24  err =  5.622498588543311e-08
it =  25  err =  3.3753841921969977e-08
it =  26  err =  2.050695590605847e-08
it =  27  err =  1.2573120984768046e-08
it =  28  err =  7.310342124602579e-09
it =  29  err =  4.290376411824819e-09
it =  30  err =  2.469419254593805e-09
it =  31  err =  1.4316972534344513e-09
it =  32  err =  8.77605482241475e-10
it =  33  err =  5.202538144553532e-10
it =  34  err =  3.177532555089666e-10
it =  35  err =  1.888099034497334e-10
it =  36  err =  1.110042908511075e-10
it =  37  err =  6.568461994601609e-11
it =  38  err =  3.867488752153068e-11
it =  39  err =  2.288029494354035e-11
it =  40  err =  1.3787282352059778e-11
it =  41  err =  8.551870860782309e-12
it =  42  err =  4.9089432018070234e-12
it =  43  err =  2.891176286251989e-12
it =  44  err =  1.7320555921750252e-12
it =  45  err =  1.0439723245773804e-12
it =  46  err =  6.985971343990278e-13
it =  47  err =  4.5870462328153e-13
it =  48  err =  2.8229382656143784e-13
it =  49  err =  1.6508749110510983e-13
it =  50  err =  9.667673584261733e-14
it =  51  err =  5.87547533037912e-14
it =  52  err =  3.487259126492939e-14
it =  53  err =  2.1353900683901444e-14
it =  54  err =  1.2821170764192051e-14
it =  55  err =  7.280305958730234e-15
it =  56  err =  4.271331523280645e-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)