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.004846766222367084
CG iteration 2, residual = 0.004063969731223288
CG iteration 3, residual = 0.0032905813899856786
CG iteration 4, residual = 0.0022607529740109223
CG iteration 5, residual = 0.001975789805749375
CG iteration 6, residual = 0.0011507895290807897
CG iteration 7, residual = 0.0008639313194826884
CG iteration 8, residual = 0.0006327449630705043
CG iteration 9, residual = 0.0004949733353699721
CG iteration 10, residual = 0.00040760541179170665
CG iteration 11, residual = 0.0003164866829728489
CG iteration 12, residual = 0.0001839395581087769
CG iteration 13, residual = 0.00013713589123245094
CG iteration 14, residual = 8.753130615568689e-05
CG iteration 15, residual = 6.404693879430606e-05
CG iteration 16, residual = 4.235156521137463e-05
CG iteration 17, residual = 2.616078059411122e-05
CG iteration 18, residual = 1.9534339579306236e-05
CG iteration 19, residual = 1.302054061616373e-05
CG iteration 20, residual = 1.7574062434873906e-05
CG iteration 21, residual = 8.201537538654727e-06
CG iteration 22, residual = 5.689206686323466e-06
CG iteration 23, residual = 3.828313132802762e-06
CG iteration 24, residual = 2.51260807337417e-06
CG iteration 25, residual = 1.913303417315332e-06
CG iteration 26, residual = 1.1882807365125641e-06
CG iteration 27, residual = 8.414430156842323e-07
CG iteration 28, residual = 5.733996791261955e-07
CG iteration 29, residual = 3.8562930457179304e-07
CG iteration 30, residual = 2.8126211896690714e-07
CG iteration 31, residual = 1.786950586169127e-07
CG iteration 32, residual = 1.2802477286708976e-07
CG iteration 33, residual = 8.441110130438809e-08
CG iteration 34, residual = 8.192589417229008e-08
CG iteration 35, residual = 5.2819970129817834e-08
CG iteration 36, residual = 4.679111969949478e-08
CG iteration 37, residual = 4.075266619229848e-08
CG iteration 38, residual = 2.174941298438666e-08
CG iteration 39, residual = 1.4464484455133477e-08
CG iteration 40, residual = 9.565616154464139e-09
CG iteration 41, residual = 7.96541941574555e-09
CG iteration 42, residual = 6.428900122437788e-09
CG iteration 43, residual = 3.634036560427776e-09
CG iteration 44, residual = 2.622700039930094e-09
CG iteration 45, residual = 1.7228638162582812e-09
CG iteration 46, residual = 1.1298894381834916e-09
CG iteration 47, residual = 7.954171158150016e-10
CG iteration 48, residual = 5.005044000433716e-10
CG iteration 49, residual = 3.5051608665289404e-10
CG iteration 50, residual = 2.2704542381127316e-10
CG iteration 51, residual = 1.5411168969079966e-10
CG iteration 52, residual = 1.0671723069314527e-10
CG iteration 53, residual = 1.5245661511011925e-10
CG iteration 54, residual = 6.138890038227151e-11
CG iteration 55, residual = 4.03211499281148e-11
CG iteration 56, residual = 2.7648925523195808e-11
CG iteration 57, residual = 1.7867535824149455e-11
CG iteration 58, residual = 1.2083458119762513e-11
CG iteration 59, residual = 7.834565590034458e-12
CG iteration 60, residual = 4.991746719957276e-12
CG iteration 61, residual = 3.665423766641319e-12
CG iteration 62, residual = 3.5008402611260288e-12
CG iteration 63, residual = 2.0806940629074273e-12
CG iteration 64, residual = 1.275137984970803e-12
CG iteration 65, residual = 8.765428977636886e-13
CG iteration 66, residual = 6.023888121153498e-13
CG iteration 67, residual = 3.8387241687017273e-13
CG iteration 68, residual = 3.9687110994626357e-13
CG iteration 69, residual = 2.883039451591226e-13
CG iteration 70, residual = 2.8303321913127094e-13
CG iteration 71, residual = 1.4002978606687492e-13
CG iteration 72, residual = 8.91809285357838e-14
CG iteration 73, residual = 5.863207657284734e-14
CG iteration 74, residual = 4.024682806808561e-14
CG iteration 75, residual = 2.675401037430638e-14
CG iteration 76, residual = 1.7722745842540177e-14
CG iteration 77, residual = 1.1619538769024652e-14
CG iteration 78, residual = 7.602190618246484e-15
CG iteration 79, residual = 4.827483862054772e-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);
[ ]:

[ ]: