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 0x7ff3b40b3b80>
[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)
WARNING: maxsteps is deprecated, use maxiter instead!
CG iteration 1, residual = 0.004832539962849042
CG iteration 2, residual = 0.0028801984652913597
CG iteration 3, residual = 0.0022636840717767896
CG iteration 4, residual = 0.0018596783443640356
CG iteration 5, residual = 0.0013706814431410249
CG iteration 6, residual = 0.0009629315328718053
CG iteration 7, residual = 0.0006784248172513843
CG iteration 8, residual = 0.0004537347998869366
CG iteration 9, residual = 0.000337398452358466
CG iteration 10, residual = 0.0002663636696667053
CG iteration 11, residual = 0.00015785865716710004
CG iteration 12, residual = 0.00010596525810194424
CG iteration 13, residual = 7.440611757876552e-05
CG iteration 14, residual = 5.0217762320984526e-05
CG iteration 15, residual = 3.565457135377309e-05
CG iteration 16, residual = 2.165277559270105e-05
CG iteration 17, residual = 1.6241088584150868e-05
CG iteration 18, residual = 1.017386729900955e-05
CG iteration 19, residual = 7.247797858001008e-06
CG iteration 20, residual = 4.782087689205537e-06
CG iteration 21, residual = 3.2977746854222806e-06
CG iteration 22, residual = 2.3147388539389554e-06
CG iteration 23, residual = 1.440107146741313e-06
CG iteration 24, residual = 1.0154963026429293e-06
CG iteration 25, residual = 6.90963507241692e-07
CG iteration 26, residual = 5.167319166465788e-07
CG iteration 27, residual = 4.97848551985073e-07
CG iteration 28, residual = 2.87795058036198e-07
CG iteration 29, residual = 1.8381145530473123e-07
CG iteration 30, residual = 1.2576252926544437e-07
CG iteration 31, residual = 8.212391174499575e-08
CG iteration 32, residual = 5.655902227170783e-08
CG iteration 33, residual = 3.718600844900678e-08
CG iteration 34, residual = 2.395621703429854e-08
CG iteration 35, residual = 1.6498257419532506e-08
CG iteration 36, residual = 1.0759738081289346e-08
CG iteration 37, residual = 7.165063589751193e-09
CG iteration 38, residual = 4.497597038041852e-09
CG iteration 39, residual = 2.9156780432837635e-09
CG iteration 40, residual = 1.9460006798112054e-09
CG iteration 41, residual = 1.401595421896183e-09
CG iteration 42, residual = 8.454842765508706e-10
CG iteration 43, residual = 5.885400862834295e-10
CG iteration 44, residual = 3.706828147779532e-10
CG iteration 45, residual = 2.4975711303061003e-10
CG iteration 46, residual = 1.911087056721211e-10
CG iteration 47, residual = 1.6413761693551003e-10
CG iteration 48, residual = 1.2823125967498335e-10
CG iteration 49, residual = 8.853339277442909e-11
CG iteration 50, residual = 5.486331424641183e-11
CG iteration 51, residual = 3.648180730367436e-11
CG iteration 52, residual = 2.364643488226334e-11
CG iteration 53, residual = 1.632122725030354e-11
CG iteration 54, residual = 1.0677549250916197e-11
CG iteration 55, residual = 7.155945595665424e-12
CG iteration 56, residual = 4.7963208152052675e-12
CG iteration 57, residual = 3.098941064794934e-12
CG iteration 58, residual = 2.0500214967226044e-12
CG iteration 59, residual = 1.3892674950800082e-12
CG iteration 60, residual = 8.940758564898918e-13
CG iteration 61, residual = 5.967414825651534e-13
CG iteration 62, residual = 3.9036835371147246e-13
CG iteration 63, residual = 2.578159077652662e-13
CG iteration 64, residual = 1.7607382847651106e-13
CG iteration 65, residual = 1.141582487571004e-13
CG iteration 66, residual = 7.660441087469153e-14
CG iteration 67, residual = 7.074310005223806e-14
CG iteration 68, residual = 5.139976624136755e-14
CG iteration 69, residual = 2.998493870197188e-14
CG iteration 70, residual = 1.9968673237555384e-14
CG iteration 71, residual = 1.4229419520287527e-14
CG iteration 72, residual = 1.0310379779260522e-14
CG iteration 73, residual = 6.316240239324829e-15
CG iteration 74, residual = 4.17644220555311e-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
[ ]: