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 = 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.004063969731501542
CG iteration 3, residual = 0.003290581389748055
CG iteration 4, residual = 0.0022607529739440226
CG iteration 5, residual = 0.0019757898056008837
CG iteration 6, residual = 0.0011507895291099671
CG iteration 7, residual = 0.0008639313194094752
CG iteration 8, residual = 0.0006327449629420533
CG iteration 9, residual = 0.0004949733354845918
CG iteration 10, residual = 0.000407605412414315
CG iteration 11, residual = 0.00031648668265075867
CG iteration 12, residual = 0.0001839395580600599
CG iteration 13, residual = 0.00013713589121295468
CG iteration 14, residual = 8.753130614940792e-05
CG iteration 15, residual = 6.40469388139117e-05
CG iteration 16, residual = 4.235156707865746e-05
CG iteration 17, residual = 2.6160940650642946e-05
CG iteration 18, residual = 1.9546762209334134e-05
CG iteration 19, residual = 1.3580852760593627e-05
CG iteration 20, residual = 1.798808636822381e-05
CG iteration 21, residual = 8.038807718846506e-06
CG iteration 22, residual = 5.687476695979545e-06
CG iteration 23, residual = 3.828297235743515e-06
CG iteration 24, residual = 2.512607983448245e-06
CG iteration 25, residual = 1.9133034156165156e-06
CG iteration 26, residual = 1.1882807373093324e-06
CG iteration 27, residual = 8.414430311849796e-07
CG iteration 28, residual = 5.733999040421527e-07
CG iteration 29, residual = 3.8563237039044054e-07
CG iteration 30, residual = 2.812965336962804e-07
CG iteration 31, residual = 1.7914275627863905e-07
CG iteration 32, residual = 1.3497785654935678e-07
CG iteration 33, residual = 1.1912809676097392e-07
CG iteration 34, residual = 7.74511350839512e-08
CG iteration 35, residual = 4.9373756139113466e-08
CG iteration 36, residual = 4.2275704064855194e-08
CG iteration 37, residual = 4.3965069064242005e-08
CG iteration 38, residual = 2.1808705528094614e-08
CG iteration 39, residual = 1.446217392478935e-08
CG iteration 40, residual = 9.528478377092738e-09
CG iteration 41, residual = 7.500741406901941e-09
CG iteration 42, residual = 6.6211804577309514e-09
CG iteration 43, residual = 3.649723574950662e-09
CG iteration 44, residual = 2.623206395249842e-09
CG iteration 45, residual = 1.7228864151696825e-09
CG iteration 46, residual = 1.1298898931680993e-09
CG iteration 47, residual = 7.954171308176308e-10
CG iteration 48, residual = 5.005044008281257e-10
CG iteration 49, residual = 3.505161311978582e-10
CG iteration 50, residual = 2.270488008941397e-10
CG iteration 51, residual = 1.5440626800501445e-10
CG iteration 52, residual = 1.2468174710905142e-10
CG iteration 53, residual = 1.2874022018142688e-10
CG iteration 54, residual = 6.068034109841038e-11
CG iteration 55, residual = 4.031804501167675e-11
CG iteration 56, residual = 2.7648912176648625e-11
CG iteration 57, residual = 1.7867635716386904e-11
CG iteration 58, residual = 1.2085567765457076e-11
CG iteration 59, residual = 7.87406428187445e-12
CG iteration 60, residual = 5.942759817867974e-12
CG iteration 61, residual = 5.4238388898342324e-12
CG iteration 62, residual = 3.05769905687807e-12
CG iteration 63, residual = 1.977053541247477e-12
CG iteration 64, residual = 1.2730049814720828e-12
CG iteration 65, residual = 8.765308201773413e-13
CG iteration 66, residual = 6.030245656808731e-13
CG iteration 67, residual = 4.1057758034037135e-13
CG iteration 68, residual = 5.168604312627611e-13
CG iteration 69, residual = 3.2580334757153545e-13
CG iteration 70, residual = 2.3846450808262284e-13
CG iteration 71, residual = 1.3877645655520924e-13
CG iteration 72, residual = 8.916468381434764e-14
CG iteration 73, residual = 5.86318571168264e-14
CG iteration 74, residual = 4.024682504021556e-14
CG iteration 75, residual = 2.6754017156928105e-14
CG iteration 76, residual = 1.7722763325835078e-14
CG iteration 77, residual = 1.1619601653904837e-14
CG iteration 78, residual = 7.60246316338744e-15
CG iteration 79, residual = 4.828790143444771e-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);
[ ]:

[ ]: