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.occ 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
box = Box( (-3,-3,-3), (3,3,3))
box.faces.name = "outer"
magnet = Cylinder((-1,0,0),X, r=0.3, h=2)
magnet.mat("magnet")
magnet.faces.col = (1,0,0)
air = box-magnet
air.mat("air")
shape = Glue([air,magnet])
geo = OCCGeometry(shape)
Draw (shape, clipping={ "z" : -1, "function":True})
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 = 32867
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, printrates=True)
CG iteration 1, residual = 0.004809678530125145
CG iteration 2, residual = 0.003322319013244241
CG iteration 3, residual = 0.0033115883177540083
CG iteration 4, residual = 0.0027467579779245728
CG iteration 5, residual = 0.00146587652705953
CG iteration 6, residual = 0.0012167024999244287
CG iteration 7, residual = 0.0008096582309573251
CG iteration 8, residual = 0.0006570076838867993
CG iteration 9, residual = 0.0004750869183442854
CG iteration 10, residual = 0.0003622421407447582
CG iteration 11, residual = 0.0002544183424853302
CG iteration 12, residual = 0.00016194989094436862
CG iteration 13, residual = 0.00011358302702052046
CG iteration 14, residual = 8.950382818008721e-05
CG iteration 15, residual = 5.394126812319095e-05
CG iteration 16, residual = 3.916041789553318e-05
CG iteration 17, residual = 2.708231866259943e-05
CG iteration 18, residual = 1.8249665057440137e-05
CG iteration 19, residual = 1.3588245352141831e-05
CG iteration 20, residual = 9.96360154447062e-06
CG iteration 21, residual = 1.2089092056026526e-05
CG iteration 22, residual = 5.636710328408904e-06
CG iteration 23, residual = 3.5935932651152056e-06
CG iteration 24, residual = 2.558906361322061e-06
CG iteration 25, residual = 1.977947754737567e-06
CG iteration 26, residual = 1.277093467549514e-06
CG iteration 27, residual = 8.74083240361175e-07
CG iteration 28, residual = 5.927228746857039e-07
CG iteration 29, residual = 4.008901703446779e-07
CG iteration 30, residual = 2.8564581969581787e-07
CG iteration 31, residual = 2.029268584658104e-07
CG iteration 32, residual = 1.3199636986770275e-07
CG iteration 33, residual = 9.850315519832965e-08
CG iteration 34, residual = 6.734296711663225e-08
CG iteration 35, residual = 4.75602267782366e-08
CG iteration 36, residual = 2.9180913843619468e-08
CG iteration 37, residual = 2.599176578891142e-08
CG iteration 38, residual = 2.365585942345862e-08
CG iteration 39, residual = 1.303991860594903e-08
CG iteration 40, residual = 8.395587887984425e-09
CG iteration 41, residual = 5.711979758537775e-09
CG iteration 42, residual = 4.18104984946842e-09
CG iteration 43, residual = 2.7949493395595322e-09
CG iteration 44, residual = 1.8595331196815859e-09
CG iteration 45, residual = 1.2654495201624762e-09
CG iteration 46, residual = 8.583836376987338e-10
CG iteration 47, residual = 5.766587297210704e-10
CG iteration 48, residual = 4.217769955792741e-10
CG iteration 49, residual = 2.761187050102319e-10
CG iteration 50, residual = 2.1705281037339853e-10
CG iteration 51, residual = 1.2859400785337043e-10
CG iteration 52, residual = 8.64106781260756e-11
CG iteration 53, residual = 7.013235661458179e-11
CG iteration 54, residual = 7.414565263777497e-11
CG iteration 55, residual = 3.502643478758853e-11
CG iteration 56, residual = 2.495200597000085e-11
CG iteration 57, residual = 2.0726807222262598e-11
CG iteration 58, residual = 1.2946389153360395e-11
CG iteration 59, residual = 8.502058162524584e-12
CG iteration 60, residual = 5.9238848752983836e-12
CG iteration 61, residual = 3.942655481923114e-12
CG iteration 62, residual = 2.695620882984124e-12
CG iteration 63, residual = 1.7387080883596372e-12
CG iteration 64, residual = 1.2104137333708435e-12
CG iteration 65, residual = 7.623361044641524e-13
CG iteration 66, residual = 4.94173155656516e-13
CG iteration 67, residual = 3.2076250675189645e-13
CG iteration 68, residual = 2.167870085909234e-13
CG iteration 69, residual = 2.7803798937306976e-13
CG iteration 70, residual = 1.4128398861050609e-13
CG iteration 71, residual = 9.428172931076452e-14
CG iteration 72, residual = 6.131662021442724e-14
CG iteration 73, residual = 4.2086028935778823e-14
CG iteration 74, residual = 2.769237354854663e-14
CG iteration 75, residual = 1.8048913752474042e-14
CG iteration 76, residual = 1.3786064194428385e-14
CG iteration 77, residual = 9.968909599710185e-15
CG iteration 78, residual = 6.438126888149863e-15
CG iteration 79, residual = 3.968445197237799e-15
[7]:
Draw (curl(gfu), mesh, "B-field", draw_surf=False, \
clipping = { "z" : -1, "function":True}, \
vectors = { "grid_size":50}, min=0, max=2e-5);
[ ]:
[ ]: