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\):

\[\DeclareMathOperator{\Grad}{grad} \DeclareMathOperator{\Curl}{curl} \DeclareMathOperator{\Div}{div} B = \mu (H + M), \quad \Div B = 0, \quad \Curl H = 0\]

Introducing a vector-potential \(A\) such that \(B = \Curl A\), and putting equations together we get

\[\Curl \mu^{-1} \Curl A = \Curl M\]

In weak form: Find \(A \in H(\Curl)\) such that

\[\int \mu^{-1} \Curl A \Curl v = \int M \Curl v \qquad \forall \, v \in H(\Curl)\]

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.004809678530125143
CG iteration 2, residual = 0.0033223190132442505
CG iteration 3, residual = 0.003311588317754026
CG iteration 4, residual = 0.002746757977924581
CG iteration 5, residual = 0.0014658765270595348
CG iteration 6, residual = 0.0012167024999244318
CG iteration 7, residual = 0.0008096582309573281
CG iteration 8, residual = 0.0006570076838868013
CG iteration 9, residual = 0.0004750869183442874
CG iteration 10, residual = 0.0003622421407447595
CG iteration 11, residual = 0.00025441834248533156
CG iteration 12, residual = 0.00016194989094436927
CG iteration 13, residual = 0.00011358302702052072
CG iteration 14, residual = 8.950382818008865e-05
CG iteration 15, residual = 5.394126812323255e-05
CG iteration 16, residual = 3.916041789777668e-05
CG iteration 17, residual = 2.7082318743444127e-05
CG iteration 18, residual = 1.82496688342686e-05
CG iteration 19, residual = 1.3588397852641488e-05
CG iteration 20, residual = 9.967377904483639e-06
CG iteration 21, residual = 1.209944533927576e-05
CG iteration 22, residual = 5.635002815842918e-06
CG iteration 23, residual = 3.593584639059988e-06
CG iteration 24, residual = 2.5589062786489067e-06
CG iteration 25, residual = 1.9779477532020695e-06
CG iteration 26, residual = 1.2770934675292937e-06
CG iteration 27, residual = 8.740832403609926e-07
CG iteration 28, residual = 5.927228746857051e-07
CG iteration 29, residual = 4.0089017034467885e-07
CG iteration 30, residual = 2.8564581969582147e-07
CG iteration 31, residual = 2.029268584659259e-07
CG iteration 32, residual = 1.319963698730612e-07
CG iteration 33, residual = 9.85031554539965e-08
CG iteration 34, residual = 6.73429761428992e-08
CG iteration 35, residual = 4.756062034938638e-08
CG iteration 36, residual = 2.9202841272564594e-08
CG iteration 37, residual = 2.708664597673188e-08
CG iteration 38, residual = 2.2998281343997115e-08
CG iteration 39, residual = 1.30233543774e-08
CG iteration 40, residual = 8.395472660932824e-09
CG iteration 41, residual = 5.711978975064796e-09
CG iteration 42, residual = 4.181049841074685e-09
CG iteration 43, residual = 2.7949493394749902e-09
CG iteration 44, residual = 1.8595331197202832e-09
CG iteration 45, residual = 1.265449520354632e-09
CG iteration 46, residual = 8.583836387660268e-10
CG iteration 47, residual = 5.766587352195382e-10
CG iteration 48, residual = 4.217770135678216e-10
CG iteration 49, residual = 2.7611877009180503e-10
CG iteration 50, residual = 2.17052979303827e-10
CG iteration 51, residual = 1.2859627633375243e-10
CG iteration 52, residual = 8.658889607066979e-11
CG iteration 53, residual = 7.812279362603221e-11
CG iteration 54, residual = 6.781572647690459e-11
CG iteration 55, residual = 3.4987340650179164e-11
CG iteration 56, residual = 2.5313273336168544e-11
CG iteration 57, residual = 2.0871879497639534e-11
CG iteration 58, residual = 1.288510900843055e-11
CG iteration 59, residual = 8.496910803410285e-12
CG iteration 60, residual = 5.923380318834461e-12
CG iteration 61, residual = 3.942613003304638e-12
CG iteration 62, residual = 2.6956168799157676e-12
CG iteration 63, residual = 1.7387078545620551e-12
CG iteration 64, residual = 1.2104137149230655e-12
CG iteration 65, residual = 7.623361116180402e-13
CG iteration 66, residual = 4.941739268307285e-13
CG iteration 67, residual = 3.208349162947707e-13
CG iteration 68, residual = 2.2276471618596654e-13
CG iteration 69, residual = 3.070773619246229e-13
CG iteration 70, residual = 1.3678304129205398e-13
CG iteration 71, residual = 9.423603167791008e-14
CG iteration 72, residual = 6.131468231581626e-14
CG iteration 73, residual = 4.2080266852672506e-14
CG iteration 74, residual = 2.7669811168073193e-14
CG iteration 75, residual = 1.7918281705070442e-14
CG iteration 76, residual = 1.3229339090091068e-14
CG iteration 77, residual = 1.0016917001315497e-14
CG iteration 78, residual = 6.476721179604968e-15
CG iteration 79, residual = 3.97145600652088e-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);
[ ]:

[ ]: