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 = 23791
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.004855037018856693
iteration 1 error = 0.0027293684352409295
iteration 2 error = 0.0019178169752064114
iteration 3 error = 0.0011868200578276746
iteration 4 error = 0.0008381444433524467
iteration 5 error = 0.0005352958317220569
iteration 6 error = 0.00041878588073819925
iteration 7 error = 0.000249667292203067
iteration 8 error = 0.0001479444156746633
iteration 9 error = 9.829281774119941e-05
iteration 10 error = 6.796526225372933e-05
iteration 11 error = 3.914655734049888e-05
iteration 12 error = 2.282461683183525e-05
iteration 13 error = 1.4024099850047308e-05
iteration 14 error = 9.109914719968633e-06
iteration 15 error = 5.472491313755589e-06
iteration 16 error = 3.332289562794918e-06
iteration 17 error = 1.972834212002664e-06
iteration 18 error = 1.187155034288811e-06
iteration 19 error = 7.339626800917817e-07
iteration 20 error = 4.3623176671614813e-07
iteration 21 error = 2.56683635918861e-07
iteration 22 error = 1.6074742788743125e-07
iteration 23 error = 1.0206447162818972e-07
iteration 24 error = 6.217499592128528e-08
iteration 25 error = 3.7265378772296614e-08
iteration 26 error = 2.1815012102851695e-08
iteration 27 error = 1.2445629160114056e-08
iteration 28 error = 7.3237410812180256e-09
iteration 29 error = 4.8530457068053855e-09
iteration 30 error = 3.1476230725435892e-09
iteration 31 error = 1.7144938912281451e-09
iteration 32 error = 1.0204511653757055e-09
iteration 33 error = 6.28515710247334e-10
iteration 34 error = 3.820108113438369e-10
iteration 35 error = 2.2615096962387014e-10
iteration 36 error = 1.333427404373648e-10
iteration 37 error = 7.592574238063816e-11
iteration 38 error = 4.462595643695611e-11
iteration 39 error = 2.5695855900998983e-11
iteration 40 error = 1.542995725711371e-11
iteration 41 error = 9.389387349130117e-12
iteration 42 error = 6.577404480574253e-12
iteration 43 error = 4.4127351188302396e-12
iteration 44 error = 2.5439829150772746e-12
iteration 45 error = 1.4755648258044578e-12
iteration 46 error = 9.042743785371341e-13
iteration 47 error = 5.394734054737619e-13
iteration 48 error = 2.989292558044478e-13
iteration 49 error = 1.7566670797841405e-13
iteration 50 error = 1.1544103553339943e-13
iteration 51 error = 8.387649737687494e-14
iteration 52 error = 4.4982603585713545e-14
iteration 53 error = 2.5981353732921825e-14
iteration 54 error = 1.5284334927115574e-14
iteration 55 error = 8.713979005543354e-15
iteration 56 error = 5.9904207437133914e-15
iteration 57 error = 3.9910358488320136e-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)