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.005042736495448089
CG iteration 3, residual = 0.003149747490985731
CG iteration 4, residual = 0.0028071729924236148
CG iteration 5, residual = 0.0018044283205905988
CG iteration 6, residual = 0.001210206441206379
CG iteration 7, residual = 0.0009603395596834496
CG iteration 8, residual = 0.0008701989504799085
CG iteration 9, residual = 0.0006325610483462708
CG iteration 10, residual = 0.00045326784080376064
CG iteration 11, residual = 0.000317616738946878
CG iteration 12, residual = 0.0002316055952750398
CG iteration 13, residual = 0.00020737408604002216
CG iteration 14, residual = 0.00013494578871856523
CG iteration 15, residual = 8.950719878195868e-05
CG iteration 16, residual = 6.491463624730127e-05
CG iteration 17, residual = 8.455640229469506e-05
CG iteration 18, residual = 4.1837115232582976e-05
CG iteration 19, residual = 2.8962982496000066e-05
CG iteration 20, residual = 1.9763311190758055e-05
CG iteration 21, residual = 1.2931661275269884e-05
CG iteration 22, residual = 9.12544258138533e-06
CG iteration 23, residual = 6.593115124201802e-06
CG iteration 24, residual = 4.346076508650894e-06
CG iteration 25, residual = 3.1836366773883166e-06
CG iteration 26, residual = 2.040807149133876e-06
CG iteration 27, residual = 1.5202325824258516e-06
CG iteration 28, residual = 9.964327864778882e-07
CG iteration 29, residual = 7.31773704332905e-07
CG iteration 30, residual = 1.1229834892298272e-06
CG iteration 31, residual = 4.76822093880901e-07
CG iteration 32, residual = 3.381816658216281e-07
CG iteration 33, residual = 3.294601528007855e-07
CG iteration 34, residual = 1.8718859240352186e-07
CG iteration 35, residual = 1.2817262824363007e-07
CG iteration 36, residual = 8.236178042531437e-08
CG iteration 37, residual = 5.8072785847919765e-08
CG iteration 38, residual = 3.836529399648158e-08
CG iteration 39, residual = 2.7371357310039538e-08
CG iteration 40, residual = 1.88477705585932e-08
CG iteration 41, residual = 1.2743085645401121e-08
CG iteration 42, residual = 9.637402036425255e-09
CG iteration 43, residual = 1.0245488785866412e-08
CG iteration 44, residual = 8.778402093517629e-09
CG iteration 45, residual = 5.497957047580652e-09
CG iteration 46, residual = 3.4205771942433558e-09
CG iteration 47, residual = 2.3527983759113126e-09
CG iteration 48, residual = 1.6197264628808092e-09
CG iteration 49, residual = 1.4323845014391073e-09
CG iteration 50, residual = 1.046756192424169e-09
CG iteration 51, residual = 6.167462015672361e-10
CG iteration 52, residual = 4.138240371562979e-10
CG iteration 53, residual = 2.7945320784357647e-10
CG iteration 54, residual = 1.912123165673685e-10
CG iteration 55, residual = 1.2891786839839103e-10
CG iteration 56, residual = 9.065699239328939e-11
CG iteration 57, residual = 1.4622080485409385e-10
CG iteration 58, residual = 6.097457098230558e-11
CG iteration 59, residual = 4.4226937971732164e-11
CG iteration 60, residual = 4.055462041402022e-11
CG iteration 61, residual = 2.195680305258772e-11
CG iteration 62, residual = 1.4625558884633547e-11
CG iteration 63, residual = 9.555413758059595e-12
CG iteration 64, residual = 6.711631565011444e-12
CG iteration 65, residual = 4.3507982100126935e-12
CG iteration 66, residual = 3.1814986997713704e-12
CG iteration 67, residual = 2.077153781811838e-12
CG iteration 68, residual = 1.3360009035811384e-12
CG iteration 69, residual = 1.0326580557993054e-12
CG iteration 70, residual = 1.3625604384147625e-12
CG iteration 71, residual = 5.404745133602227e-13
CG iteration 72, residual = 3.774814819448205e-13
CG iteration 73, residual = 2.615925866412988e-13
CG iteration 74, residual = 1.964710142430213e-13
CG iteration 75, residual = 1.707833968707878e-13
CG iteration 76, residual = 9.135190544420854e-14
CG iteration 77, residual = 6.51777395567606e-14
CG iteration 78, residual = 4.6419357357291975e-14
CG iteration 79, residual = 4.417025224475848e-14
CG iteration 80, residual = 2.5075804795065262e-14
CG iteration 81, residual = 1.602713294624545e-14
CG iteration 82, residual = 1.1671462963662311e-14
CG iteration 83, residual = 1.467867956353133e-14
CG iteration 84, residual = 9.130109103147672e-15
CG iteration 85, residual = 6.3936659713649175e-15
CG iteration 86, residual = 5.512785609356204e-15
CG iteration 87, residual = 3.4396620842010227e-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);
[ ]:

[ ]: