This page was generated from wta/coil.ipynb.

Magnetostatics

[1]:
from netgen.occ import *
from ngsolve import *
from ngsolve.webgui import Draw
from netgen.webgui import Draw as DrawGeo
import math

model of the coil:

[2]:
cyl = Cylinder((0,0,0), Z, r=0.01, h=0.03).faces[0]
heli = Edge(Segment((0,0), (12*math.pi, 0.03)), cyl)
ps = heli.start
vs = heli.start_tangent
pe = heli.end
ve = heli.end_tangent

e1 = Segment((0,0,-0.03), (0,0,-0.01))
c1 = BezierCurve( [(0,0,-0.01), (0,0,0), ps-vs, ps])
e2 = Segment((0,0,0.04), (0,0,0.06))
c2 = BezierCurve( [pe, pe+ve, (0,0,0.03), (0,0,0.04)])
spiral = Wire([e1, c1, heli, c2, e2])
circ = Face(Wire([Circle((0,0,-0.03), Z, 0.001)]))
coil = Pipe(spiral, circ)

coil.faces.maxh=0.2
coil.faces.name="coilbnd"
coil.faces.Max(Z).name="in"
coil.faces.Min(Z).name="out"
coil.mat("coil")
crosssection = coil.faces.Max(Z).mass
[3]:
DrawGeo (coil);
[4]:
box = Box((-0.04,-0.04,-0.03), (0.04,0.04,0.06))
box.faces.name = "outer"
air = box-coil
air.mat("air");

mesh-generation of coil and air-box:

[5]:
geo = OCCGeometry(Glue([coil,air]))
with TaskManager():
    mesh = Mesh(geo.GenerateMesh(meshsize.coarse, maxh=0.01)).Curve(3)

Draw (mesh, clipping={"y":1, "z":0, "dist":0.012});

checking mesh data materials and boundaries:

[6]:
mesh.ne, mesh.nv, mesh.GetMaterials(), mesh.GetBoundaries()
[6]:
(143106,
 24380,
 ('coil', 'air'),
 ('out',
  'coilbnd',
  'coilbnd',
  'coilbnd',
  'coilbnd',
  'coilbnd',
  'in',
  'outer',
  'outer',
  'outer',
  'outer',
  'outer',
  'outer'))

Solve a potential problem to determine current density in wire:

on the domain \(\Omega_{\text{coil}}\): \begin{eqnarray*} j & = & \sigma \nabla \Phi \\ \operatorname{div} j & = & 0 \end{eqnarray*} port boundary conditions: \begin{eqnarray*} \Phi & = & 0 \qquad \qquad \text{on } \Gamma_{\text{out}}, \\ j_n & = & \frac{1}{|S|} \quad \qquad \text{on } \Gamma_{\text{in}}, \end{eqnarray*} and \(j_n=0\) else

[7]:
fespot = H1(mesh, order=3, definedon="coil", dirichlet="out")
phi,psi = fespot.TnT()
sigma = 58.7e6
bfa = BilinearForm(sigma*grad(phi)*grad(psi)*dx).Assemble()
inv = bfa.mat.Inverse(freedofs=fespot.FreeDofs(), inverse="sparsecholesky")
lff = LinearForm(1/crosssection*psi*ds("in")).Assemble()
gfphi = GridFunction(fespot)
gfphi.vec.data = inv * lff.vec
[8]:
Draw (gfphi, draw_vol=False, clipping={"y":1, "z":0, "dist":0.012});

Solve magnetostatic problem:

current source is current from potential equation:

\[\int \mu^{-1} \operatorname{curl} u \cdot \operatorname{curl} v \, dx = \int j \cdot v \, dx\]
[9]:
fes = HCurl(mesh, order=2, nograds=True)
print ("HCurl dofs:", fes.ndof)
u,v = fes.TnT()
mu = 4*math.pi*1e-7
a = BilinearForm(1/mu*curl(u)*curl(v)*dx+1e-6/mu*u*v*dx)
pre = Preconditioner(a, "bddc")
f = LinearForm(sigma*grad(gfphi)*v*dx("coil"))
with TaskManager():
    a.Assemble()
    f.Assemble()
HCurl dofs: 741967
[10]:
from ngsolve.krylovspace import CGSolver
inv = CGSolver(a.mat, pre, printrates=True)
gfu = GridFunction(fes)
with TaskManager():
    gfu.vec.data = inv * f.vec
CG iteration 1, residual = 23.456360961780817
CG iteration 2, residual = 0.07358662709979243
CG iteration 3, residual = 0.015569080694485034
CG iteration 4, residual = 0.008258992917359598
CG iteration 5, residual = 0.004964089157337705
CG iteration 6, residual = 0.003024989175809241
CG iteration 7, residual = 0.0018217235413854557
CG iteration 8, residual = 0.0010894049318672307
CG iteration 9, residual = 0.0007371555849158308
CG iteration 10, residual = 0.00048222228766432496
CG iteration 11, residual = 0.0003553057013706158
CG iteration 12, residual = 0.0002646023719653133
CG iteration 13, residual = 0.00021137354793235758
CG iteration 14, residual = 0.0001680015300483049
CG iteration 15, residual = 0.00013290851994139058
CG iteration 16, residual = 0.00010491021410418639
CG iteration 17, residual = 8.03292433777739e-05
CG iteration 18, residual = 5.7935937544777016e-05
CG iteration 19, residual = 3.935742813959517e-05
CG iteration 20, residual = 2.5519855951361304e-05
CG iteration 21, residual = 1.7874947330350458e-05
CG iteration 22, residual = 1.1600340279539776e-05
CG iteration 23, residual = 7.831456827343547e-06
CG iteration 24, residual = 5.378854297698119e-06
CG iteration 25, residual = 3.5546170164860835e-06
CG iteration 26, residual = 2.319185939321161e-06
CG iteration 27, residual = 1.5674686961720885e-06
CG iteration 28, residual = 1.0428391406918824e-06
CG iteration 29, residual = 7.545623038959704e-07
CG iteration 30, residual = 5.551543247735974e-07
CG iteration 31, residual = 4.134976705533413e-07
CG iteration 32, residual = 3.0206405854810665e-07
CG iteration 33, residual = 2.2601567027506402e-07
CG iteration 34, residual = 1.707756974524792e-07
CG iteration 35, residual = 1.2775060292097993e-07
CG iteration 36, residual = 9.564477987983552e-08
CG iteration 37, residual = 7.457576099549836e-08
CG iteration 38, residual = 5.669867109636814e-08
CG iteration 39, residual = 4.118957135081684e-08
CG iteration 40, residual = 2.987652391167933e-08
CG iteration 41, residual = 2.080806954340328e-08
CG iteration 42, residual = 1.4685704228419295e-08
CG iteration 43, residual = 1.0739121025018784e-08
CG iteration 44, residual = 8.014573564384242e-09
CG iteration 45, residual = 5.812393523791446e-09
CG iteration 46, residual = 4.3200960810082094e-09
CG iteration 47, residual = 3.082840170935855e-09
CG iteration 48, residual = 2.287398773465948e-09
CG iteration 49, residual = 1.7353997017735513e-09
CG iteration 50, residual = 1.3044757793080283e-09
CG iteration 51, residual = 9.425240315019854e-10
CG iteration 52, residual = 6.780798641156394e-10
CG iteration 53, residual = 4.964893114910025e-10
CG iteration 54, residual = 3.589564982250396e-10
CG iteration 55, residual = 2.7340249327339066e-10
CG iteration 56, residual = 2.010455449975614e-10
CG iteration 57, residual = 1.3982168599071236e-10
CG iteration 58, residual = 1.0518828315265392e-10
CG iteration 59, residual = 8.173748892187305e-11
CG iteration 60, residual = 5.861377653580455e-11
CG iteration 61, residual = 4.194578626153289e-11
CG iteration 62, residual = 3.048293839790532e-11
CG iteration 63, residual = 2.2147654704009247e-11
[11]:
Draw (curl(gfu), mesh, draw_surf=False, \
      min=0, max=3e-4, clipping = { "y":1, "z" : 0, "function":False}, vectors = { "grid_size":100});
[ ]: