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.004770751268765439
iteration 1 error = 0.0024813274399324496
iteration 2 error = 0.0018918095766121933
iteration 3 error = 0.0012587090378891432
iteration 4 error = 0.0011364536233189856
iteration 5 error = 0.0008195457478849483
iteration 6 error = 0.0006060825336657866
iteration 7 error = 0.00045341288612874134
iteration 8 error = 0.00035785131847072924
iteration 9 error = 0.00023195184860900968
iteration 10 error = 0.00015491751841894308
iteration 11 error = 9.808156355275195e-05
iteration 12 error = 7.798886269681494e-05
iteration 13 error = 4.4312695788086745e-05
iteration 14 error = 3.282881496334988e-05
iteration 15 error = 2.2920087554377008e-05
iteration 16 error = 1.5012379279908464e-05
iteration 17 error = 1.029058627655539e-05
iteration 18 error = 6.797403084573732e-06
iteration 19 error = 5.324080434202822e-06
iteration 20 error = 3.093510781029749e-06
iteration 21 error = 2.273250129978486e-06
iteration 22 error = 1.5139195317352765e-06
iteration 23 error = 9.976809810274072e-07
iteration 24 error = 6.972869052375702e-07
iteration 25 error = 4.710001624204142e-07
iteration 26 error = 3.0898824460187017e-07
iteration 27 error = 2.12012130655872e-07
iteration 28 error = 1.4360554950978657e-07
iteration 29 error = 9.438532197572706e-08
iteration 30 error = 6.136548682445997e-08
iteration 31 error = 4.3555245111856965e-08
iteration 32 error = 2.9553593431748453e-08
iteration 33 error = 1.922027552414756e-08
iteration 34 error = 1.2599840426637963e-08
iteration 35 error = 8.495118540108491e-09
iteration 36 error = 6.288414564260635e-09
iteration 37 error = 5.186181699234099e-09
iteration 38 error = 3.4149114484062986e-09
iteration 39 error = 2.162174034429246e-09
iteration 40 error = 1.3842580873483715e-09
iteration 41 error = 9.570031865659943e-10
iteration 42 error = 6.402025394737825e-10
iteration 43 error = 4.4181610234789745e-10
iteration 44 error = 2.925666658074692e-10
iteration 45 error = 1.921475136555857e-10
iteration 46 error = 1.3783807444081773e-10
iteration 47 error = 1.0895713498539733e-10
iteration 48 error = 7.008132998254199e-11
iteration 49 error = 4.628832486895283e-11
iteration 50 error = 2.798560006875936e-11
iteration 51 error = 2.071501989744877e-11
iteration 52 error = 1.3207485491828499e-11
iteration 53 error = 8.90747723811717e-12
iteration 54 error = 5.3922633353490745e-12
iteration 55 error = 3.82958550820803e-12
iteration 56 error = 2.3572446498091386e-12
iteration 57 error = 1.561108914292894e-12
iteration 58 error = 1.0214140488849917e-12
iteration 59 error = 6.905126472884541e-13
iteration 60 error = 4.4083517145598173e-13
iteration 61 error = 2.9091940632831437e-13
iteration 62 error = 1.8805396230179367e-13
iteration 63 error = 1.3345054093655074e-13
iteration 64 error = 1.2312321679740248e-13
iteration 65 error = 6.601554195751265e-14
iteration 66 error = 4.5485434987956036e-14
iteration 67 error = 2.844956490926703e-14
iteration 68 error = 1.854149476825619e-14
iteration 69 error = 1.2778207966338202e-14
iteration 70 error = 7.886134426736732e-15
iteration 71 error = 5.371107187364038e-15
iteration 72 error = 3.760018959950903e-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)