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.0033223190132442492
CG iteration 3, residual = 0.0033115883177540157
CG iteration 4, residual = 0.00274675797792458
CG iteration 5, residual = 0.0014658765270595341
CG iteration 6, residual = 0.0012167024999244316
CG iteration 7, residual = 0.0008096582309573266
CG iteration 8, residual = 0.0006570076838868005
CG iteration 9, residual = 0.000475086918344287
CG iteration 10, residual = 0.00036224214074475903
CG iteration 11, residual = 0.0002544183424853308
CG iteration 12, residual = 0.00016194989094436878
CG iteration 13, residual = 0.00011358302702051957
CG iteration 14, residual = 8.95038281800596e-05
CG iteration 15, residual = 5.3941268122132255e-05
CG iteration 16, residual = 3.916041783812964e-05
CG iteration 17, residual = 2.708231659395375e-05
CG iteration 18, residual = 1.8249568415699526e-05
CG iteration 19, residual = 1.3584341992335776e-05
CG iteration 20, residual = 9.865797328223227e-06
CG iteration 21, residual = 1.1717898658226638e-05
CG iteration 22, residual = 5.694098471406397e-06
CG iteration 23, residual = 3.593886920834201e-06
CG iteration 24, residual = 2.5589091760807176e-06
CG iteration 25, residual = 1.9779478070166387e-06
CG iteration 26, residual = 1.2770934682382533e-06
CG iteration 27, residual = 8.740832403675801e-07
CG iteration 28, residual = 5.927228746857623e-07
CG iteration 29, residual = 4.00890170344678e-07
CG iteration 30, residual = 2.856458196958077e-07
CG iteration 31, residual = 2.0292685846545208e-07
CG iteration 32, residual = 1.3199636985095956e-07
CG iteration 33, residual = 9.850315439939665e-08
CG iteration 34, residual = 6.734293891043426e-08
CG iteration 35, residual = 4.75589968675574e-08
CG iteration 36, residual = 2.9112219827966917e-08
CG iteration 37, residual = 2.1743872712007513e-08
CG iteration 38, residual = 2.6338740874032187e-08
CG iteration 39, residual = 1.3369873480006615e-08
CG iteration 40, residual = 8.397968936764506e-09
CG iteration 41, residual = 5.71199595431689e-09
CG iteration 42, residual = 4.181050023060027e-09
CG iteration 43, residual = 2.7949493415055176e-09
CG iteration 44, residual = 1.85953311993934e-09
CG iteration 45, residual = 1.2654495213556228e-09
CG iteration 46, residual = 8.583836443253641e-10
CG iteration 47, residual = 5.766587638558959e-10
CG iteration 48, residual = 4.217771070533169e-10
CG iteration 49, residual = 2.76119098576979e-10
CG iteration 50, residual = 2.1705347801330806e-10
CG iteration 51, residual = 1.2859337764573906e-10
CG iteration 52, residual = 8.62478069499179e-11
CG iteration 53, residual = 5.98846158421347e-11
CG iteration 54, residual = 7.570495284578722e-11
CG iteration 55, residual = 3.6940483222420186e-11
CG iteration 56, residual = 2.6896637747628242e-11
CG iteration 57, residual = 2.0926826732675046e-11
CG iteration 58, residual = 1.2729199207393305e-11
CG iteration 59, residual = 8.48530776962353e-12
CG iteration 60, residual = 5.92225979307696e-12
CG iteration 61, residual = 3.942518859594153e-12
CG iteration 62, residual = 2.695608010318349e-12
CG iteration 63, residual = 1.7387073365484449e-12
CG iteration 64, residual = 1.2104136736248744e-12
CG iteration 65, residual = 7.623360953366181e-13
CG iteration 66, residual = 4.94172640950329e-13
CG iteration 67, residual = 3.2071440052959655e-13
CG iteration 68, residual = 2.1267062868871536e-13
CG iteration 69, residual = 1.5969897524178984e-13
CG iteration 70, residual = 1.9241927879092297e-13
CG iteration 71, residual = 9.596240175248145e-14
CG iteration 72, residual = 6.132528870034374e-14
CG iteration 73, residual = 4.207245136372843e-14
CG iteration 74, residual = 2.7638618382865565e-14
CG iteration 75, residual = 1.7731286613846187e-14
CG iteration 76, residual = 1.2142714506756057e-14
CG iteration 77, residual = 9.608898203456885e-15
CG iteration 78, residual = 6.65147721210909e-15
CG iteration 79, residual = 3.991460839739057e-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);
[ ]:
[ ]: