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 ngsolve.webgui import Draw
from netgen.csg 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
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)
[2]:
<ngsolve.comp.Mesh at 0x7f0d4525c360>
[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 = 25068
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)
WARNING: maxsteps is deprecated, use maxiter instead!
CG iteration 1, residual = 0.0048325399628490456
CG iteration 2, residual = 0.0028801984652913493
CG iteration 3, residual = 0.0022636840717767817
CG iteration 4, residual = 0.0018596783443640282
CG iteration 5, residual = 0.0013706814431410205
CG iteration 6, residual = 0.0009629315328718021
CG iteration 7, residual = 0.0006784248172513831
CG iteration 8, residual = 0.00045373479988693595
CG iteration 9, residual = 0.00033739845235846664
CG iteration 10, residual = 0.00026636366966670483
CG iteration 11, residual = 0.0001578586571670998
CG iteration 12, residual = 0.00010596525810194424
CG iteration 13, residual = 7.440611757876551e-05
CG iteration 14, residual = 5.021776232098456e-05
CG iteration 15, residual = 3.5654571353773085e-05
CG iteration 16, residual = 2.1652775592701054e-05
CG iteration 17, residual = 1.624108858415096e-05
CG iteration 18, residual = 1.0173867299011198e-05
CG iteration 19, residual = 7.247797858031916e-06
CG iteration 20, residual = 4.7820876898995745e-06
CG iteration 21, residual = 3.2977746990720467e-06
CG iteration 22, residual = 2.314739131021156e-06
CG iteration 23, residual = 1.4401135174593936e-06
CG iteration 24, residual = 1.0156692744746578e-06
CG iteration 25, residual = 6.943663671593801e-07
CG iteration 26, residual = 5.60205009422576e-07
CG iteration 27, residual = 4.873246323834464e-07
CG iteration 28, residual = 2.834339597189828e-07
CG iteration 29, residual = 1.8372748823973257e-07
CG iteration 30, residual = 1.2576074145450456e-07
CG iteration 31, residual = 8.212387903172255e-08
CG iteration 32, residual = 5.655902159990351e-08
CG iteration 33, residual = 3.7186008437427496e-08
CG iteration 34, residual = 2.3956217034140503e-08
CG iteration 35, residual = 1.649825741954352e-08
CG iteration 36, residual = 1.0759738081355497e-08
CG iteration 37, residual = 7.16506359022578e-09
CG iteration 38, residual = 4.497597041051191e-09
CG iteration 39, residual = 2.9156780764484104e-09
CG iteration 40, residual = 1.946000989702628e-09
CG iteration 41, residual = 1.4015970373316643e-09
CG iteration 42, residual = 8.454927601870392e-10
CG iteration 43, residual = 5.886122205394698e-10
CG iteration 44, residual = 3.7119789980860563e-10
CG iteration 45, residual = 2.5506242491834047e-10
CG iteration 46, residual = 2.2081302274063233e-10
CG iteration 47, residual = 1.7818556763699376e-10
CG iteration 48, residual = 1.310330793750153e-10
CG iteration 49, residual = 8.618066083500679e-11
CG iteration 50, residual = 5.432591771444806e-11
CG iteration 51, residual = 3.644418530879367e-11
CG iteration 52, residual = 2.364366419611503e-11
CG iteration 53, residual = 1.6321065493437588e-11
CG iteration 54, residual = 1.0677533127223656e-11
CG iteration 55, residual = 7.155944541936287e-12
CG iteration 56, residual = 4.7963207050367265e-12
CG iteration 57, residual = 3.098941060108693e-12
CG iteration 58, residual = 2.0500215060812245e-12
CG iteration 59, residual = 1.3892675227012812e-12
CG iteration 60, residual = 8.940759449286533e-13
CG iteration 61, residual = 5.967417279404007e-13
CG iteration 62, residual = 3.9036915425913924e-13
CG iteration 63, residual = 2.578187071195877e-13
CG iteration 64, residual = 1.7608681449249252e-13
CG iteration 65, residual = 1.1437174582179819e-13
CG iteration 66, residual = 8.220275845728704e-14
CG iteration 67, residual = 8.613376955958004e-14
CG iteration 68, residual = 4.6817117474089865e-14
CG iteration 69, residual = 3.00772839094338e-14
CG iteration 70, residual = 2.0408341316525554e-14
CG iteration 71, residual = 1.4867469352209805e-14
CG iteration 72, residual = 1.0277312538739666e-14
CG iteration 73, residual = 6.278554627068153e-15
CG iteration 74, residual = 4.175264377281126e-15
[7]:
# the vector potential is not supposed to look nice
Draw (gfu, mesh, "vector-potential", draw_surf=False, clipping=True)
Draw (curl(gfu), mesh, "B-field", draw_surf=False, clipping=True)
Draw (1/(mu0*mur)*curl(gfu)-mag, mesh, "H-field", draw_surf=False, clipping=True)
[7]:
BaseWebGuiScene
[ ]: