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

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.004809678530125144
CG iteration 2, residual = 0.0033223190132442384
CG iteration 3, residual = 0.0033115883177540087
CG iteration 4, residual = 0.0027467579779245728
CG iteration 5, residual = 0.0014658765270595291
CG iteration 6, residual = 0.0012167024999244274
CG iteration 7, residual = 0.0008096582309573248
CG iteration 8, residual = 0.0006570076838867976
CG iteration 9, residual = 0.00047508691834428484
CG iteration 10, residual = 0.0003622421407447581
CG iteration 11, residual = 0.00025441834248533026
CG iteration 12, residual = 0.00016194989094436848
CG iteration 13, residual = 0.00011358302702051934
CG iteration 14, residual = 8.950382818007835e-05
CG iteration 15, residual = 5.3941268122853256e-05
CG iteration 16, residual = 3.91604178772213e-05
CG iteration 17, residual = 2.7082318002701106e-05
CG iteration 18, residual = 1.824963422880539e-05
CG iteration 19, residual = 1.358700042755724e-05
CG iteration 20, residual = 9.932648748222109e-06
CG iteration 21, residual = 1.199442026163102e-05
CG iteration 22, residual = 5.651897017353856e-06
CG iteration 23, residual = 3.5936702662126547e-06
CG iteration 24, residual = 2.558907099334125e-06
CG iteration 25, residual = 1.9779477684447895e-06
CG iteration 26, residual = 1.277093467730097e-06
CG iteration 27, residual = 8.74083240362855e-07
CG iteration 28, residual = 5.927228746857206e-07
CG iteration 29, residual = 4.0089017034467864e-07
CG iteration 30, residual = 2.8564581969581485e-07
CG iteration 31, residual = 2.029268584657071e-07
CG iteration 32, residual = 1.3199636986288536e-07
CG iteration 33, residual = 9.850315496846103e-08
CG iteration 34, residual = 6.734295900115594e-08
CG iteration 35, residual = 4.755987291491994e-08
CG iteration 36, residual = 2.9161176120466222e-08
CG iteration 37, residual = 2.491268825267502e-08
CG iteration 38, residual = 2.4420264131733007e-08
CG iteration 39, residual = 1.3064260116111161e-08
CG iteration 40, residual = 8.395757940121985e-09
CG iteration 41, residual = 5.711980914838391e-09
CG iteration 42, residual = 4.181049861860478e-09
CG iteration 43, residual = 2.7949493396945736e-09
CG iteration 44, residual = 1.8595331196791575e-09
CG iteration 45, residual = 1.2654495201458923e-09
CG iteration 46, residual = 8.583836376065832e-10
CG iteration 47, residual = 5.766587292461291e-10
CG iteration 48, residual = 4.2177699401519373e-10
CG iteration 49, residual = 2.7611869884666256e-10
CG iteration 50, residual = 2.1705277602536106e-10
CG iteration 51, residual = 1.2859304908174663e-10
CG iteration 52, residual = 8.632931214976317e-11
CG iteration 53, residual = 6.57217074795215e-11
CG iteration 54, residual = 7.908737145571526e-11
CG iteration 55, residual = 3.5157855659066073e-11
CG iteration 56, residual = 2.492062952426885e-11
CG iteration 57, residual = 2.0709586307058002e-11
CG iteration 58, residual = 1.295265013576492e-11
CG iteration 59, residual = 8.502606300961224e-12
CG iteration 60, residual = 5.9239388812022774e-12
CG iteration 61, residual = 3.942660031891841e-12
CG iteration 62, residual = 2.6956213117960028e-12
CG iteration 63, residual = 1.7387081133943132e-12
CG iteration 64, residual = 1.2104137352767245e-12
CG iteration 65, residual = 7.623361014336461e-13
CG iteration 66, residual = 4.94172880365677e-13
CG iteration 67, residual = 3.207368061618341e-13
CG iteration 68, residual = 2.146032485306518e-13
CG iteration 69, residual = 2.4031089167819587e-13
CG iteration 70, residual = 1.482888582510695e-13
CG iteration 71, residual = 9.43608697273049e-14
CG iteration 72, residual = 6.131439538358818e-14
CG iteration 73, residual = 4.2075998069618303e-14
CG iteration 74, residual = 2.7652997076156222e-14
CG iteration 75, residual = 1.7818483727030582e-14
CG iteration 76, residual = 1.2700904639474365e-14
CG iteration 77, residual = 9.960222799258856e-15
CG iteration 78, residual = 6.536072955472234e-15
CG iteration 79, residual = 3.976860202550116e-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);
[ ]:

[ ]: