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 = 33152
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.004846766222292576
CG iteration 2, residual = 0.0040639697315015535
CG iteration 3, residual = 0.003290581389748042
CG iteration 4, residual = 0.002260752973944027
CG iteration 5, residual = 0.0019757898056008724
CG iteration 6, residual = 0.0011507895291099635
CG iteration 7, residual = 0.0008639313194094742
CG iteration 8, residual = 0.0006327449629420519
CG iteration 9, residual = 0.0004949733354845924
CG iteration 10, residual = 0.00040760541241431484
CG iteration 11, residual = 0.0003164866826507574
CG iteration 12, residual = 0.00018393955806005895
CG iteration 13, residual = 0.00013713589121294687
CG iteration 14, residual = 8.753130614903451e-05
CG iteration 15, residual = 6.404693879077445e-05
CG iteration 16, residual = 4.235156603880363e-05
CG iteration 17, residual = 2.6160851585889872e-05
CG iteration 18, residual = 1.953985148853433e-05
CG iteration 19, residual = 1.3274429462562e-05
CG iteration 20, residual = 1.807506432019514e-05
CG iteration 21, residual = 8.096218322633755e-06
CG iteration 22, residual = 5.688075482269196e-06
CG iteration 23, residual = 3.828302736606345e-06
CG iteration 24, residual = 2.512608014418754e-06
CG iteration 25, residual = 1.9133034159980976e-06
CG iteration 26, residual = 1.1882807370649346e-06
CG iteration 27, residual = 8.414430273125002e-07
CG iteration 28, residual = 5.733998486020003e-07
CG iteration 29, residual = 3.8563161445226834e-07
CG iteration 30, residual = 2.8128804922057385e-07
CG iteration 31, residual = 1.7903264909374636e-07
CG iteration 32, residual = 1.3336135486635135e-07
CG iteration 33, residual = 1.157497715545079e-07
CG iteration 34, residual = 7.859128325743426e-08
CG iteration 35, residual = 4.9300053017296735e-08
CG iteration 36, residual = 3.8138752187723044e-08
CG iteration 37, residual = 4.676975196162582e-08
CG iteration 38, residual = 2.2171123364103118e-08
CG iteration 39, residual = 1.4465678089806094e-08
CG iteration 40, residual = 9.548079569467562e-09
CG iteration 41, residual = 7.757887875976946e-09
CG iteration 42, residual = 6.5141087411272525e-09
CG iteration 43, residual = 3.6398645675411454e-09
CG iteration 44, residual = 2.622886829395068e-09
CG iteration 45, residual = 1.7228721505528968e-09
CG iteration 46, residual = 1.1298896049497874e-09
CG iteration 47, residual = 7.954171210878555e-10
CG iteration 48, residual = 5.005043999958202e-10
CG iteration 49, residual = 3.505160812737958e-10
CG iteration 50, residual = 2.270450219562435e-10
CG iteration 51, residual = 1.5407658494048545e-10
CG iteration 52, residual = 1.0420159579312062e-10
CG iteration 53, residual = 1.5067939886783798e-10
CG iteration 54, residual = 6.20047694146377e-11
CG iteration 55, residual = 4.032391530282409e-11
CG iteration 56, residual = 2.76489474596565e-11
CG iteration 57, residual = 1.7867643088094327e-11
CG iteration 58, residual = 1.2085719421074547e-11
CG iteration 59, residual = 7.876886307810658e-12
CG iteration 60, residual = 5.99814878854235e-12
CG iteration 61, residual = 5.396805971549421e-12
CG iteration 62, residual = 3.0552168641231883e-12
CG iteration 63, residual = 1.9769921868818155e-12
CG iteration 64, residual = 1.2730011732667527e-12
CG iteration 65, residual = 8.764662276002234e-13
CG iteration 66, residual = 6.018715775829037e-13
CG iteration 67, residual = 3.7299205236739034e-13
CG iteration 68, residual = 3.2761561587167917e-13
CG iteration 69, residual = 3.086959282267892e-13
CG iteration 70, residual = 1.8596893018177687e-13
CG iteration 71, residual = 1.6599137023894148e-13
CG iteration 72, residual = 9.039428906733414e-14
CG iteration 73, residual = 5.8651738780523e-14
CG iteration 74, residual = 4.024731209272248e-14
CG iteration 75, residual = 2.6754025903216562e-14
CG iteration 76, residual = 1.7722758516246538e-14
CG iteration 77, residual = 1.1619583550134709e-14
CG iteration 78, residual = 7.602385397051636e-15
CG iteration 79, residual = 4.828420840630908e-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);
[ ]:
[ ]: