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.003322319013244254
CG iteration 3, residual = 0.0033115883177540235
CG iteration 4, residual = 0.0027467579779245884
CG iteration 5, residual = 0.0014658765270595376
CG iteration 6, residual = 0.0012167024999244344
CG iteration 7, residual = 0.0008096582309573288
CG iteration 8, residual = 0.0006570076838868025
CG iteration 9, residual = 0.00047508691834428814
CG iteration 10, residual = 0.0003622421407447603
CG iteration 11, residual = 0.00025441834248533243
CG iteration 12, residual = 0.00016194989094437008
CG iteration 13, residual = 0.00011358302702051981
CG iteration 14, residual = 8.950382818004574e-05
CG iteration 15, residual = 5.394126812159272e-05
CG iteration 16, residual = 3.9160417808859495e-05
CG iteration 17, residual = 2.7082315539137376e-05
CG iteration 18, residual = 1.8249519137072344e-05
CG iteration 19, residual = 1.3582350763722765e-05
CG iteration 20, residual = 9.81503814171092e-06
CG iteration 21, residual = 1.142050084451806e-05
CG iteration 22, residual = 5.739115607454143e-06
CG iteration 23, residual = 3.5941224679546945e-06
CG iteration 24, residual = 2.5589114343140366e-06
CG iteration 25, residual = 1.977947848959404e-06
CG iteration 26, residual = 1.2770934687908235e-06
CG iteration 27, residual = 8.740832403727209e-07
CG iteration 28, residual = 5.927228746858118e-07
CG iteration 29, residual = 4.008901703446805e-07
CG iteration 30, residual = 2.856458196958118e-07
CG iteration 31, residual = 2.0292685846552723e-07
CG iteration 32, residual = 1.3199636985439951e-07
CG iteration 33, residual = 9.850315456347249e-08
CG iteration 34, residual = 6.734294470306049e-08
CG iteration 35, residual = 4.755924945563088e-08
CG iteration 36, residual = 2.9126348762299984e-08
CG iteration 37, residual = 2.2743036727385148e-08
CG iteration 38, residual = 2.6174045056855043e-08
CG iteration 39, residual = 1.3182711691092847e-08
CG iteration 40, residual = 8.39659794515357e-09
CG iteration 41, residual = 5.7119866274733565e-09
CG iteration 42, residual = 4.181049923072965e-09
CG iteration 43, residual = 2.7949493403364955e-09
CG iteration 44, residual = 1.8595331195321567e-09
CG iteration 45, residual = 1.2654495194043614e-09
CG iteration 46, residual = 8.583836334875037e-10
CG iteration 47, residual = 5.766587080274636e-10
CG iteration 48, residual = 4.217769246938908e-10
CG iteration 49, residual = 2.761184527613205e-10
CG iteration 50, residual = 2.1705230864320657e-10
CG iteration 51, residual = 1.2859141950022875e-10
CG iteration 52, residual = 8.625566999427506e-11
CG iteration 53, residual = 6.165481666185553e-11
CG iteration 54, residual = 8.222037028801157e-11
CG iteration 55, residual = 3.536262832510192e-11
CG iteration 56, residual = 2.3267443233279652e-11
CG iteration 57, residual = 1.7504845872153298e-11
CG iteration 58, residual = 1.3525931283928205e-11
CG iteration 59, residual = 8.647843592581343e-12
CG iteration 60, residual = 5.9405929113891226e-12
CG iteration 61, residual = 3.9440937112944866e-12
CG iteration 62, residual = 2.6957567912553127e-12
CG iteration 63, residual = 1.7387160281934691e-12
CG iteration 64, residual = 1.2104143628848349e-12
CG iteration 65, residual = 7.62336134688047e-13
CG iteration 66, residual = 4.941727043614916e-13
CG iteration 67, residual = 3.2072021337709897e-13
CG iteration 68, residual = 2.131747363714954e-13
CG iteration 69, residual = 1.89393132645741e-13
CG iteration 70, residual = 1.6804845952407852e-13
CG iteration 71, residual = 9.4685059634303e-14
CG iteration 72, residual = 6.131559417231642e-14
CG iteration 73, residual = 4.207181959641362e-14
CG iteration 74, residual = 2.7636406658385283e-14
CG iteration 75, residual = 1.771772996965591e-14
CG iteration 76, residual = 1.2046264664322196e-14
CG iteration 77, residual = 9.485331263512134e-15
CG iteration 78, residual = 6.6804312402189264e-15
CG iteration 79, residual = 3.9966198059204554e-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);
[ ]:

[ ]: