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

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.00477075126876544
iteration 1 error = 0.0024813274399324453
iteration 2 error = 0.001891809576612191
iteration 3 error = 0.0012587090378891413
iteration 4 error = 0.001136453623318983
iteration 5 error = 0.0008195457478849468
iteration 6 error = 0.0006060825336657854
iteration 7 error = 0.0004534128861287399
iteration 8 error = 0.00035785131847072886
iteration 9 error = 0.00023195184860900933
iteration 10 error = 0.00015491751841894291
iteration 11 error = 9.808156355275184e-05
iteration 12 error = 7.798886269681482e-05
iteration 13 error = 4.431269578808669e-05
iteration 14 error = 3.282881496334983e-05
iteration 15 error = 2.292008755437698e-05
iteration 16 error = 1.5012379279908444e-05
iteration 17 error = 1.0290586276555379e-05
iteration 18 error = 6.797403084573727e-06
iteration 19 error = 5.32408043420282e-06
iteration 20 error = 3.0935107810297475e-06
iteration 21 error = 2.2732501299784854e-06
iteration 22 error = 1.513919531735276e-06
iteration 23 error = 9.97680981027407e-07
iteration 24 error = 6.972869052375712e-07
iteration 25 error = 4.710001624204149e-07
iteration 26 error = 3.089882446018695e-07
iteration 27 error = 2.1201213065586627e-07
iteration 28 error = 1.4360554950973686e-07
iteration 29 error = 9.438532197533975e-08
iteration 30 error = 6.136548681935135e-08
iteration 31 error = 4.355524506643756e-08
iteration 32 error = 2.9553593068140377e-08
iteration 33 error = 1.9220272899373628e-08
iteration 34 error = 1.2599807004166964e-08
iteration 35 error = 8.494750167552192e-09
iteration 36 error = 6.286103375548955e-09
iteration 37 error = 5.183946170794216e-09
iteration 38 error = 3.4157067658145383e-09
iteration 39 error = 2.162227236378485e-09
iteration 40 error = 1.3842597985048054e-09
iteration 41 error = 9.570035371748851e-10
iteration 42 error = 6.402042970245908e-10
iteration 43 error = 4.4182326930262824e-10
iteration 44 error = 2.926001367032175e-10
iteration 45 error = 1.9234435803006636e-10
iteration 46 error = 1.39041462879003e-10
iteration 47 error = 1.0930029811166763e-10
iteration 48 error = 6.993359899888951e-11
iteration 49 error = 4.6262890739911466e-11
iteration 50 error = 2.7984555829562866e-11
iteration 51 error = 2.071495355696642e-11
iteration 52 error = 1.3207477509240725e-11
iteration 53 error = 8.907476692979874e-12
iteration 54 error = 5.392263312210664e-12
iteration 55 error = 3.829585507045024e-12
iteration 56 error = 2.3572446496231288e-12
iteration 57 error = 1.56110891223033e-12
iteration 58 error = 1.0214140192811456e-12
iteration 59 error = 6.905122516714876e-13
iteration 60 error = 4.408306174492103e-13
iteration 61 error = 2.908526394928918e-13
iteration 62 error = 1.8708288446484247e-13
iteration 63 error = 1.18488294454598e-13
iteration 64 error = 9.7948983305966e-14
iteration 65 error = 7.306076249242146e-14
iteration 66 error = 4.5968971499907486e-14
iteration 67 error = 2.846450916525336e-14
iteration 68 error = 1.8543026469140923e-14
iteration 69 error = 1.2782122277042395e-14
iteration 70 error = 7.90163796471695e-15
iteration 71 error = 5.449188038664005e-15
iteration 72 error = 3.917698008089792e-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)