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 0x7fcf5d1b2c70>
[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.0048325399628490456
CG iteration 2, residual = 0.0028801984652913554
CG iteration 3, residual = 0.0022636840717767865
CG iteration 4, residual = 0.0018596783443640308
CG iteration 5, residual = 0.0013706814431410227
CG iteration 6, residual = 0.0009629315328718032
CG iteration 7, residual = 0.0006784248172513837
CG iteration 8, residual = 0.00045373479988693643
CG iteration 9, residual = 0.0003373984523584661
CG iteration 10, residual = 0.00026636366966670543
CG iteration 11, residual = 0.00015785865716710015
CG iteration 12, residual = 0.00010596525810194436
CG iteration 13, residual = 7.44061175787656e-05
CG iteration 14, residual = 5.021776232098458e-05
CG iteration 15, residual = 3.565457135377311e-05
CG iteration 16, residual = 2.1652775592701058e-05
CG iteration 17, residual = 1.6241088584150993e-05
CG iteration 18, residual = 1.017386729901183e-05
CG iteration 19, residual = 7.247797858043955e-06
CG iteration 20, residual = 4.782087690169919e-06
CG iteration 21, residual = 3.29777470438918e-06
CG iteration 22, residual = 2.3147392389559875e-06
CG iteration 23, residual = 1.4401159990987946e-06
CG iteration 24, residual = 1.0157366363760883e-06
CG iteration 25, residual = 6.956810136230163e-07
CG iteration 26, residual = 5.74237137795663e-07
CG iteration 27, residual = 4.817251754068111e-07
CG iteration 28, residual = 2.827641003814129e-07
CG iteration 29, residual = 1.8371509385496837e-07
CG iteration 30, residual = 1.2576047785542574e-07
CG iteration 31, residual = 8.212387420856848e-08
CG iteration 32, residual = 5.655902150085426e-08
CG iteration 33, residual = 3.718600843572007e-08
CG iteration 34, residual = 2.3956217034117113e-08
CG iteration 35, residual = 1.649825741954535e-08
CG iteration 36, residual = 1.0759738081366624e-08
CG iteration 37, residual = 7.1650635903059945e-09
CG iteration 38, residual = 4.497597041560091e-09
CG iteration 39, residual = 2.915678082057253e-09
CG iteration 40, residual = 1.946001042119997e-09
CG iteration 41, residual = 1.4015973107772708e-09
CG iteration 42, residual = 8.454942009987504e-10
CG iteration 43, residual = 5.886246227171942e-10
CG iteration 44, residual = 3.712905553482017e-10
CG iteration 45, residual = 2.561099409236985e-10
CG iteration 46, residual = 2.2578223583369505e-10
CG iteration 47, residual = 1.7876135162532347e-10
CG iteration 48, residual = 1.3106897620734052e-10
CG iteration 49, residual = 8.59623092714088e-11
CG iteration 50, residual = 5.429090483744898e-11
CG iteration 51, residual = 3.644176051524264e-11
CG iteration 52, residual = 2.3643484535969752e-11
CG iteration 53, residual = 1.632105498140577e-11
CG iteration 54, residual = 1.0677532078512073e-11
CG iteration 55, residual = 7.155944473052101e-12
CG iteration 56, residual = 4.796320697298884e-12
CG iteration 57, residual = 3.0989410582893747e-12
CG iteration 58, residual = 2.0500215005961733e-12
CG iteration 59, residual = 1.3892675071206807e-12
CG iteration 60, residual = 8.940758950937319e-13
CG iteration 61, residual = 5.967415897835504e-13
CG iteration 62, residual = 3.9036870680033935e-13
CG iteration 63, residual = 2.5781723950111103e-13
CG iteration 64, residual = 1.7608230134619278e-13
CG iteration 65, residual = 1.1433010679554052e-13
CG iteration 66, residual = 8.131912720031088e-14
CG iteration 67, residual = 8.58922967909135e-14
CG iteration 68, residual = 4.69623854475367e-14
CG iteration 69, residual = 2.9970911581909587e-14
CG iteration 70, residual = 2.017154159328763e-14
CG iteration 71, residual = 1.456924608195764e-14
CG iteration 72, residual = 1.0306665728430678e-14
CG iteration 73, residual = 6.2944285808616154e-15
CG iteration 74, residual = 4.175754806222453e-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
[ ]:

[ ]: