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.004801904878959065
iteration 1 error = 0.0025206764324378827
iteration 2 error = 0.0020551864069441655
iteration 3 error = 0.0021126409485256493
iteration 4 error = 0.0011634336768156861
iteration 5 error = 0.0008257478480355326
iteration 6 error = 0.0006489422159427328
iteration 7 error = 0.00044904321604994144
iteration 8 error = 0.00036444673646794705
iteration 9 error = 0.000220361400021187
iteration 10 error = 0.00017327061351427568
iteration 11 error = 0.00010853931819067255
iteration 12 error = 7.085992463910578e-05
iteration 13 error = 5.1703327118476115e-05
iteration 14 error = 3.2741974697312774e-05
iteration 15 error = 2.3498971820604713e-05
iteration 16 error = 1.4876027757753357e-05
iteration 17 error = 1.0121351561770796e-05
iteration 18 error = 7.011616090982956e-06
iteration 19 error = 4.674823034053791e-06
iteration 20 error = 3.096491172463652e-06
iteration 21 error = 2.1485821580399096e-06
iteration 22 error = 1.4413996112764504e-06
iteration 23 error = 9.760512003510296e-07
iteration 24 error = 7.015646876660583e-07
iteration 25 error = 6.919151630725699e-07
iteration 26 error = 3.831029835023217e-07
iteration 27 error = 2.3506733663218052e-07
iteration 28 error = 1.6521040030133825e-07
iteration 29 error = 1.0719777249250855e-07
iteration 30 error = 7.12659703694528e-08
iteration 31 error = 5.1028363116894886e-08
iteration 32 error = 3.177702387231607e-08
iteration 33 error = 2.222551259876341e-08
iteration 34 error = 1.4397290742277198e-08
iteration 35 error = 9.498340038055631e-09
iteration 36 error = 6.209554707548018e-09
iteration 37 error = 4.196986667356919e-09
iteration 38 error = 2.806960988842398e-09
iteration 39 error = 1.829402396959425e-09
iteration 40 error = 1.197387066112526e-09
iteration 41 error = 7.552591638846266e-10
iteration 42 error = 5.353415322673602e-10
iteration 43 error = 3.655741141342465e-10
iteration 44 error = 3.5446274568648325e-10
iteration 45 error = 2.3071937615826502e-10
iteration 46 error = 1.4585866157161644e-10
iteration 47 error = 1.2198808681967328e-10
iteration 48 error = 7.644523982797499e-11
iteration 49 error = 4.913754809618186e-11
iteration 50 error = 2.950780796133023e-11
iteration 51 error = 2.2553257526734677e-11
iteration 52 error = 1.3711701450174411e-11
iteration 53 error = 9.064501632008665e-12
iteration 54 error = 5.948248008333114e-12
iteration 55 error = 4.098994595033413e-12
iteration 56 error = 2.5050716761575012e-12
iteration 57 error = 1.6349370686554411e-12
iteration 58 error = 1.1707784881759781e-12
iteration 59 error = 7.566339537526967e-13
iteration 60 error = 4.924144945024119e-13
iteration 61 error = 3.2244468366540966e-13
iteration 62 error = 2.2283319190038358e-13
iteration 63 error = 1.3646053755584538e-13
iteration 64 error = 1.2933114815608787e-13
iteration 65 error = 8.749025217961355e-14
iteration 66 error = 5.0259238660638144e-14
iteration 67 error = 3.431385544822773e-14
iteration 68 error = 2.254255160153578e-14
iteration 69 error = 1.4854696963372665e-14
iteration 70 error = 9.320336263306305e-15
iteration 71 error = 5.975183767711616e-15
iteration 72 error = 3.714099248639716e-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)