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

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.004841340224441176
CG iteration 2, residual = 0.005042736495448076
CG iteration 3, residual = 0.0031497474909857313
CG iteration 4, residual = 0.0028071729924236143
CG iteration 5, residual = 0.0018044283205905904
CG iteration 6, residual = 0.0012102064412063717
CG iteration 7, residual = 0.0009603395596834449
CG iteration 8, residual = 0.0008701989504798975
CG iteration 9, residual = 0.0006325610483462679
CG iteration 10, residual = 0.0004532678408037573
CG iteration 11, residual = 0.0003176167389469271
CG iteration 12, residual = 0.00023160559528073835
CG iteration 13, residual = 0.00020737408621508292
CG iteration 14, residual = 0.00013494579444833606
CG iteration 15, residual = 8.950767293589883e-05
CG iteration 16, residual = 6.496946495430704e-05
CG iteration 17, residual = 8.576981125209072e-05
CG iteration 18, residual = 4.168154528717489e-05
CG iteration 19, residual = 2.896228291847354e-05
CG iteration 20, residual = 1.9763307938682017e-05
CG iteration 21, residual = 1.2931661264989331e-05
CG iteration 22, residual = 9.125442581354188e-06
CG iteration 23, residual = 6.593115124196742e-06
CG iteration 24, residual = 4.346076508580304e-06
CG iteration 25, residual = 3.183636676658752e-06
CG iteration 26, residual = 2.0408071379092687e-06
CG iteration 27, residual = 1.520232603493176e-06
CG iteration 28, residual = 9.964537130012956e-07
CG iteration 29, residual = 7.344144752910433e-07
CG iteration 30, residual = 1.1374944630944586e-06
CG iteration 31, residual = 4.7480448951960487e-07
CG iteration 32, residual = 3.3521257176406486e-07
CG iteration 33, residual = 3.3073457527679086e-07
CG iteration 34, residual = 1.8744512024432851e-07
CG iteration 35, residual = 1.2818291630128365e-07
CG iteration 36, residual = 8.236198392602643e-08
CG iteration 37, residual = 5.807279330884868e-08
CG iteration 38, residual = 3.836531819591171e-08
CG iteration 39, residual = 2.7371690774227348e-08
CG iteration 40, residual = 1.8851428976036967e-08
CG iteration 41, residual = 1.2787965711708865e-08
CG iteration 42, residual = 1.001185854222392e-08
CG iteration 43, residual = 9.964355944073217e-09
CG iteration 44, residual = 9.625170452234161e-09
CG iteration 45, residual = 5.336601899292746e-09
CG iteration 46, residual = 3.4154142552997586e-09
CG iteration 47, residual = 2.4108684359982347e-09
CG iteration 48, residual = 2.1159800840119797e-09
CG iteration 49, residual = 1.46475303759538e-09
CG iteration 50, residual = 9.646512907666098e-10
CG iteration 51, residual = 6.139726098524082e-10
CG iteration 52, residual = 4.1375534096212836e-10
CG iteration 53, residual = 2.7945130686628115e-10
CG iteration 54, residual = 1.912122555202489e-10
CG iteration 55, residual = 1.2891783389552736e-10
CG iteration 56, residual = 9.068482984685906e-11
CG iteration 57, residual = 1.4627100827304846e-10
CG iteration 58, residual = 6.092246290060429e-11
CG iteration 59, residual = 4.3772696290532126e-11
CG iteration 60, residual = 4.0825813424684236e-11
CG iteration 61, residual = 2.197145063355673e-11
CG iteration 62, residual = 1.4625858393181818e-11
CG iteration 63, residual = 9.55541865556648e-12
CG iteration 64, residual = 6.71163167909754e-12
CG iteration 65, residual = 4.350798211952307e-12
CG iteration 66, residual = 3.1814986974978467e-12
CG iteration 67, residual = 2.077153452042909e-12
CG iteration 68, residual = 1.3359319855222575e-12
CG iteration 69, residual = 1.021548004292903e-12
CG iteration 70, residual = 1.386663551414497e-12
CG iteration 71, residual = 5.406510033161867e-13
CG iteration 72, residual = 3.779364885583652e-13
CG iteration 73, residual = 2.675106645450843e-13
CG iteration 74, residual = 2.433056510550929e-13
CG iteration 75, residual = 1.5295553897517435e-13
CG iteration 76, residual = 9.100518339196898e-14
CG iteration 77, residual = 6.657855369983222e-14
CG iteration 78, residual = 5.957191485017937e-14
CG iteration 79, residual = 3.904961936834696e-14
CG iteration 80, residual = 2.477046698464265e-14
CG iteration 81, residual = 1.6161609003529393e-14
CG iteration 82, residual = 1.2142628909030473e-14
CG iteration 83, residual = 1.0895191075885318e-14
CG iteration 84, residual = 1.146977390867298e-14
CG iteration 85, residual = 6.4706773148276774e-15
CG iteration 86, residual = 5.710346826413518e-15
CG iteration 87, residual = 3.3909730500159356e-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);
[ ]:

[ ]: