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.004809678530125144
CG iteration 2, residual = 0.0033223190132442527
CG iteration 3, residual = 0.0033115883177540282
CG iteration 4, residual = 0.002746757977924581
CG iteration 5, residual = 0.0014658765270595356
CG iteration 6, residual = 0.0012167024999244324
CG iteration 7, residual = 0.0008096582309573278
CG iteration 8, residual = 0.000657007683886801
CG iteration 9, residual = 0.00047508691834428674
CG iteration 10, residual = 0.0003622421407447592
CG iteration 11, residual = 0.00025441834248533124
CG iteration 12, residual = 0.0001619498909443692
CG iteration 13, residual = 0.00011358302702052124
CG iteration 14, residual = 8.950382818010456e-05
CG iteration 15, residual = 5.3941268123842773e-05
CG iteration 16, residual = 3.916041793086991e-05
CG iteration 17, residual = 2.708231993602935e-05
CG iteration 18, residual = 1.8249724548314423e-05
CG iteration 19, residual = 1.3590647072081406e-05
CG iteration 20, residual = 1.0022692097301634e-05
CG iteration 21, residual = 1.2225512971630618e-05
CG iteration 22, residual = 5.612959160271574e-06
CG iteration 23, residual = 3.593473846519095e-06
CG iteration 24, residual = 2.5589052168485445e-06
CG iteration 25, residual = 1.9779477334810707e-06
CG iteration 26, residual = 1.2770934672694703e-06
CG iteration 27, residual = 8.740832403585674e-07
CG iteration 28, residual = 5.927228746856776e-07
CG iteration 29, residual = 4.008901703446766e-07
CG iteration 30, residual = 2.856458196958214e-07
CG iteration 31, residual = 2.0292685846595422e-07
CG iteration 32, residual = 1.3199636987443592e-07
CG iteration 33, residual = 9.850315551961416e-08
CG iteration 34, residual = 6.734297845951856e-08
CG iteration 35, residual = 4.756072135969078e-08
CG iteration 36, residual = 2.9208464705298893e-08
CG iteration 37, residual = 2.735172333094666e-08
CG iteration 38, residual = 2.2855142780999904e-08
CG iteration 39, residual = 1.3020140654899339e-08
CG iteration 40, residual = 8.395450351029657e-09
CG iteration 41, residual = 5.711978823373224e-09
CG iteration 42, residual = 4.181049839445579e-09
CG iteration 43, residual = 2.7949493394486036e-09
CG iteration 44, residual = 1.8595331196745397e-09
CG iteration 45, residual = 1.2654495201319464e-09
CG iteration 46, residual = 8.583836375292805e-10
CG iteration 47, residual = 5.766587288490987e-10
CG iteration 48, residual = 4.217769927700998e-10
CG iteration 49, residual = 2.7611869696381963e-10
CG iteration 50, residual = 2.1705286641823273e-10
CG iteration 51, residual = 1.2859684625158194e-10
CG iteration 52, residual = 8.665805790118842e-11
CG iteration 53, residual = 8.087901529435022e-11
CG iteration 54, residual = 6.623790939834542e-11
CG iteration 55, residual = 3.490914949309879e-11
CG iteration 56, residual = 2.489161136014723e-11
CG iteration 57, residual = 2.0694360665551432e-11
CG iteration 58, residual = 1.2958043698978588e-11
CG iteration 59, residual = 8.50308214118398e-12
CG iteration 60, residual = 5.923985808282691e-12
CG iteration 61, residual = 3.942663985985598e-12
CG iteration 62, residual = 2.695621684456591e-12
CG iteration 63, residual = 1.7387081351597967e-12
CG iteration 64, residual = 1.2104137372337175e-12
CG iteration 65, residual = 7.623361192250387e-13
CG iteration 66, residual = 4.941745306716762e-13
CG iteration 67, residual = 3.208916301853916e-13
CG iteration 68, residual = 2.272745860288322e-13
CG iteration 69, residual = 3.0791582254563115e-13
CG iteration 70, residual = 1.3573476080037426e-13
CG iteration 71, residual = 9.42254452416223e-14
CG iteration 72, residual = 6.131244531037756e-14
CG iteration 73, residual = 4.2072515394935235e-14
CG iteration 74, residual = 2.763930132077061e-14
CG iteration 75, residual = 1.7735652734792958e-14
CG iteration 76, residual = 1.2174309928228326e-14
CG iteration 77, residual = 9.64634197142272e-15
CG iteration 78, residual = 6.643690963006274e-15
CG iteration 79, residual = 3.989844239757128e-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);
[ ]:

[ ]: