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.004809678530125145
CG iteration 2, residual = 0.003322319013244243
CG iteration 3, residual = 0.003311588317754009
CG iteration 4, residual = 0.002746757977924576
CG iteration 5, residual = 0.0014658765270595317
CG iteration 6, residual = 0.001216702499924429
CG iteration 7, residual = 0.000809658230957325
CG iteration 8, residual = 0.0006570076838867991
CG iteration 9, residual = 0.0004750869183442853
CG iteration 10, residual = 0.00036224214074475833
CG iteration 11, residual = 0.00025441834248533064
CG iteration 12, residual = 0.0001619498909443687
CG iteration 13, residual = 0.00011358302702051922
CG iteration 14, residual = 8.950382818006825e-05
CG iteration 15, residual = 5.394126812247181e-05
CG iteration 16, residual = 3.916041785654374e-05
CG iteration 17, residual = 2.708231725754564e-05
CG iteration 18, residual = 1.824959941708822e-05
CG iteration 19, residual = 1.3585594383092305e-05
CG iteration 20, residual = 9.897420169931647e-06
CG iteration 21, residual = 1.1862519518942149e-05
CG iteration 22, residual = 5.672249675176412e-06
CG iteration 23, residual = 3.5937742573560138e-06
CG iteration 24, residual = 2.5589080960998443e-06
CG iteration 25, residual = 1.977947786957905e-06
CG iteration 26, residual = 1.2770934679739879e-06
CG iteration 27, residual = 8.740832403651163e-07
CG iteration 28, residual = 5.927228746857361e-07
CG iteration 29, residual = 4.0089017034467604e-07
CG iteration 30, residual = 2.856458196958184e-07
CG iteration 31, residual = 2.0292685846589974e-07
CG iteration 32, residual = 1.3199636987195287e-07
CG iteration 33, residual = 9.850315540116874e-08
CG iteration 34, residual = 6.734297427784079e-08
CG iteration 35, residual = 4.756053902801611e-08
CG iteration 36, residual = 2.9198312697085005e-08
CG iteration 37, residual = 2.6868702154784086e-08
CG iteration 38, residual = 2.3120432849520918e-08
CG iteration 39, residual = 1.3026197474420954e-08
CG iteration 40, residual = 8.395492410361875e-09
CG iteration 41, residual = 5.711979109343909e-09
CG iteration 42, residual = 4.1810498425022584e-09
CG iteration 43, residual = 2.794949339461186e-09
CG iteration 44, residual = 1.8595331195625792e-09
CG iteration 45, residual = 1.265449519583777e-09
CG iteration 46, residual = 8.583836344844582e-10
CG iteration 47, residual = 5.766587131640494e-10
CG iteration 48, residual = 4.2177694151955625e-10
CG iteration 49, residual = 2.7611851466713986e-10
CG iteration 50, residual = 2.1705250681804394e-10
CG iteration 51, residual = 1.285950970961264e-10
CG iteration 52, residual = 8.655706785450058e-11
CG iteration 53, residual = 7.713457336664312e-11
CG iteration 54, residual = 6.839789731206933e-11
CG iteration 55, residual = 3.477938978160038e-11
CG iteration 56, residual = 2.3705852575695383e-11
CG iteration 57, residual = 1.9090848972834223e-11
CG iteration 58, residual = 1.335519527548304e-11
CG iteration 59, residual = 8.553536510474682e-12
CG iteration 60, residual = 5.929207225773354e-12
CG iteration 61, residual = 3.943106905920753e-12
CG iteration 62, residual = 2.695663462934199e-12
CG iteration 63, residual = 1.7387105754137911e-12
CG iteration 64, residual = 1.2104139306427563e-12
CG iteration 65, residual = 7.623361233406957e-13
CG iteration 66, residual = 4.941739079787258e-13
CG iteration 67, residual = 3.208332319416464e-13
CG iteration 68, residual = 2.226294122189969e-13
CG iteration 69, residual = 3.0690592832133193e-13
CG iteration 70, residual = 1.3682804314400704e-13
CG iteration 71, residual = 9.423558910368108e-14
CG iteration 72, residual = 6.131182773847482e-14
CG iteration 73, residual = 4.207002668758155e-14
CG iteration 74, residual = 2.7629436320843212e-14
CG iteration 75, residual = 1.7674924165778345e-14
CG iteration 76, residual = 1.1723300787131031e-14
CG iteration 77, residual = 8.797610857770427e-15
CG iteration 78, residual = 6.7658108708319185e-15
CG iteration 79, residual = 4.039532199219642e-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);
[ ]:

[ ]: