This page was generated from jupyter-files/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.

In [1]:
from ngsolve import *
from netgen.csg import *
import netgen.gui
%gui tk

Geometric model and meshing of a bar magnet:

In [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").transp())
geo.Add (magnet.mat("magnet").maxh(1))
geo.Draw()
mesh = Mesh(geo.GenerateMesh(maxh=2, curvaturesafety=1))
mesh.Curve(3)
In [3]:
mesh.GetMaterials(), mesh.GetBoundaries()
Out[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’)
In [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 += SymbolicBFI(1/(mu0*mur) * curl(u) * curl(v) + 1e-8/(mu0*mur)*u*v)
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 = 21439

Assemble system and setup preconditioner using task-parallelization:

In [5]:
with TaskManager():
    a.Assemble()
    f.Assemble()

Finally, declare GridFunction and solve by preconditioned CG iteration:

In [6]:
gfu = GridFunction(fes)
with TaskManager():
    solvers.CG(sol=gfu.vec, rhs=f.vec, mat=a.mat, pre=c.mat)
it =  0  err =  0.004834556218517527
it =  1  err =  0.0028938795257410117
it =  2  err =  0.0018764053599106428
it =  3  err =  0.0012167024420446035
it =  4  err =  0.0008154078108554166
it =  5  err =  0.0005512980608231248
it =  6  err =  0.0003783687160655943
it =  7  err =  0.00024850704793636
it =  8  err =  0.00014751591939837823
it =  9  err =  0.00010046059754497959
it =  10  err =  6.0599386314131916e-05
it =  11  err =  3.7306746707688986e-05
it =  12  err =  2.2865488298664312e-05
it =  13  err =  1.3402569890317977e-05
it =  14  err =  8.744091600065984e-06
it =  15  err =  5.5933636969128035e-06
it =  16  err =  3.192193939595891e-06
it =  17  err =  1.911582068812383e-06
it =  18  err =  1.1537923931961965e-06
it =  19  err =  6.66288841368631e-07
it =  20  err =  4.1644321801017405e-07
it =  21  err =  2.7183534856130845e-07
it =  22  err =  1.6639462226981296e-07
it =  23  err =  1.023184651335171e-07
it =  24  err =  5.894888743781949e-08
it =  25  err =  3.5542134277234365e-08
it =  26  err =  2.1714186678656858e-08
it =  27  err =  1.3111960150866855e-08
it =  28  err =  7.622354075626605e-09
it =  29  err =  4.451915516950855e-09
it =  30  err =  2.5667653447970765e-09
it =  31  err =  1.5254455666689282e-09
it =  32  err =  9.506467120575732e-10
it =  33  err =  5.620910321157875e-10
it =  34  err =  3.4593450793876496e-10
it =  35  err =  1.9766816375141168e-10
it =  36  err =  1.174537677413651e-10
it =  37  err =  7.079331072431792e-11
it =  38  err =  4.15449701360247e-11
it =  39  err =  2.4756971398668627e-11
it =  40  err =  1.5341064928142025e-11
it =  41  err =  9.156455648648875e-12
it =  42  err =  5.280191027621488e-12
it =  43  err =  3.2879905796687907e-12
it =  44  err =  1.937322886980886e-12
it =  45  err =  1.1575787499694853e-12
it =  46  err =  6.775420772822517e-13
it =  47  err =  4.017300590188354e-13
it =  48  err =  2.740211128678103e-13
it =  49  err =  1.864275852353925e-13
it =  50  err =  1.063271698211812e-13
it =  51  err =  6.335621398663346e-14
it =  52  err =  3.7404270848382554e-14
it =  53  err =  2.3138588501848596e-14
it =  54  err =  1.4444778054166232e-14
it =  55  err =  9.297641777021612e-15
it =  56  err =  6.14275145884241e-15
In [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)