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\):
Introducing a vector-potential \(A\) such that \(B = \Curl A\), and putting equations together we get
In weak form: Find \(A \in H(\Curl)\) such that
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 netgen.csg import *
import netgen.gui
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
geo = CSGeometry()
geo.Add (air.mat("air"), transparent=True)
geo.Add (magnet.mat("magnet").maxh(1), col=(0.3,0.3,0.1))
geo.Draw()
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 = CoefficientFunction( [1000 if mat== "magnet" else 1
for mat in mesh.GetMaterials()])
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 = CoefficientFunction((1,0,0)) * \
CoefficientFunction( [1 if mat == "magnet" else 0 for mat in mesh.GetMaterials()])
f += SymbolicLFI(mag*curl(v), definedon=mesh.Materials("magnet"))
ndof = 25442
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)
iteration 0 error = 0.004770751268765441
iteration 1 error = 0.002481327439932446
iteration 2 error = 0.001891809576612191
iteration 3 error = 0.0012587090378891413
iteration 4 error = 0.001136453623318983
iteration 5 error = 0.0008195457478849471
iteration 6 error = 0.0006060825336657849
iteration 7 error = 0.0004534128861287396
iteration 8 error = 0.0003578513184707288
iteration 9 error = 0.00023195184860900925
iteration 10 error = 0.00015491751841894275
iteration 11 error = 9.808156355275174e-05
iteration 12 error = 7.798886269681473e-05
iteration 13 error = 4.431269578808666e-05
iteration 14 error = 3.282881496334983e-05
iteration 15 error = 2.2920087554376978e-05
iteration 16 error = 1.5012379279908444e-05
iteration 17 error = 1.0290586276555383e-05
iteration 18 error = 6.797403084573732e-06
iteration 19 error = 5.324080434202823e-06
iteration 20 error = 3.093510781029754e-06
iteration 21 error = 2.27325012997849e-06
iteration 22 error = 1.5139195317352793e-06
iteration 23 error = 9.976809810274086e-07
iteration 24 error = 6.972869052375705e-07
iteration 25 error = 4.710001624204075e-07
iteration 26 error = 3.089882446017927e-07
iteration 27 error = 2.120121306552677e-07
iteration 28 error = 1.4360554950468746e-07
iteration 29 error = 9.438532193607466e-08
iteration 30 error = 6.13654863017183e-08
iteration 31 error = 4.3555240464504774e-08
iteration 32 error = 2.9553556227435373e-08
iteration 33 error = 1.922000695267886e-08
iteration 34 error = 1.2596419357016705e-08
iteration 35 error = 8.45716810361544e-09
iteration 36 error = 6.0288851026242664e-09
iteration 37 error = 4.681078727636483e-09
iteration 38 error = 3.571787167884541e-09
iteration 39 error = 2.1775545427642136e-09
iteration 40 error = 1.3847515545400096e-09
iteration 41 error = 9.5702307809913e-10
iteration 42 error = 6.401902122346198e-10
iteration 43 error = 4.4176223752100043e-10
iteration 44 error = 2.923142824134638e-10
iteration 45 error = 1.9063764595174233e-10
iteration 46 error = 1.2685031323783118e-10
iteration 47 error = 9.25616768749433e-11
iteration 48 error = 7.000424763695312e-11
iteration 49 error = 4.801207725381548e-11
iteration 50 error = 2.8073252219298783e-11
iteration 51 error = 2.072067385839785e-11
iteration 52 error = 1.3208168011459023e-11
iteration 53 error = 8.907523861270863e-12
iteration 54 error = 5.392265314270885e-12
iteration 55 error = 3.829585606793108e-12
iteration 56 error = 2.3572446553249417e-12
iteration 57 error = 1.5611089125073317e-12
iteration 58 error = 1.0214140200560482e-12
iteration 59 error = 6.905122619611176e-13
iteration 60 error = 4.408307361769625e-13
iteration 61 error = 2.908543818417076e-13
iteration 62 error = 1.8710839040593245e-13
iteration 63 error = 1.189436558642572e-13
iteration 64 error = 1.0177286297701935e-13
iteration 65 error = 7.201730730935995e-14
iteration 66 error = 4.585180959999352e-14
iteration 67 error = 2.84604895897364e-14
iteration 68 error = 1.854090448054995e-14
iteration 69 error = 1.277535390589116e-14
iteration 70 error = 7.874642446788625e-15
iteration 71 error = 5.309775678344668e-15
iteration 72 error = 3.586305940995534e-15
[7]:
# the vector potential is not supposed to look nice
Draw (gfu, mesh, "vector-potential", draw_surf=False)
Draw (curl(gfu), mesh, "B-field", draw_surf=False)
Draw (1/(mu0*mur)*curl(gfu)-mag, mesh, "H-field", draw_surf=False)