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 ngsolve.webgui import Draw
from netgen.csg import *

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)
[2]:
<ngsolve.comp.Mesh at 0x7f2f1ba58e50>
[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 = 25068

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)
CG iteration 1, residual = 0.004832539962849043
CG iteration 2, residual = 0.002880198465291355
CG iteration 3, residual = 0.0022636840717767852
CG iteration 4, residual = 0.0018596783443640293
CG iteration 5, residual = 0.0013706814431410216
CG iteration 6, residual = 0.0009629315328718029
CG iteration 7, residual = 0.0006784248172513836
CG iteration 8, residual = 0.0004537347998869367
CG iteration 9, residual = 0.0003373984523584697
CG iteration 10, residual = 0.00026636366966670505
CG iteration 11, residual = 0.00015785865716709974
CG iteration 12, residual = 0.00010596525810194435
CG iteration 13, residual = 7.440611757876558e-05
CG iteration 14, residual = 5.0217762320984635e-05
CG iteration 15, residual = 3.565457135377311e-05
CG iteration 16, residual = 2.1652775592701065e-05
CG iteration 17, residual = 1.6241088584150902e-05
CG iteration 18, residual = 1.0173867299009934e-05
CG iteration 19, residual = 7.24779785800793e-06
CG iteration 20, residual = 4.782087689360765e-06
CG iteration 21, residual = 3.29777468847513e-06
CG iteration 22, residual = 2.3147389159099555e-06
CG iteration 23, residual = 1.4401085715934553e-06
CG iteration 24, residual = 1.0155349944492385e-06
CG iteration 25, residual = 6.917281013524945e-07
CG iteration 26, residual = 5.275355831372497e-07
CG iteration 27, residual = 4.974775101189182e-07
CG iteration 28, residual = 2.8614865085129467e-07
CG iteration 29, residual = 1.837790684575086e-07
CG iteration 30, residual = 1.2576183927106167e-07
CG iteration 31, residual = 8.212389911927664e-08
CG iteration 32, residual = 5.655902201242402e-08
CG iteration 33, residual = 3.7186008444538384e-08
CG iteration 34, residual = 2.395621703424371e-08
CG iteration 35, residual = 1.6498257419582005e-08
CG iteration 36, residual = 1.0759738081536739e-08
CG iteration 37, residual = 7.165063591525583e-09
CG iteration 38, residual = 4.497597049293491e-09
CG iteration 39, residual = 2.915678167278719e-09
CG iteration 40, residual = 1.9460018382445743e-09
CG iteration 41, residual = 1.4016014562368145e-09
CG iteration 42, residual = 8.455158583231661e-10
CG iteration 43, residual = 5.888050459800812e-10
CG iteration 44, residual = 3.72468194199161e-10
CG iteration 45, residual = 2.64269025466415e-10
CG iteration 46, residual = 2.247041813931647e-10
CG iteration 47, residual = 1.7287317345711777e-10
CG iteration 48, residual = 1.4007896639773846e-10
CG iteration 49, residual = 8.512279473019448e-11
CG iteration 50, residual = 5.4099565834011057e-11
CG iteration 51, residual = 3.64277886156205e-11
CG iteration 52, residual = 2.3642432722910823e-11
CG iteration 53, residual = 1.6320993211923992e-11
CG iteration 54, residual = 1.0677525909959887e-11
CG iteration 55, residual = 7.1559440691513135e-12
CG iteration 56, residual = 4.796320654071837e-12
CG iteration 57, residual = 3.098941053682878e-12
CG iteration 58, residual = 2.0500214928599125e-12
CG iteration 59, residual = 1.3892674854168662e-12
CG iteration 60, residual = 8.940758256903124e-13
CG iteration 61, residual = 5.967413972297896e-13
CG iteration 62, residual = 3.9036807868964244e-13
CG iteration 63, residual = 2.5781504624508836e-13
CG iteration 64, residual = 1.7607220119455501e-13
CG iteration 65, residual = 1.1416505704637415e-13
CG iteration 66, residual = 7.696217324529728e-14
CG iteration 67, residual = 7.404857761605009e-14
CG iteration 68, residual = 5.026879238615375e-14
CG iteration 69, residual = 2.9882744779073735e-14
CG iteration 70, residual = 1.9789697148190823e-14
CG iteration 71, residual = 1.385029025125802e-14
CG iteration 72, residual = 1.0259418415452686e-14
CG iteration 73, residual = 6.347911070384764e-15
CG iteration 74, residual = 4.1774670519387576e-15
[7]:
# the vector potential is not supposed to look nice
Draw (gfu, mesh, "vector-potential", draw_surf=False, clipping=True)

Draw (curl(gfu), mesh, "B-field", draw_surf=False, clipping=True)
Draw (1/(mu0*mur)*curl(gfu)-mag, mesh, "H-field", draw_surf=False, clipping=True)
[7]:
BaseWebGuiScene
[ ]:

[ ]: