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.004809678530125144
CG iteration 2, residual = 0.003322319013244251
CG iteration 3, residual = 0.003311588317754024
CG iteration 4, residual = 0.002746757977924581
CG iteration 5, residual = 0.0014658765270595348
CG iteration 6, residual = 0.0012167024999244322
CG iteration 7, residual = 0.0008096582309573271
CG iteration 8, residual = 0.0006570076838868006
CG iteration 9, residual = 0.0004750869183442866
CG iteration 10, residual = 0.0003622421407447591
CG iteration 11, residual = 0.0002544183424853311
CG iteration 12, residual = 0.0001619498909443691
CG iteration 13, residual = 0.00011358302702052068
CG iteration 14, residual = 8.950382818010079e-05
CG iteration 15, residual = 5.394126812369696e-05
CG iteration 16, residual = 3.9160417922962924e-05
CG iteration 17, residual = 2.7082319651085978e-05
CG iteration 18, residual = 1.82497112366325e-05
CG iteration 19, residual = 1.3590109736917365e-05
CG iteration 20, residual = 1.000954246646836e-05
CG iteration 21, residual = 1.2199556387211128e-05
CG iteration 22, residual = 5.617743871009024e-06
CG iteration 23, residual = 3.5934978054460547e-06
CG iteration 24, residual = 2.5589054464556036e-06
CG iteration 25, residual = 1.977947737745601e-06
CG iteration 26, residual = 1.2770934673256553e-06
CG iteration 27, residual = 8.740832403590926e-07
CG iteration 28, residual = 5.92722874685685e-07
CG iteration 29, residual = 4.0089017034467827e-07
CG iteration 30, residual = 2.8564581969582374e-07
CG iteration 31, residual = 2.0292685846602097e-07
CG iteration 32, residual = 1.319963698775578e-07
CG iteration 33, residual = 9.85031556685903e-08
CG iteration 34, residual = 6.734298371909184e-08
CG iteration 35, residual = 4.7560950688458965e-08
CG iteration 36, residual = 2.922122542146771e-08
CG iteration 37, residual = 2.793152057140147e-08
CG iteration 38, residual = 2.256154568206488e-08
CG iteration 39, residual = 1.3013912792235378e-08
CG iteration 40, residual = 8.395407159153344e-09
CG iteration 41, residual = 5.711978529702247e-09
CG iteration 42, residual = 4.181049836292295e-09
CG iteration 43, residual = 2.7949493393987395e-09
CG iteration 44, residual = 1.8595331195918353e-09
CG iteration 45, residual = 1.2654495197290788e-09
CG iteration 46, residual = 8.583836352915739e-10
CG iteration 47, residual = 5.766587173225386e-10
CG iteration 48, residual = 4.217769551403914e-10
CG iteration 49, residual = 2.761185647288493e-10
CG iteration 50, residual = 2.170526651845165e-10
CG iteration 51, residual = 1.2859799422356132e-10
CG iteration 52, residual = 8.67930069222062e-11
CG iteration 53, residual = 8.563472952334675e-11
CG iteration 54, residual = 6.401579525193446e-11
CG iteration 55, residual = 3.4776559389365203e-11
CG iteration 56, residual = 2.4042635040091006e-11
CG iteration 57, residual = 1.9808466653816863e-11
CG iteration 58, residual = 1.320011450011249e-11
CG iteration 59, residual = 8.528885298512687e-12
CG iteration 60, residual = 5.926593138543212e-12
CG iteration 61, residual = 3.942884422021349e-12
CG iteration 62, residual = 2.6956424685262465e-12
CG iteration 63, residual = 1.7387093490918876e-12
CG iteration 64, residual = 1.2104138336109068e-12
CG iteration 65, residual = 7.623361337242638e-13
CG iteration 66, residual = 4.941753815515835e-13
CG iteration 67, residual = 3.209714106764555e-13
CG iteration 68, residual = 2.333815165013505e-13
CG iteration 69, residual = 3.020211759353333e-13
CG iteration 70, residual = 1.3502404389380167e-13
CG iteration 71, residual = 9.421856120851726e-14
CG iteration 72, residual = 6.131149089424588e-14
CG iteration 73, residual = 4.206926694251458e-14
CG iteration 74, residual = 2.7626466126513066e-14
CG iteration 75, residual = 1.765668840130954e-14
CG iteration 76, residual = 1.1576815002459953e-14
CG iteration 77, residual = 8.25338171094234e-15
CG iteration 78, residual = 6.5823159194537296e-15
CG iteration 79, residual = 4.126329519302202e-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);
[ ]:
[ ]: