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 = mesh.MaterialCF({"magnet" : 1000}, default=1)

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 = mesh.MaterialCF({"magnet" : (1,0,0)}, default=(0,0,0))
f += mag*curl(v) * dx("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.0025206764324378845
iteration 2 error = 0.002055186406944164
iteration 3 error = 0.0021126409485256498
iteration 4 error = 0.0011634336768156844
iteration 5 error = 0.0008257478480355314
iteration 6 error = 0.0006489422159427318
iteration 7 error = 0.00044904321604994106
iteration 8 error = 0.0003644467364679465
iteration 9 error = 0.0002203614000211869
iteration 10 error = 0.00017327061351427568
iteration 11 error = 0.00010853931819067254
iteration 12 error = 7.085992463910571e-05
iteration 13 error = 5.1703327118476074e-05
iteration 14 error = 3.274197469731272e-05
iteration 15 error = 2.3498971820604595e-05
iteration 16 error = 1.487602775775188e-05
iteration 17 error = 1.0121351561733626e-05
iteration 18 error = 7.011616090030651e-06
iteration 19 error = 4.674823014212931e-06
iteration 20 error = 3.0964907292966043e-06
iteration 21 error = 2.148573019610739e-06
iteration 22 error = 1.4412349059811364e-06
iteration 23 error = 9.725910194891257e-07
iteration 24 error = 6.048980912364284e-07
iteration 25 error = 4.479485932080253e-07
iteration 26 error = 4.6668061666087147e-07
iteration 27 error = 2.440893455057485e-07
iteration 28 error = 1.6535541752363835e-07
iteration 29 error = 1.0720033565637214e-07
iteration 30 error = 7.126601233023247e-08
iteration 31 error = 5.1028364125004926e-08
iteration 32 error = 3.1777023888268163e-08
iteration 33 error = 2.222551259904319e-08
iteration 34 error = 1.4397290742253317e-08
iteration 35 error = 9.498340037812984e-09
iteration 36 error = 6.2095547060367535e-09
iteration 37 error = 4.196986654057495e-09
iteration 38 error = 2.806960930577301e-09
iteration 39 error = 1.8294019698344646e-09
iteration 40 error = 1.197384148155887e-09
iteration 41 error = 7.552196796748971e-10
iteration 42 error = 5.347427686855083e-10
iteration 43 error = 3.5556446643883766e-10
iteration 44 error = 2.3694695243882644e-10
iteration 45 error = 1.9475844116900805e-10
iteration 46 error = 1.6238656815194184e-10
iteration 47 error = 1.0364427392110557e-10
iteration 48 error = 7.717981861273891e-11
iteration 49 error = 5.0765220346469365e-11
iteration 50 error = 2.95580284862959e-11
iteration 51 error = 2.2558395947271914e-11
iteration 52 error = 1.371224446567973e-11
iteration 53 error = 9.064527759564304e-12
iteration 54 error = 5.948249257829294e-12
iteration 55 error = 4.0989946902029855e-12
iteration 56 error = 2.5050716802785445e-12
iteration 57 error = 1.6349370688098089e-12
iteration 58 error = 1.1707784881198017e-12
iteration 59 error = 7.566339519141841e-13
iteration 60 error = 4.924144395113912e-13
iteration 61 error = 3.2244288023348155e-13
iteration 62 error = 2.2278510884364117e-13
iteration 63 error = 1.3500312474197775e-13
iteration 64 error = 8.651582875526798e-14
iteration 65 error = 8.238246803471227e-14
iteration 66 error = 5.679653422998336e-14
iteration 67 error = 3.442486185193374e-14
iteration 68 error = 2.254464297531787e-14
iteration 69 error = 1.485541248943598e-14
iteration 70 error = 9.323211241173782e-15
iteration 71 error = 5.98827851804261e-15
iteration 72 error = 3.796346277776626e-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)