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.0048096785301251456
CG iteration 2, residual = 0.003322319013244242
CG iteration 3, residual = 0.0033115883177540105
CG iteration 4, residual = 0.002746757977924575
CG iteration 5, residual = 0.0014658765270595313
CG iteration 6, residual = 0.0012167024999244277
CG iteration 7, residual = 0.0008096582309573243
CG iteration 8, residual = 0.000657007683886798
CG iteration 9, residual = 0.00047508691834428484
CG iteration 10, residual = 0.0003622421407447577
CG iteration 11, residual = 0.0002544183424853302
CG iteration 12, residual = 0.00016194989094436865
CG iteration 13, residual = 0.0001135830270205203
CG iteration 14, residual = 8.95038281800762e-05
CG iteration 15, residual = 5.394126812277027e-05
CG iteration 16, residual = 3.916041787272139e-05
CG iteration 17, residual = 2.7082317840531823e-05
CG iteration 18, residual = 1.824962665269008e-05
CG iteration 19, residual = 1.3586694453239409e-05
CG iteration 20, residual = 9.925007020263908e-06
CG iteration 21, residual = 1.1968145608492963e-05
CG iteration 22, residual = 5.656007731807217e-06
CG iteration 23, residual = 3.59369119589875e-06
CG iteration 24, residual = 2.5589072999408923e-06
CG iteration 25, residual = 1.9779477721706918e-06
CG iteration 26, residual = 1.2770934677791828e-06
CG iteration 27, residual = 8.74083240363311e-07
CG iteration 28, residual = 5.927228746857242e-07
CG iteration 29, residual = 4.0089017034467954e-07
CG iteration 30, residual = 2.8564581969581755e-07
CG iteration 31, residual = 2.0292685846573618e-07
CG iteration 32, residual = 1.319963698641989e-07
CG iteration 33, residual = 9.850315503109854e-08
CG iteration 34, residual = 6.734296121252867e-08
CG iteration 35, residual = 4.75599693390351e-08
CG iteration 36, residual = 2.916655656887256e-08
CG iteration 37, residual = 2.5216455114836116e-08
CG iteration 38, residual = 2.4192871094006452e-08
CG iteration 39, residual = 1.3056287590018261e-08
CG iteration 40, residual = 8.395702148158744e-09
CG iteration 41, residual = 5.711980535462875e-09
CG iteration 42, residual = 4.18104985778961e-09
CG iteration 43, residual = 2.794949339637109e-09
CG iteration 44, residual = 1.8595331196096027e-09
CG iteration 45, residual = 1.2654495198076316e-09
CG iteration 46, residual = 8.583836357277361e-10
CG iteration 47, residual = 5.766587195680032e-10
CG iteration 48, residual = 4.21776962417377e-10
CG iteration 49, residual = 2.7611858768799715e-10
CG iteration 50, residual = 2.1705260236090234e-10
CG iteration 51, residual = 1.2859383149896257e-10
CG iteration 52, residual = 8.642763345896675e-11
CG iteration 53, residual = 7.116059911067164e-11
CG iteration 54, residual = 7.312162153021118e-11
CG iteration 55, residual = 3.490819945039616e-11
CG iteration 56, residual = 2.4218272296030902e-11
CG iteration 57, residual = 2.007634876177969e-11
CG iteration 58, residual = 1.3135223703013364e-11
CG iteration 59, residual = 8.520921959614959e-12
CG iteration 60, residual = 5.925775017434715e-12
CG iteration 61, residual = 3.942815096809531e-12
CG iteration 62, residual = 2.6956359302845585e-12
CG iteration 63, residual = 1.7387089671975686e-12
CG iteration 64, residual = 1.210413803027969e-12
CG iteration 65, residual = 7.623361103354899e-13
CG iteration 66, residual = 4.941733659970352e-13
CG iteration 67, residual = 3.2078243705083564e-13
CG iteration 68, residual = 2.1845938363358796e-13
CG iteration 69, residual = 2.9265236269145245e-13
CG iteration 70, residual = 1.3917135785898355e-13
CG iteration 71, residual = 9.42585219340617e-14
CG iteration 72, residual = 6.131173994340097e-14
CG iteration 73, residual = 4.206912343187545e-14
CG iteration 74, residual = 2.7625882783016854e-14
CG iteration 75, residual = 1.7653035220735684e-14
CG iteration 76, residual = 1.1546271341885469e-14
CG iteration 77, residual = 8.107554902988951e-15
CG iteration 78, residual = 6.411198117609457e-15
CG iteration 79, residual = 4.173583661506341e-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);
[ ]:
[ ]: