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.003322319013244239
CG iteration 3, residual = 0.003311588317754007
CG iteration 4, residual = 0.0027467579779245745
CG iteration 5, residual = 0.0014658765270595296
CG iteration 6, residual = 0.0012167024999244268
CG iteration 7, residual = 0.0008096582309573242
CG iteration 8, residual = 0.0006570076838867986
CG iteration 9, residual = 0.00047508691834428543
CG iteration 10, residual = 0.00036224214074475757
CG iteration 11, residual = 0.00025441834248533015
CG iteration 12, residual = 0.00016194989094436854
CG iteration 13, residual = 0.00011358302702051946
CG iteration 14, residual = 8.950382818006961e-05
CG iteration 15, residual = 5.394126812252011e-05
CG iteration 16, residual = 3.916041785915708e-05
CG iteration 17, residual = 2.708231735171425e-05
CG iteration 18, residual = 1.8249603816404518e-05
CG iteration 19, residual = 1.3585772087464202e-05
CG iteration 20, residual = 9.901888497645214e-06
CG iteration 21, residual = 1.1880839567839292e-05
CG iteration 22, residual = 5.669457727842172e-06
CG iteration 23, residual = 3.593759937623263e-06
CG iteration 24, residual = 2.558907958839026e-06
CG iteration 25, residual = 1.977947784408537e-06
CG iteration 26, residual = 1.2770934679404039e-06
CG iteration 27, residual = 8.740832403648063e-07
CG iteration 28, residual = 5.927228746857347e-07
CG iteration 29, residual = 4.0089017034467715e-07
CG iteration 30, residual = 2.856458196958177e-07
CG iteration 31, residual = 2.0292685846583545e-07
CG iteration 32, residual = 1.319963698689081e-07
CG iteration 33, residual = 9.850315525587333e-08
CG iteration 34, residual = 6.734296914821343e-08
CG iteration 35, residual = 4.756031536156286e-08
CG iteration 36, residual = 2.9185851480466015e-08
CG iteration 37, residual = 2.6247277546711307e-08
CG iteration 38, residual = 2.3492272995949813e-08
CG iteration 39, residual = 1.3035488764833918e-08
CG iteration 40, residual = 8.395557033393525e-09
CG iteration 41, residual = 5.711979548739659e-09
CG iteration 42, residual = 4.1810498472084595e-09
CG iteration 43, residual = 2.794949339505362e-09
CG iteration 44, residual = 1.8595331195235108e-09
CG iteration 45, residual = 1.2654495193911474e-09
CG iteration 46, residual = 8.583836334145359e-10
CG iteration 47, residual = 5.766587076527566e-10
CG iteration 48, residual = 4.2177692351924746e-10
CG iteration 49, residual = 2.7611845100638133e-10
CG iteration 50, residual = 2.1705239475457787e-10
CG iteration 51, residual = 1.285950341926577e-10
CG iteration 52, residual = 8.65690783071117e-11
CG iteration 53, residual = 7.770869336953199e-11
CG iteration 54, residual = 6.80194819101193e-11
CG iteration 55, residual = 3.472078767948872e-11
CG iteration 56, residual = 2.323041073679045e-11
CG iteration 57, residual = 1.7339140993725202e-11
CG iteration 58, residual = 1.3511600576261893e-11
CG iteration 59, residual = 8.663977388327058e-12
CG iteration 60, residual = 5.942816176753264e-12
CG iteration 61, residual = 3.9442898796857385e-12
CG iteration 62, residual = 2.6957753850472994e-12
CG iteration 63, residual = 1.7387171148035777e-12
CG iteration 64, residual = 1.2104144492576381e-12
CG iteration 65, residual = 7.623361546652123e-13
CG iteration 66, residual = 4.941741149527372e-13
CG iteration 67, residual = 3.20852535831465e-13
CG iteration 68, residual = 2.2418213105360847e-13
CG iteration 69, residual = 3.08234096801784e-13
CG iteration 70, residual = 1.3636894357907712e-13
CG iteration 71, residual = 9.423119165842108e-14
CG iteration 72, residual = 6.13117281261854e-14
CG iteration 73, residual = 4.2069779919132366e-14
CG iteration 74, residual = 2.7628450476114503e-14
CG iteration 75, residual = 1.7668774240604515e-14
CG iteration 76, residual = 1.1673849512006252e-14
CG iteration 77, residual = 8.63354014860851e-15
CG iteration 78, residual = 6.747779694446461e-15
CG iteration 79, residual = 4.0565144680953554e-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);
[ ]:
[ ]: