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

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.00481037704725234
iteration 1 error = 0.002376049542145329
iteration 2 error = 0.002009132992263431
iteration 3 error = 0.0018347748421582207
iteration 4 error = 0.0011702983943112164
iteration 5 error = 0.0007994173970521496
iteration 6 error = 0.0006047754213553879
iteration 7 error = 0.000448473017956756
iteration 8 error = 0.0003371739560894319
iteration 9 error = 0.0002183015509288613
iteration 10 error = 0.00016470550887057354
iteration 11 error = 0.00010632701022515559
iteration 12 error = 6.970620403330876e-05
iteration 13 error = 4.956032015510294e-05
iteration 14 error = 3.212912718883845e-05
iteration 15 error = 2.250030288281255e-05
iteration 16 error = 1.4624085914637864e-05
iteration 17 error = 9.992986502221412e-06
iteration 18 error = 6.8146244248731805e-06
iteration 19 error = 4.660924327744275e-06
iteration 20 error = 2.9828720383834402e-06
iteration 21 error = 2.0715246602185236e-06
iteration 22 error = 1.4187981026274697e-06
iteration 23 error = 9.377814057586366e-07
iteration 24 error = 6.484364795299737e-07
iteration 25 error = 6.751264440067065e-07
iteration 26 error = 4.0200632951431836e-07
iteration 27 error = 2.4550117026716015e-07
iteration 28 error = 1.731366175408803e-07
iteration 29 error = 1.1155191218858394e-07
iteration 30 error = 7.144838753127387e-08
iteration 31 error = 4.782664516594915e-08
iteration 32 error = 3.116814110480835e-08
iteration 33 error = 2.0859579896309427e-08
iteration 34 error = 1.4126543325276659e-08
iteration 35 error = 8.73706845793194e-09
iteration 36 error = 5.818889511146841e-09
iteration 37 error = 3.964936383845389e-09
iteration 38 error = 2.724584550353857e-09
iteration 39 error = 1.6443926888048447e-09
iteration 40 error = 1.090458361946667e-09
iteration 41 error = 7.480108578806244e-10
iteration 42 error = 4.966951136791269e-10
iteration 43 error = 3.4890829754479225e-10
iteration 44 error = 3.071509100625328e-10
iteration 45 error = 2.3814630193199996e-10
iteration 46 error = 1.7016410627015397e-10
iteration 47 error = 1.1429795921698263e-10
iteration 48 error = 6.7747844078969e-11
iteration 49 error = 4.8104844218206994e-11
iteration 50 error = 3.0259251572191896e-11
iteration 51 error = 2.1363406067840313e-11
iteration 52 error = 1.403310784542784e-11
iteration 53 error = 8.783192511318458e-12
iteration 54 error = 5.976802992245573e-12
iteration 55 error = 4.048720285273533e-12
iteration 56 error = 2.6176677617093826e-12
iteration 57 error = 1.675822735195238e-12
iteration 58 error = 1.1760797136460367e-12
iteration 59 error = 7.80804128814466e-13
iteration 60 error = 4.882351973503057e-13
iteration 61 error = 3.200465034524917e-13
iteration 62 error = 2.2386326136235992e-13
iteration 63 error = 1.4443122758147973e-13
iteration 64 error = 1.2980121118747188e-13
iteration 65 error = 9.370258834429704e-14
iteration 66 error = 5.2811724613658416e-14
iteration 67 error = 3.6114198393148824e-14
iteration 68 error = 2.3896271076565704e-14
iteration 69 error = 1.6247213611635304e-14
iteration 70 error = 1.0620749421115862e-14
iteration 71 error = 6.5470582673904535e-15
iteration 72 error = 4.250449540262065e-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)