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 = 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 = 25045
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.004810377047252341
iteration 1 error = 0.0023760495421453337
iteration 2 error = 0.0020091329922634353
iteration 3 error = 0.0018347748421582242
iteration 4 error = 0.0011702983943112177
iteration 5 error = 0.0007994173970521511
iteration 6 error = 0.0006047754213553882
iteration 7 error = 0.0004484730179567566
iteration 8 error = 0.00033717395608943193
iteration 9 error = 0.00021830155092886124
iteration 10 error = 0.00016470550887057335
iteration 11 error = 0.00010632701022515557
iteration 12 error = 6.970620403330886e-05
iteration 13 error = 4.9560320155102976e-05
iteration 14 error = 3.2129127188838485e-05
iteration 15 error = 2.2500302882812562e-05
iteration 16 error = 1.4624085914637495e-05
iteration 17 error = 9.992986502210904e-06
iteration 18 error = 6.8146244246295535e-06
iteration 19 error = 4.660924322297562e-06
iteration 20 error = 2.9828719291827927e-06
iteration 21 error = 2.071521893949705e-06
iteration 22 error = 1.4187542219683928e-06
iteration 23 error = 9.367188301894539e-07
iteration 24 error = 6.236666144351196e-07
iteration 25 error = 4.731556889502696e-07
iteration 26 error = 4.680578138503136e-07
iteration 27 error = 2.5162951585000246e-07
iteration 28 error = 1.7324839707902926e-07
iteration 29 error = 1.1155422114024779e-07
iteration 30 error = 7.144841979008892e-08
iteration 31 error = 4.7826645592949286e-08
iteration 32 error = 3.116814111090041e-08
iteration 33 error = 2.0859579896370473e-08
iteration 34 error = 1.4126543325152402e-08
iteration 35 error = 8.737068456831768e-09
iteration 36 error = 5.818889500393854e-09
iteration 37 error = 3.964936296984502e-09
iteration 38 error = 2.7245842040946543e-09
iteration 39 error = 1.6443902164812437e-09
iteration 40 error = 1.0904314859880575e-09
iteration 41 error = 7.478234533904745e-10
iteration 42 error = 4.955051428643892e-10
iteration 43 error = 3.3932299424693045e-10
iteration 44 error = 2.6894543286645653e-10
iteration 45 error = 1.843391540195826e-10
iteration 46 error = 1.3014582835786054e-10
iteration 47 error = 1.2936993028524715e-10
iteration 48 error = 7.139927779159655e-11
iteration 49 error = 4.8181094875929546e-11
iteration 50 error = 3.026094351359654e-11
iteration 51 error = 2.136345870445469e-11
iteration 52 error = 1.4033113336067073e-11
iteration 53 error = 8.783192723278603e-12
iteration 54 error = 5.976803003112857e-12
iteration 55 error = 4.048720286106233e-12
iteration 56 error = 2.617667761759695e-12
iteration 57 error = 1.6758227351954219e-12
iteration 58 error = 1.1760797135391065e-12
iteration 59 error = 7.80804126494363e-13
iteration 60 error = 4.882351296887458e-13
iteration 61 error = 3.200437844043975e-13
iteration 62 error = 2.23793044491635e-13
iteration 63 error = 1.4263266677986444e-13
iteration 64 error = 8.99066773690261e-14
iteration 65 error = 5.911960368574551e-14
iteration 66 error = 4.5568653943696326e-14
iteration 67 error = 4.5424442495257585e-14
iteration 68 error = 2.441542786588041e-14
iteration 69 error = 1.625748789182941e-14
iteration 70 error = 1.0620979998821664e-14
iteration 71 error = 6.547288162939496e-15
iteration 72 error = 4.251824211003897e-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)