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

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.004846766222292574
CG iteration 2, residual = 0.004063969731501557
CG iteration 3, residual = 0.0032905813897480497
CG iteration 4, residual = 0.0022607529739440304
CG iteration 5, residual = 0.001975789805600883
CG iteration 6, residual = 0.001150789529109967
CG iteration 7, residual = 0.0008639313194094767
CG iteration 8, residual = 0.0006327449629420533
CG iteration 9, residual = 0.0004949733354845948
CG iteration 10, residual = 0.000407605412414322
CG iteration 11, residual = 0.0003164866826507561
CG iteration 12, residual = 0.00018393955806005982
CG iteration 13, residual = 0.00013713589121295788
CG iteration 14, residual = 8.753130614956616e-05
CG iteration 15, residual = 6.404693882372028e-05
CG iteration 16, residual = 4.235156751948365e-05
CG iteration 17, residual = 2.6160978407818126e-05
CG iteration 18, residual = 1.9549690624864868e-05
CG iteration 19, residual = 1.3706932985390065e-05
CG iteration 20, residual = 1.786835019062794e-05
CG iteration 21, residual = 8.024090476884843e-06
CG iteration 22, residual = 5.6873251602057515e-06
CG iteration 23, residual = 3.828295843896282e-06
CG iteration 24, residual = 2.512607975609427e-06
CG iteration 25, residual = 1.9133034154962675e-06
CG iteration 26, residual = 1.1882807370919078e-06
CG iteration 27, residual = 8.414430277905431e-07
CG iteration 28, residual = 5.733998554464014e-07
CG iteration 29, residual = 3.8563170778110635e-07
CG iteration 30, residual = 2.8128909681632414e-07
CG iteration 31, residual = 1.790462546330986e-07
CG iteration 32, residual = 1.3356424957586778e-07
CG iteration 33, residual = 1.1625360789667777e-07
CG iteration 34, residual = 7.842498482703725e-08
CG iteration 35, residual = 4.951106056608611e-08
CG iteration 36, residual = 4.540189835676575e-08
CG iteration 37, residual = 4.153057854609812e-08
CG iteration 38, residual = 2.1739658754304404e-08
CG iteration 39, residual = 1.4458259388836865e-08
CG iteration 40, residual = 9.468150258472693e-09
CG iteration 41, residual = 6.480348639216918e-09
CG iteration 42, residual = 6.111683861087161e-09
CG iteration 43, residual = 3.887823541942877e-09
CG iteration 44, residual = 2.6327130294972733e-09
CG iteration 45, residual = 1.7233169259306949e-09
CG iteration 46, residual = 1.1298985968368898e-09
CG iteration 47, residual = 7.954174244088379e-10
CG iteration 48, residual = 5.005044072262324e-10
CG iteration 49, residual = 3.5051614913088827e-10
CG iteration 50, residual = 2.2705014624846867e-10
CG iteration 51, residual = 1.5452337634830617e-10
CG iteration 52, residual = 1.30619161222193e-10
CG iteration 53, residual = 1.2340492843512664e-10
CG iteration 54, residual = 6.063096969841642e-11
CG iteration 55, residual = 4.031783151871738e-11
CG iteration 56, residual = 2.764890893920315e-11
CG iteration 57, residual = 1.7867597501233748e-11
CG iteration 58, residual = 1.208476180702244e-11
CG iteration 59, residual = 7.85902786623752e-12
CG iteration 60, residual = 5.62299733553768e-12
CG iteration 61, residual = 5.562910921677541e-12
CG iteration 62, residual = 3.0799609584333397e-12
CG iteration 63, residual = 1.977614575199879e-12
CG iteration 64, residual = 1.2730127093023197e-12
CG iteration 65, residual = 8.764731567000056e-13
CG iteration 66, residual = 6.020965028732679e-13
CG iteration 67, residual = 3.884800076640244e-13
CG iteration 68, residual = 5.503968512956826e-13
CG iteration 69, residual = 2.5340979022382327e-13
CG iteration 70, residual = 2.4083875409584796e-13
CG iteration 71, residual = 1.467726228185055e-13
CG iteration 72, residual = 8.93536089653212e-14
CG iteration 73, residual = 5.863484394383812e-14
CG iteration 74, residual = 4.0246896860594264e-14
CG iteration 75, residual = 2.6754014194805164e-14
CG iteration 76, residual = 1.7722751672555513e-14
CG iteration 77, residual = 1.1619559466601413e-14
CG iteration 78, residual = 7.602279748808699e-15
CG iteration 79, residual = 4.827908386850748e-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);
[ ]:

[ ]: