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.004809678530125143
CG iteration 2, residual = 0.0033223190132442406
CG iteration 3, residual = 0.0033115883177539988
CG iteration 4, residual = 0.0027467579779245866
CG iteration 5, residual = 0.0014658765270595335
CG iteration 6, residual = 0.0012167024999244303
CG iteration 7, residual = 0.0008096582309573253
CG iteration 8, residual = 0.0006570076838867988
CG iteration 9, residual = 0.00047508691834428516
CG iteration 10, residual = 0.00036224214074475817
CG iteration 11, residual = 0.0002544183424853305
CG iteration 12, residual = 0.00016194989094436856
CG iteration 13, residual = 0.00011358302702051901
CG iteration 14, residual = 8.950382818005822e-05
CG iteration 15, residual = 5.3941268122085824e-05
CG iteration 16, residual = 3.9160417835622944e-05
CG iteration 17, residual = 2.7082316503629277e-05
CG iteration 18, residual = 1.8249564195961574e-05
CG iteration 19, residual = 1.3584171505943444e-05
CG iteration 20, residual = 9.861474603202478e-06
CG iteration 21, residual = 1.169592719127149e-05
CG iteration 22, residual = 5.6974037777398235e-06
CG iteration 23, residual = 3.5939040582154726e-06
CG iteration 24, residual = 2.5589093403661247e-06
CG iteration 25, residual = 1.9779478100679428e-06
CG iteration 26, residual = 1.2770934682784486e-06
CG iteration 27, residual = 8.74083240367951e-07
CG iteration 28, residual = 5.927228746857655e-07
CG iteration 29, residual = 4.008901703446785e-07
CG iteration 30, residual = 2.856458196958117e-07
CG iteration 31, residual = 2.029268584655523e-07
CG iteration 32, residual = 1.3199636985563356e-07
CG iteration 33, residual = 9.850315462240981e-08
CG iteration 34, residual = 6.734294678385936e-08
CG iteration 35, residual = 4.755934018850992e-08
CG iteration 36, residual = 2.9131421380542793e-08
CG iteration 37, residual = 2.308337789665299e-08
CG iteration 38, residual = 2.5920863845176665e-08
CG iteration 39, residual = 1.3151244389605855e-08
CG iteration 40, residual = 8.396372748587321e-09
CG iteration 41, residual = 5.7119850958366605e-09
CG iteration 42, residual = 4.181049906672536e-09
CG iteration 43, residual = 2.7949493401940434e-09
CG iteration 44, residual = 1.8595331197302832e-09
CG iteration 45, residual = 1.2654495203786272e-09
CG iteration 46, residual = 8.583836388991428e-10
CG iteration 47, residual = 5.766587359042361e-10
CG iteration 48, residual = 4.217770157549865e-10
CG iteration 49, residual = 2.761187754272466e-10
CG iteration 50, residual = 2.1705289947937238e-10
CG iteration 51, residual = 1.2859266462568247e-10
CG iteration 52, residual = 8.627489350287253e-11
CG iteration 53, residual = 6.228470623639583e-11
CG iteration 54, residual = 8.218874928210366e-11
CG iteration 55, residual = 3.5541259840710894e-11
CG iteration 56, residual = 2.5361068255544396e-11
CG iteration 57, residual = 2.0884367952010116e-11
CG iteration 58, residual = 1.2878677095261773e-11
CG iteration 59, residual = 8.496391898424228e-12
CG iteration 60, residual = 5.923329711329595e-12
CG iteration 61, residual = 3.942608745636641e-12
CG iteration 62, residual = 2.6956164787177432e-12
CG iteration 63, residual = 1.7387078311195144e-12
CG iteration 64, residual = 1.2104137128461129e-12
CG iteration 65, residual = 7.623360981257901e-13
CG iteration 66, residual = 4.941726981817531e-13
CG iteration 67, residual = 3.207197816706313e-13
CG iteration 68, residual = 2.1313745606694316e-13
CG iteration 69, residual = 1.8750669017928274e-13
CG iteration 70, residual = 1.6926040158149583e-13
CG iteration 71, residual = 9.471233260551058e-14
CG iteration 72, residual = 6.131553726140268e-14
CG iteration 73, residual = 4.207091109216714e-14
CG iteration 74, residual = 2.763285497764286e-14
CG iteration 75, residual = 1.7696140262899673e-14
CG iteration 76, residual = 1.1888372322366841e-14
CG iteration 77, residual = 9.218117600383106e-15
CG iteration 78, residual = 6.735187487242989e-15
CG iteration 79, residual = 4.009784296830559e-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);
[ ]:

[ ]: