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.0048096785301251456
CG iteration 2, residual = 0.003322319013244248
CG iteration 3, residual = 0.0033115883177540213
CG iteration 4, residual = 0.0027467579779245793
CG iteration 5, residual = 0.001465876527059535
CG iteration 6, residual = 0.0012167024999244322
CG iteration 7, residual = 0.0008096582309573275
CG iteration 8, residual = 0.000657007683886801
CG iteration 9, residual = 0.00047508691834428717
CG iteration 10, residual = 0.0003622421407447591
CG iteration 11, residual = 0.00025441834248533145
CG iteration 12, residual = 0.00016194989094436922
CG iteration 13, residual = 0.00011358302702051953
CG iteration 14, residual = 8.950382818005711e-05
CG iteration 15, residual = 5.394126812203761e-05
CG iteration 16, residual = 3.9160417832988265e-05
CG iteration 17, residual = 2.7082316408668464e-05
CG iteration 18, residual = 1.824955975961164e-05
CG iteration 19, residual = 1.358399226334221e-05
CG iteration 20, residual = 9.856925200297016e-06
CG iteration 21, residual = 1.1672185035047181e-05
CG iteration 22, residual = 5.700975271296948e-06
CG iteration 23, residual = 3.593922603522579e-06
CG iteration 24, residual = 2.5589095181508826e-06
CG iteration 25, residual = 1.9779478133699813e-06
CG iteration 26, residual = 1.2770934683219537e-06
CG iteration 27, residual = 8.740832403683578e-07
CG iteration 28, residual = 5.927228746857693e-07
CG iteration 29, residual = 4.0089017034467906e-07
CG iteration 30, residual = 2.8564581969581347e-07
CG iteration 31, residual = 2.0292685846560852e-07
CG iteration 32, residual = 1.3199636985823652e-07
CG iteration 33, residual = 9.850315474659601e-08
CG iteration 34, residual = 6.734295116823182e-08
CG iteration 35, residual = 4.7559531367207864e-08
CG iteration 36, residual = 2.914210501446377e-08
CG iteration 37, residual = 2.377145434664581e-08
CG iteration 38, residual = 2.5348151084433113e-08
CG iteration 39, residual = 1.3107284110683976e-08
CG iteration 40, residual = 8.396060633655318e-09
CG iteration 41, residual = 5.711982973205267e-09
CG iteration 42, residual = 4.181049883920803e-09
CG iteration 43, residual = 2.7949493399371733e-09
CG iteration 44, residual = 1.859533119686744e-09
CG iteration 45, residual = 1.2654495201745907e-09
CG iteration 46, residual = 8.583836377659406e-10
CG iteration 47, residual = 5.766587300670479e-10
CG iteration 48, residual = 4.217769966963711e-10
CG iteration 49, residual = 2.761187083262066e-10
CG iteration 50, residual = 2.1705279260536034e-10
CG iteration 51, residual = 1.2859305412127208e-10
CG iteration 52, residual = 8.632714172544615e-11
CG iteration 53, residual = 6.557977075345888e-11
CG iteration 54, residual = 7.925620931593986e-11
CG iteration 55, residual = 3.517320747544392e-11
CG iteration 56, residual = 2.4976530586137566e-11
CG iteration 57, residual = 2.0738870360107556e-11
CG iteration 58, residual = 1.2941933267669771e-11
CG iteration 59, residual = 8.501670625023827e-12
CG iteration 60, residual = 5.923846724434238e-12
CG iteration 61, residual = 3.942652268102095e-12
CG iteration 62, residual = 2.695620580097348e-12
CG iteration 63, residual = 1.738708070654579e-12
CG iteration 64, residual = 1.210413731865133e-12
CG iteration 65, residual = 7.623361011777225e-13
CG iteration 66, residual = 4.941728841979672e-13
CG iteration 67, residual = 3.2073724182702684e-13
CG iteration 68, residual = 2.14641105500546e-13
CG iteration 69, residual = 2.4126632130548694e-13
CG iteration 70, residual = 1.4806449034575843e-13
CG iteration 71, residual = 9.435758660599259e-14
CG iteration 72, residual = 6.131270799515811e-14
CG iteration 73, residual = 4.2070027794156755e-14
CG iteration 74, residual = 2.76294431959026e-14
CG iteration 75, residual = 1.7675110987001787e-14
CG iteration 76, residual = 1.1725768538363476e-14
CG iteration 77, residual = 8.809172930686377e-15
CG iteration 78, residual = 6.770477627427718e-15
CG iteration 79, residual = 4.038959802087778e-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);
[ ]:

[ ]: