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.occ 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
box = Box( (-3,-3,-3), (3,3,3))
box.faces.name = "outer"

magnet = Cylinder((-1,0,0),X, r=0.3, h=2)
magnet.mat("magnet")
magnet.faces.col = (1,0,0)

air = box-magnet
air.mat("air")
shape = Glue([air,magnet])
geo = OCCGeometry(shape)

Draw (shape, clipping={ "z" : -1, "function":True})

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

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, printrates=True)
CG iteration 1, residual = 0.0048467662223670846
CG iteration 2, residual = 0.004063969731223297
CG iteration 3, residual = 0.003290581389985664
CG iteration 4, residual = 0.0022607529740109175
CG iteration 5, residual = 0.001975789805749372
CG iteration 6, residual = 0.0011507895290807876
CG iteration 7, residual = 0.0008639313194826863
CG iteration 8, residual = 0.0006327449630705042
CG iteration 9, residual = 0.0004949733353699826
CG iteration 10, residual = 0.0004076054117917417
CG iteration 11, residual = 0.0003164866829728323
CG iteration 12, residual = 0.00018393955810877544
CG iteration 13, residual = 0.00013713589123244755
CG iteration 14, residual = 8.753130615553811e-05
CG iteration 15, residual = 6.404693878509319e-05
CG iteration 16, residual = 4.235156479732741e-05
CG iteration 17, residual = 2.616074513020146e-05
CG iteration 18, residual = 1.953158540405145e-05
CG iteration 19, residual = 1.2890330394317926e-05
CG iteration 20, residual = 1.679374195783725e-05
CG iteration 21, residual = 8.318118514290396e-06
CG iteration 22, residual = 5.6905108796282555e-06
CG iteration 23, residual = 3.8283251266826155e-06
CG iteration 24, residual = 2.5126081409070526e-06
CG iteration 25, residual = 1.9133034181940437e-06
CG iteration 26, residual = 1.1882807365298684e-06
CG iteration 27, residual = 8.414430158583573e-07
CG iteration 28, residual = 5.733996816170292e-07
CG iteration 29, residual = 3.8562933851698383e-07
CG iteration 30, residual = 2.812624998735838e-07
CG iteration 31, residual = 1.787000267279629e-07
CG iteration 32, residual = 1.2810773917514924e-07
CG iteration 33, residual = 8.54219770362267e-08
CG iteration 34, residual = 8.468097573352891e-08
CG iteration 35, residual = 5.194744854138652e-08
CG iteration 36, residual = 4.4601347388207606e-08
CG iteration 37, residual = 4.218667864584577e-08
CG iteration 38, residual = 2.1776026931855853e-08
CG iteration 39, residual = 1.4486485099081498e-08
CG iteration 40, residual = 9.92480156091136e-09
CG iteration 41, residual = 9.900343942140837e-09
CG iteration 42, residual = 5.825143409888559e-09
CG iteration 43, residual = 3.6107486588435233e-09
CG iteration 44, residual = 2.6219691508038524e-09
CG iteration 45, residual = 1.7228312591187435e-09
CG iteration 46, residual = 1.129888780428514e-09
CG iteration 47, residual = 7.954170936295282e-10
CG iteration 48, residual = 5.0050439963471e-10
CG iteration 49, residual = 3.5051609124720297e-10
CG iteration 50, residual = 2.27045772532501e-10
CG iteration 51, residual = 1.541421479783508e-10
CG iteration 52, residual = 1.0882371646986881e-10
CG iteration 53, residual = 1.5022625857223836e-10
CG iteration 54, residual = 6.11515832130514e-11
CG iteration 55, residual = 4.032010097634139e-11
CG iteration 56, residual = 2.7648919269150712e-11
CG iteration 57, residual = 1.7867535486913367e-11
CG iteration 58, residual = 1.2083451687514519e-11
CG iteration 59, residual = 7.834444459578419e-12
CG iteration 60, residual = 4.988192875742276e-12
CG iteration 61, residual = 3.5922341055591435e-12
CG iteration 62, residual = 3.306240592148215e-12
CG iteration 63, residual = 2.1349757691874792e-12
CG iteration 64, residual = 1.2765156449749133e-12
CG iteration 65, residual = 8.767869766936132e-13
CG iteration 66, residual = 6.056796878119146e-13
CG iteration 67, residual = 4.412175555488807e-13
CG iteration 68, residual = 4.0604198614202363e-13
CG iteration 69, residual = 2.29877820617353e-13
CG iteration 70, residual = 2.8859672995975756e-13
CG iteration 71, residual = 1.465167791308156e-13
CG iteration 72, residual = 8.919495451644322e-14
CG iteration 73, residual = 5.863182603920526e-14
CG iteration 74, residual = 4.024682020830959e-14
CG iteration 75, residual = 2.6754014524711405e-14
CG iteration 76, residual = 1.772275748433823e-14
CG iteration 77, residual = 1.1619582418086307e-14
CG iteration 78, residual = 7.602387030552521e-15
CG iteration 79, residual = 4.828459378341719e-15
[7]:
Draw (curl(gfu), mesh, "B-field", draw_surf=False, \
      clipping = { "z" : -1, "function":True}, \
      vectors = { "grid_size":50}, min=0, max=2e-5);
[ ]:

[ ]: