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]:
(107150,
18459,
('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: 556236
[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 = 22.997114412048084
CG iteration 2, residual = 0.09922187848371194
CG iteration 3, residual = 0.020412708461064053
CG iteration 4, residual = 0.01238877483952836
CG iteration 5, residual = 0.006063383666771182
CG iteration 6, residual = 0.0034272345659258157
CG iteration 7, residual = 0.002205213486797272
CG iteration 8, residual = 0.0013707875616166994
CG iteration 9, residual = 0.0010107977218527951
CG iteration 10, residual = 0.0007888390081900122
CG iteration 11, residual = 0.0006245761949393266
CG iteration 12, residual = 0.0005214004287838638
CG iteration 13, residual = 0.0004258155122100788
CG iteration 14, residual = 0.0003004772010953787
CG iteration 15, residual = 0.00021631450078339409
CG iteration 16, residual = 0.0001524860601835929
CG iteration 17, residual = 9.951340092841463e-05
CG iteration 18, residual = 6.0371583769903635e-05
CG iteration 19, residual = 4.242970923830777e-05
CG iteration 20, residual = 3.094770618498057e-05
CG iteration 21, residual = 2.1406775662635987e-05
CG iteration 22, residual = 1.3572670784205e-05
CG iteration 23, residual = 8.757076239506368e-06
CG iteration 24, residual = 5.654575033255604e-06
CG iteration 25, residual = 3.821768713815917e-06
CG iteration 26, residual = 2.550143412798042e-06
CG iteration 27, residual = 1.8994478106639868e-06
CG iteration 28, residual = 1.4010126024324198e-06
CG iteration 29, residual = 1.0361268124921403e-06
CG iteration 30, residual = 7.540083494145983e-07
CG iteration 31, residual = 6.065134104776754e-07
CG iteration 32, residual = 4.766804328771368e-07
CG iteration 33, residual = 3.513674872404517e-07
CG iteration 34, residual = 2.581111846843345e-07
CG iteration 35, residual = 1.8980505663473784e-07
CG iteration 36, residual = 1.5907913186915514e-07
CG iteration 37, residual = 1.1260129606446747e-07
CG iteration 38, residual = 7.544915762277095e-08
CG iteration 39, residual = 4.955714687524776e-08
CG iteration 40, residual = 3.211509914845666e-08
CG iteration 41, residual = 2.0617706400639565e-08
CG iteration 42, residual = 1.4845268353175168e-08
CG iteration 43, residual = 1.0441271041735683e-08
CG iteration 44, residual = 6.899215895738054e-09
CG iteration 45, residual = 4.5151770783977536e-09
CG iteration 46, residual = 3.683532328913326e-09
CG iteration 47, residual = 2.962706563375356e-09
CG iteration 48, residual = 2.1657790657609903e-09
CG iteration 49, residual = 1.6799370788069685e-09
CG iteration 50, residual = 1.2089505758409931e-09
CG iteration 51, residual = 9.195325585388838e-10
CG iteration 52, residual = 6.898557221341236e-10
CG iteration 53, residual = 5.110602219041415e-10
CG iteration 54, residual = 3.5918114963535376e-10
CG iteration 55, residual = 2.753744648002812e-10
CG iteration 56, residual = 1.9872440598448202e-10
CG iteration 57, residual = 1.3302729712469486e-10
CG iteration 58, residual = 9.112874611764849e-11
CG iteration 59, residual = 7.389510909252902e-11
CG iteration 60, residual = 5.4812747804558675e-11
CG iteration 61, residual = 3.2163331128080684e-11
CG iteration 62, residual = 2.2868864296625466e-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});
[ ]: