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:

B=μ(H+M),divB=0,curlH=0

Introducing a vector-potential A such that B=curlA, and putting equations together we get

curlμ1curlA=curlM

In weak form: Find AH(curl) such that

μ1curlAcurlv=McurlvvH(curl)

Usually, the permeability μ is given as μ=μrμ0, with μ0=4π107 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 L2 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 converged in 79 iterations to residual 3.671626039004356e-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);
[ ]:

[ ]: