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 = 37684
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='\r')
CG iteration 1, residual = 0.004850973085913294
CG iteration 2, residual = 0.005070311950390225 CG iteration 3, residual = 0.004061148250617966 CG iteration 4, residual = 0.0027848857704494343 CG iteration 5, residual = 0.0021869996822521747 CG iteration 6, residual = 0.0013947920765419103 CG iteration 7, residual = 0.000991067909052922 CG iteration 8, residual = 0.0007882253802487303 CG iteration 9, residual = 0.000649062902879689 CG iteration 10, residual = 0.0005617065918108688 CG iteration 11, residual = 0.00035196344676141275 CG iteration 12, residual = 0.00022872339864623404 CG iteration 13, residual = 0.00016539023229782408 CG iteration 14, residual = 0.00012878164665465646 CG iteration 15, residual = 9.28716269596758e-05 CG iteration 16, residual = 5.482801476858597e-05 CG iteration 17, residual = 5.767787533886232e-05 CG iteration 18, residual = 4.225659177147738e-05 CG iteration 19, residual = 2.546438006603943e-05 CG iteration 20, residual = 1.7523713578956545e-05 CG iteration 21, residual = 1.1633918805191769e-05 CG iteration 22, residual = 7.830536631448518e-06 CG iteration 23, residual = 5.55087703593193e-06 CG iteration 24, residual = 3.7356806435012735e-06 CG iteration 25, residual = 2.5276590244645992e-06 CG iteration 26, residual = 1.7668766351352113e-06 CG iteration 27, residual = 1.1776767494697926e-06 CG iteration 28, residual = 7.943753716959176e-07 CG iteration 29, residual = 5.81995856292263e-07 CG iteration 30, residual = 4.556030490488636e-07 CG iteration 31, residual = 5.969727942688124e-07 CG iteration 32, residual = 3.53413336987771e-07 CG iteration 33, residual = 2.3165170433186799e-07 CG iteration 34, residual = 1.523579904274309e-07 CG iteration 35, residual = 1.0054466510145232e-07 CG iteration 36, residual = 6.714545175168269e-08 CG iteration 37, residual = 4.963279352688796e-08 CG iteration 38, residual = 3.077821007387553e-08 CG iteration 39, residual = 2.2910434477570848e-08 CG iteration 40, residual = 1.4896517137350837e-08 CG iteration 41, residual = 1.0767665973717847e-08 CG iteration 42, residual = 9.559394555030711e-09 CG iteration 43, residual = 6.946145871709278e-09 CG iteration 44, residual = 4.2565830452215405e-09 CG iteration 45, residual = 3.7076262948446738e-09 CG iteration 46, residual = 3.995964482657597e-09 CG iteration 47, residual = 1.9905656697202943e-09 CG iteration 48, residual = 1.3307392471866413e-09 CG iteration 49, residual = 9.182664330218499e-10 CG iteration 50, residual = 6.413685294769359e-10 CG iteration 51, residual = 4.2320693272587554e-10 CG iteration 52, residual = 2.982873067494714e-10
</pre>
CG iteration 2, residual = 0.005070311950390225 CG iteration 3, residual = 0.004061148250617966 CG iteration 4, residual = 0.0027848857704494343 CG iteration 5, residual = 0.0021869996822521747 CG iteration 6, residual = 0.0013947920765419103 CG iteration 7, residual = 0.000991067909052922 CG iteration 8, residual = 0.0007882253802487303 CG iteration 9, residual = 0.000649062902879689 CG iteration 10, residual = 0.0005617065918108688 CG iteration 11, residual = 0.00035196344676141275 CG iteration 12, residual = 0.00022872339864623404 CG iteration 13, residual = 0.00016539023229782408 CG iteration 14, residual = 0.00012878164665465646 CG iteration 15, residual = 9.28716269596758e-05 CG iteration 16, residual = 5.482801476858597e-05 CG iteration 17, residual = 5.767787533886232e-05 CG iteration 18, residual = 4.225659177147738e-05 CG iteration 19, residual = 2.546438006603943e-05 CG iteration 20, residual = 1.7523713578956545e-05 CG iteration 21, residual = 1.1633918805191769e-05 CG iteration 22, residual = 7.830536631448518e-06 CG iteration 23, residual = 5.55087703593193e-06 CG iteration 24, residual = 3.7356806435012735e-06 CG iteration 25, residual = 2.5276590244645992e-06 CG iteration 26, residual = 1.7668766351352113e-06 CG iteration 27, residual = 1.1776767494697926e-06 CG iteration 28, residual = 7.943753716959176e-07 CG iteration 29, residual = 5.81995856292263e-07 CG iteration 30, residual = 4.556030490488636e-07 CG iteration 31, residual = 5.969727942688124e-07 CG iteration 32, residual = 3.53413336987771e-07 CG iteration 33, residual = 2.3165170433186799e-07 CG iteration 34, residual = 1.523579904274309e-07 CG iteration 35, residual = 1.0054466510145232e-07 CG iteration 36, residual = 6.714545175168269e-08 CG iteration 37, residual = 4.963279352688796e-08 CG iteration 38, residual = 3.077821007387553e-08 CG iteration 39, residual = 2.2910434477570848e-08 CG iteration 40, residual = 1.4896517137350837e-08 CG iteration 41, residual = 1.0767665973717847e-08 CG iteration 42, residual = 9.559394555030711e-09 CG iteration 43, residual = 6.946145871709278e-09 CG iteration 44, residual = 4.2565830452215405e-09 CG iteration 45, residual = 3.7076262948446738e-09 CG iteration 46, residual = 3.995964482657597e-09 CG iteration 47, residual = 1.9905656697202943e-09 CG iteration 48, residual = 1.3307392471866413e-09 CG iteration 49, residual = 9.182664330218499e-10 CG iteration 50, residual = 6.413685294769359e-10 CG iteration 51, residual = 4.2320693272587554e-10 CG iteration 52, residual = 2.982873067494714e-10
end{sphinxVerbatim}
[2KCG iteration 2, residual = 0.005070311950390225 [2KCG iteration 3, residual = 0.004061148250617966 [2KCG iteration 4, residual = 0.0027848857704494343 [2KCG iteration 5, residual = 0.0021869996822521747 [2KCG iteration 6, residual = 0.0013947920765419103 [2KCG iteration 7, residual = 0.000991067909052922 [2KCG iteration 8, residual = 0.0007882253802487303 [2KCG iteration 9, residual = 0.000649062902879689 [2KCG iteration 10, residual = 0.0005617065918108688 [2KCG iteration 11, residual = 0.00035196344676141275 [2KCG iteration 12, residual = 0.00022872339864623404 [2KCG iteration 13, residual = 0.00016539023229782408 [2KCG iteration 14, residual = 0.00012878164665465646 [2KCG iteration 15, residual = 9.28716269596758e-05 [2KCG iteration 16, residual = 5.482801476858597e-05 [2KCG iteration 17, residual = 5.767787533886232e-05 [2KCG iteration 18, residual = 4.225659177147738e-05 [2KCG iteration 19, residual = 2.546438006603943e-05 [2KCG iteration 20, residual = 1.7523713578956545e-05 [2KCG iteration 21, residual = 1.1633918805191769e-05 [2KCG iteration 22, residual = 7.830536631448518e-06 [2KCG iteration 23, residual = 5.55087703593193e-06 [2KCG iteration 24, residual = 3.7356806435012735e-06 [2KCG iteration 25, residual = 2.5276590244645992e-06 [2KCG iteration 26, residual = 1.7668766351352113e-06 [2KCG iteration 27, residual = 1.1776767494697926e-06 [2KCG iteration 28, residual = 7.943753716959176e-07 [2KCG iteration 29, residual = 5.81995856292263e-07 [2KCG iteration 30, residual = 4.556030490488636e-07 [2KCG iteration 31, residual = 5.969727942688124e-07 [2KCG iteration 32, residual = 3.53413336987771e-07 [2KCG iteration 33, residual = 2.3165170433186799e-07 [2KCG iteration 34, residual = 1.523579904274309e-07 [2KCG iteration 35, residual = 1.0054466510145232e-07 [2KCG iteration 36, residual = 6.714545175168269e-08 [2KCG iteration 37, residual = 4.963279352688796e-08 [2KCG iteration 38, residual = 3.077821007387553e-08 [2KCG iteration 39, residual = 2.2910434477570848e-08 [2KCG iteration 40, residual = 1.4896517137350837e-08 [2KCG iteration 41, residual = 1.0767665973717847e-08 [2KCG iteration 42, residual = 9.559394555030711e-09 [2KCG iteration 43, residual = 6.946145871709278e-09 [2KCG iteration 44, residual = 4.2565830452215405e-09 [2KCG iteration 45, residual = 3.7076262948446738e-09 [2KCG iteration 46, residual = 3.995964482657597e-09 [2KCG iteration 47, residual = 1.9905656697202943e-09 [2KCG iteration 48, residual = 1.3307392471866413e-09 [2KCG iteration 49, residual = 9.182664330218499e-10 [2KCG iteration 50, residual = 6.413685294769359e-10 [2KCG iteration 51, residual = 4.2320693272587554e-10 [2KCG iteration 52, residual = 2.982873067494714e-10
CG iteration 53, residual = 1.9967848129516596e-10
CG iteration 54, residual = 1.331274500597547e-10 CG iteration 55, residual = 9.828100664655117e-11 CG iteration 56, residual = 6.84191762470924e-11 CG iteration 57, residual = 5.873539291511925e-11 CG iteration 58, residual = 3.695489118817188e-11 CG iteration 59, residual = 3.892786485799725e-11 CG iteration 60, residual = 4.016901801956522e-11 CG iteration 61, residual = 2.2423777959621626e-11 CG iteration 62, residual = 1.336705637699512e-11 CG iteration 63, residual = 9.67059864790468e-12 CG iteration 64, residual = 6.2120778447300835e-12 CG iteration 65, residual = 4.138466154392459e-12 CG iteration 66, residual = 2.7898181330352406e-12 CG iteration 67, residual = 2.0273009055889288e-12 CG iteration 68, residual = 1.2820786422091424e-12 CG iteration 69, residual = 9.032615755382615e-13 CG iteration 70, residual = 6.032276948557004e-13 CG iteration 71, residual = 3.954367512322008e-13 CG iteration 72, residual = 3.23881065971075e-13 CG iteration 73, residual = 4.1624078595039915e-13 CG iteration 74, residual = 2.1565106912329267e-13 CG iteration 75, residual = 1.7923895176345238e-13 CG iteration 76, residual = 1.1345744799045737e-13 CG iteration 77, residual = 7.2777689317211e-14 CG iteration 78, residual = 4.639753321863095e-14 CG iteration 79, residual = 3.482767062908518e-14 CG iteration 80, residual = 2.2254815344036962e-14 CG iteration 81, residual = 1.5389646784977875e-14 CG iteration 82, residual = 1.0607980872501223e-14 CG iteration 83, residual = 7.158278070991633e-15 CG iteration 84, residual = 5.516885156970535e-15 CG iteration 85, residual = 4.916092936040603e-15 CG iteration 86, residual = 3.4917874699350262e-15 CG converged in 86 iterations to residual 3.4917874699350262e-15
</pre>
CG iteration 54, residual = 1.331274500597547e-10 CG iteration 55, residual = 9.828100664655117e-11 CG iteration 56, residual = 6.84191762470924e-11 CG iteration 57, residual = 5.873539291511925e-11 CG iteration 58, residual = 3.695489118817188e-11 CG iteration 59, residual = 3.892786485799725e-11 CG iteration 60, residual = 4.016901801956522e-11 CG iteration 61, residual = 2.2423777959621626e-11 CG iteration 62, residual = 1.336705637699512e-11 CG iteration 63, residual = 9.67059864790468e-12 CG iteration 64, residual = 6.2120778447300835e-12 CG iteration 65, residual = 4.138466154392459e-12 CG iteration 66, residual = 2.7898181330352406e-12 CG iteration 67, residual = 2.0273009055889288e-12 CG iteration 68, residual = 1.2820786422091424e-12 CG iteration 69, residual = 9.032615755382615e-13 CG iteration 70, residual = 6.032276948557004e-13 CG iteration 71, residual = 3.954367512322008e-13 CG iteration 72, residual = 3.23881065971075e-13 CG iteration 73, residual = 4.1624078595039915e-13 CG iteration 74, residual = 2.1565106912329267e-13 CG iteration 75, residual = 1.7923895176345238e-13 CG iteration 76, residual = 1.1345744799045737e-13 CG iteration 77, residual = 7.2777689317211e-14 CG iteration 78, residual = 4.639753321863095e-14 CG iteration 79, residual = 3.482767062908518e-14 CG iteration 80, residual = 2.2254815344036962e-14 CG iteration 81, residual = 1.5389646784977875e-14 CG iteration 82, residual = 1.0607980872501223e-14 CG iteration 83, residual = 7.158278070991633e-15 CG iteration 84, residual = 5.516885156970535e-15 CG iteration 85, residual = 4.916092936040603e-15 CG iteration 86, residual = 3.4917874699350262e-15 CG converged in 86 iterations to residual 3.4917874699350262e-15
end{sphinxVerbatim}
[2KCG iteration 54, residual = 1.331274500597547e-10 [2KCG iteration 55, residual = 9.828100664655117e-11 [2KCG iteration 56, residual = 6.84191762470924e-11 [2KCG iteration 57, residual = 5.873539291511925e-11 [2KCG iteration 58, residual = 3.695489118817188e-11 [2KCG iteration 59, residual = 3.892786485799725e-11 [2KCG iteration 60, residual = 4.016901801956522e-11 [2KCG iteration 61, residual = 2.2423777959621626e-11 [2KCG iteration 62, residual = 1.336705637699512e-11 [2KCG iteration 63, residual = 9.67059864790468e-12 [2KCG iteration 64, residual = 6.2120778447300835e-12 [2KCG iteration 65, residual = 4.138466154392459e-12 [2KCG iteration 66, residual = 2.7898181330352406e-12 [2KCG iteration 67, residual = 2.0273009055889288e-12 [2KCG iteration 68, residual = 1.2820786422091424e-12 [2KCG iteration 69, residual = 9.032615755382615e-13 [2KCG iteration 70, residual = 6.032276948557004e-13 [2KCG iteration 71, residual = 3.954367512322008e-13 [2KCG iteration 72, residual = 3.23881065971075e-13 [2KCG iteration 73, residual = 4.1624078595039915e-13 [2KCG iteration 74, residual = 2.1565106912329267e-13 [2KCG iteration 75, residual = 1.7923895176345238e-13 [2KCG iteration 76, residual = 1.1345744799045737e-13 [2KCG iteration 77, residual = 7.2777689317211e-14 [2KCG iteration 78, residual = 4.639753321863095e-14 [2KCG iteration 79, residual = 3.482767062908518e-14 [2KCG iteration 80, residual = 2.2254815344036962e-14 [2KCG iteration 81, residual = 1.5389646784977875e-14 [2KCG iteration 82, residual = 1.0607980872501223e-14 [2KCG iteration 83, residual = 7.158278070991633e-15 [2KCG iteration 84, residual = 5.516885156970535e-15 [2KCG iteration 85, residual = 4.916092936040603e-15 [2KCG iteration 86, residual = 3.4917874699350262e-15 [2KCG converged in 86 iterations to residual 3.4917874699350262e-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);
[ ]:
[ ]: