Reference¶
inductance¶
Inductance.
Coils¶
Coil inductance calculations.
Defines a coil class to keep track of coil parameters.
Benchmarking against LDX values, which come from old Mathematica routines and other tests.
- Filaments are defined as an numpy 3 element vector
r, z, and n.
- class inductance.coils.Coil(r, z, dr, dz, nt=1, at=None, nr=0, nz=0, theta=0, **kwargs)¶
Rectangular coil object to keep track of coil parameters.
- Fr_filament(C2)¶
Radial force of two coils by filamentation.
- Parameters:
C2 (Coil)
- Return type:
float
- Fr_self()¶
Radial force of coil on itself.
- Return type:
float
- Fz_filament(C2)¶
Vertical force of two coils by filamentation.
- Parameters:
C2 (Coil)
- Return type:
float
- L_Lyle4()¶
Inductance by Lyle’s formula, 4th order.
- L_Lyle6()¶
Inductance by Lyle’s formula, 6th order.
- L_Lyle6A()¶
Inductance by Lyle’s formula, 6th order, appendix.
- L_Maxwell()¶
Inductance by Maxwell’s formula.
- L_best()¶
Inductance by best formula.
- L_filament(nr=0, nz=0)¶
Inductance by filamentation.
- L_long_solenoid_butterworth()¶
Inductance by Butterworth’s formula.
- L_lorentz()¶
Inductance by Lorentz’s formula.
- M_filament(C2)¶
Mutual inductance of two coils by filamentation.
- Parameters:
C2 (Coil)
- Return type:
float
- dLdR_Lyle6()¶
Derivative of inductance by Lyle’s formula, 6th order.
- filamentize(nr, nz)¶
Create an array of filaments to represent the coil.
- classmethod from_bounds(r1, r2, z1, z2, nt=1, at=1, nr=0, nz=0)¶
Create a coil from bounds instead of center and width.
- classmethod from_dict(d)¶
Create a coil from a dictionary.
- class inductance.coils.CompositeCoil(coils, **kwargs)¶
A coil made of multiple rectangular coils.
- Parameters:
coils (list[Coil])
- L_best()¶
Inductance of composite coils by best formula.
- class inductance.coils.Conductor(shape=Shape.round, r=0.0, dr=0.0, dz=0.0, r_i=0.0, rho=0.0)¶
Conductor object to keep track of conductor parameters.
- Parameters:
shape (Shape | str)
r (float)
dr (float)
dz (float)
r_i (float)
rho (float)
- class inductance.coils.Shape(value)¶
Enum for conductor shapes.
- inductance.coils.coilset_Fr_greens(coils)¶
Calculate the radial force matrix of a set of coils.
- Parameters:
coils (list[Coil]) – set of coils
- Returns:
2D array of radial forces per amp**2
- Return type:
np.ndarray
Self inductance¶
Self nductance calculations for coils.
author: Darren Garnier <garnier@mit.edu>
basic equations come from analytic approximations of old. tested in real life with the LDX Fcoil / Ccoil / Lcoil system.
One from Maxwell himself, and better ones from:
- Lyle, T. R. “On the Self-inductance of Circular Coils of
Rectangular Section”. Roy. Soc. London. A. V213 (1914) pp 421-435. https://doi.org/10.1098/rsta.1914.0009
Unfortunately, Lyle doesn’t work that well with large dz/R coils. Other approximations are also included.
This code now uses numba to do just-in-time compiliation and parallel execution to greatly increase speed. Also requires numba-scipy for elliptical functions. numba-scipy can be fragile and sometimes needs to be “updated” before installation to the newest version of numba.
- inductance.self.L_hollow_round(r, a, n)¶
Self inductance of a round conductor coil with skin current.
- Parameters:
r (float) – coil centerline radius
a (float) – coil conductor radius
n (int) – number of turns
- Returns:
coil self inductance in Henrys
- Return type:
float
- inductance.self.L_long_solenoid_butterworth(r, dr, dz, n)¶
Self inductance of a long solenoid by Butterworth’s formula.
As written in Grover, Bulletin of the Bureau of Standards, Vol. 14 pg. 558 https://nvlpubs.nist.gov/nistpubs/bulletin/14/nbsbulletinv14n4p537_A2b.pdf
Original S Butterworth 1914 Proc. Phys. Soc. London 27 371
Applies when dz > 2*r.
- Parameters:
r (float) – coil centerline radius
dr (float) – coil radial width
dz (float) – coil height
n (int) – number of turns
- Returns:
coil self inductance in Henrys
- Return type:
float
- inductance.self.L_lorentz(r, _dr, dz, n)¶
Self inductance of a thin wall solenoid by Lorentz’s formula.
Given in: Rosa and Grover, Formulas and Tables for the Calculation of Mutual and Self-Inductance https://nvlpubs.nist.gov/nistpubs/bulletin/08/nbsbulletinv8n1p1_A2b.pdf or https://www.jstor.org/stable/pdf/24521000.pdf (Formula 72 on page 118)
- Parameters:
r (float) – coil centerline radius
_dr (float) – coil radial width (ignored)
dz (float) – coil height
n (int) – number of turns
- Returns:
coil self inductance in Henrys
- Return type:
float
- inductance.self.L_lyle4(r, dr, dz, n)¶
Self inductance of a rectangular coil via Lyle to 4th order, Eq3.
- Parameters:
r (float) – coil centerline radius
dr (float) – coil radial width
dz (float) – coil height
n (int) – number of turns
- Returns:
coil self inductance in Henrys
- Return type:
float
- inductance.self.L_lyle4_eq4(r, dr, dz, n)¶
Self inductance of a rectangular coil via Lyle to 4th order, Eq4.
this doesn’t give quite the same answer as eq3 above. and it doesn’t seem to work as well.
- Parameters:
r (float) – coil centerline radius
dr (float) – coil radial width
dz (float) – coil height
n (int) – number of turns
- Returns:
coil self inductance in Henrys
- Return type:
float
- inductance.self.L_lyle6(r, dr, dz, n)¶
Self inductance of a rectangular coil via Lyle to 6th order.
- Parameters:
r (float) – coil centerline radius
dr (float) – coil radial width
dz (float) – coil height
n (int) – number of turns
- Returns:
coil self inductance in Henrys
- Return type:
float
- inductance.self.L_lyle6_appendix(r, dr, dz, n)¶
Self inductance of a rectangular coil via Lyle to 6th order, appendix.
- Parameters:
r (float) – coil centerline radius
dr (float) – coil radial width
dz (float) – coil height
n (int) – number of turns
- Returns:
coil self inductance in Henrys
- Return type:
float
- inductance.self.L_lyle_sectioning(r, dr, dz, nt, nr, nz)¶
Self inductance by sectioning, Lyle’s method, and Lyle’s 4th order.
- Parameters:
r (float) – coil centerline radius
dr (float) – coil radial width
dz (float) – coil height
nt (float) – number of turns
nr (int) – number of radial sections
nz (int) – number of axial sections
- Returns:
coil self inductance in Henrys
- Return type:
float
- inductance.self.L_maxwell(r, dr, dz, n)¶
Self inductance of a rectangular coil with constant current density by Maxwell.
- Parameters:
r (float) – coil centerline radius
dr (float) – coil radial width
dz (float) – coil height
n (int) – number of turns
- Returns:
coil self inductance in Henrys
- Return type:
float
- inductance.self.L_round(r, a, n)¶
Self inductance of a round conductor coil with constant current density.
- Parameters:
r (float) – coil centerline radius
a (float) – coil conductor radius
n (int) – number of turns
- Returns:
coil self inductance in Henrys
- Return type:
float
- inductance.self.dLdR_lyle6(r, dr, dz, n)¶
Radial derivative of self inductance of a rectangular coil via Lyle to 6th order.
- Parameters:
r (float) – coil centerline radius
dr (float) – coil radial width
dz (float) – coil height
n (int) – number of turns
- Returns:
radial derivative of inductance in Henrys/meter
- Return type:
float
- inductance.self.self_inductance_by_filaments(f, conductor='round', a=0.01, dr=0.01, dz=0.01)¶
Self inductance of filament set.
- Parameters:
f (array) – first filament array
conductor (str, optional) – conductor shape. Defaults to “round”.
a (float, optional) – conductor radius. Defaults to 0.01.
dr (float, optional) – conductor radial width. Defaults to 0.01
dz (float, optional) – conductor vertical height. Defaults to 0.01
- Returns:
self inductance of filament set in Henries
- Return type:
float
Mutual inductance¶
Mutual inductance calculations for coils.
author: Darren Garnier <garnier@mit.edu>
- inductance.mutual.lyle_equivalent_filaments(r, z, dr, dz, nt, fils)¶
Compute the equivalent filament locations for Lyle’s method.
Using Lyle’s Method of Equivalent Filaments.
- originally from:
Lyle, Phil. Mag., 3, p. 310; 1902.
- reproduced in:
Rosa and Grover, “Formulas for the Mutual Inductance of Coaxial Circular Coils of Rectangular Section,” Bull. Natl. Bur. Stand., vol 8, no. 1, p. 38-39; 1911.
- Parameters:
r (float) – inner radius of coil
z (float) – inner height of coil
dr (float) – radial width of coil
dz (float) – height of coil
nt (float) – number of turns in coil
fils (numpy.ndarray) – array of filaments 2 x (r, z, n)
- inductance.mutual.lyle_equivalent_subcoil_filaments(subcoils)¶
Compute the equivalent filament locations for set of subcoils.
- inductance.mutual.mutual_lyles_method(r1, z1, dr1, dz1, nt1, r2, z2, dr2, dz2, nt2)¶
Mutual inductance of two coils by Lyle’s method.
Using Lyle’s Method of Equivalent Filaments.
originally from: Lyle, Phil. Mag., 3, p. 310; 1902.
reproduced: Rosa and Grover, “Formulas for the Mutual Inductance of Coaxial Circular Coils of Rectangular Section,” Bull. Natl. Bur. Stand., vol 8, no. 1, p. 38-39; 1911.
- Parameters:
r1 (float) – inner radius of coil 1
z1 (float) – inner height of coil 1
dr1 (float) – radial width of coil 1
dz1 (float) – height of coil 1
nt1 (int) – number of turns in coil 1
r2 (float) – inner radius of coil 2
z2 (float) – inner height of coil 2
dr2 (float) – radial width of coil 2
dz2 (float) – height of coil 2
nt2 (int) – number of turns in coil 2
- Returns:
mutual inductance of the two coils
- Return type:
float
- inductance.mutual.mutual_rayleigh(r1, z1, dr1, dz1, n1, r2, z2, dr2, dz2, n2)¶
Mutual inductance of two coils by Rayleigh’s Quadrature Method.
reproduced in : Rosa and Grover, “Formulas for the Mutual Inductance of Coaxial Circular Coils of Rectangular Section,” Bull. Natl. Bur. Stand., vol 8, no. 1, p. 34-35; 1911.
- Parameters:
r1 (float) – inner radius of coil 1
z1 (float) – inner height of coil 1
dr1 (float) – radial width of coil 1
dz1 (float) – height of coil 1
n1 (int) – number of turns in coil 1
r2 (float) – inner radius of coil 2
z2 (float) – inner height of coil 2
dr2 (float) – radial width of coil 2
dz2 (float) – height of coil 2
n2 (int) – number of turns in coil 2
- Returns:
mutual inductance of the two coils
- Return type:
float
- inductance.mutual.mutual_sectioning_lyle(r1, z1, dr1, dz1, nt1, nr1, nz1, r2, z2, dr2, dz2, nt2, nr2, nz2)¶
Mutual inductance by sectioning of two coils by Lyle’s Method.
Section cois into subcoils and compute mutual inductance of each set of subcoil using Lyle’s method of equivalent filaments.
- Parameters:
r1 (float) – inner radius of coil 1
z1 (float) – inner height of coil 1
dr1 (float) – radial width of coil 1
dz1 (float) – height of coil 1
nt1 (float) – number of turns in coil 1
nr1 (int) – number of radial sections in coil 1
nz1 (int) – number of vertical sections in coil 1
r2 (float) – inner radius of coil 2
z2 (float) – inner height of coil 2
dr2 (float) – radial width of coil 2
dz2 (float) – height of coil 2
nt2 (float) – number of turns in coil 2
nr2 (int) – number of radial sections in coil 2
nz2 (int) – number of vertical sections in coil 2
- Returns:
mutual inductance of the two coils
- Return type:
float
elliptics¶
Precise and fast elliptic integrals.
Author: Fukushima, T. <Toshio.Fukushima@nao.ac.jp>
Reference: Toshio Fukushima, “Precise and fast computation of a general incomplete elliptic integral of second kind by half and double argument transformations” Journal of Computational and Applied Mathematics 235 (2011) 4140–4148 DOI: https://doi.org/10.1090/S0025-5718-2011-02455-5
With Numba, this runs ~10x faster than scipy.special if you need both k and e so.. don’t need numba_scipy which lags numba and scipy often.
- inductance.elliptics.celbd(mc)¶
Complete elliptic integrals of second kind, B(m) and D(m).
- Parameters:
mc (float) – mc = 1-k^2, where k is the elliptic modulus.
- Returns:
(B(m), D(m)), the associate complete elliptic integrals of second kind
- inductance.elliptics.ellipe(k2)¶
Complete elliptic integral of the second kind.
- Parameters:
k2 (float) – K^2 the square of the elliptic modulus.
- Returns:
Elliptic integral E(k^2).
- Return type:
(float)
- inductance.elliptics.ellipk(k2)¶
Complete elliptic integral of the first kind.
- Parameters:
k2 (float) – K^2 the square of the elliptic modulus.
- Returns:
Elliptic integral K(k^2).
- Return type:
(float)
- inductance.elliptics.ellipke(k2)¶
Complete elliptic integrals of first and second kind.
- Parameters:
k2 (float) – K^2 the square of the elliptic modulus.
- Returns:
Elliptic integrals K(k^2) and E(k^2).
- Return type:
(float, float)
filaments¶
Filamentary inductance calculations.
- inductance.filaments.AGreen(r, z, a, b)¶
Psi at position r, z, due to a filament with radius a, and z postion b.
- Parameters:
r (float) – radial position of a point
z (float) – vertical position of a point
a (float) – radius of a filament
b (float) – vertical position of a filament
- Returns:
Psi at r, z in Weber / Amp
- Return type:
float
- inductance.filaments.ASegment(pts, xyz, uvw)¶
Psi at positions pts, due to a linear current segment uvw located at xyz.
Derived from https://doi.org/10.1016/j.cpc.2023.108692, Eq. 20
- Parameters:
pts (array) – coordinates to calculate vector potential at, N x (x,y,z)
xyz (array) – coordinates of base of current element (x,y,z)
uvw (array) – displacement vector for length of element
- Returns:
Psi at pts in Weber / Amp, N x (Psix,Psiy,Psiz)
- Return type:
array
- inductance.filaments.BrGreen(r, z, a, b)¶
Br at position r, z, due to a filament with radius a, and z postion b.
- Parameters:
r (float) – radial position of a point
z (float) – vertical position of a point
a (float) – radius of a filament
b (float) – vertical position of a filament
- Returns:
Br at r, z in Tesla / Amp
- Return type:
float
- inductance.filaments.BzGreen(r, z, a, b)¶
Bz at position r, z, due to a filament with radius a, and z postion b.
- Parameters:
r (float) – radial position of a point
z (float) – vertical position of a point
a (float) – radius of a filament
b (float) – vertical position of a filament
- Returns:
Bz at r, z in Tesla / Amp
- Return type:
float
- inductance.filaments.L_approx_path_rect(pts, w, h, n, ds=1)¶
Approximate self inductance of a path of points with a rectangular cross section.
take a path of points n x 3, with a cross section b x c and approximate self inductance using Dengler
- inductance.filaments.M_filsol_path(fils, pts, nt, ds=0)¶
Mutual inductance between a set of axisymmetric filaments and a path from pts.
- inductance.filaments.M_path_path(pts1, pts2, ds1=0, ds2=0)¶
Mutual inductance between two pts paths.
- inductance.filaments.filament_coil(r, z, dr, dz, nt, nr, nz, theta=0)¶
Create an array of filaments, each with its own radius, height, and amperage.
r : Major radius of coil center. z : Vertical center of coil. dr : Radial width of coil. dz : Height of coil. nt : number of turns in coil nr : Number of radial slices nz : Number of vertical slices theta : Rotation of coil about phi axis
Returns: Array of shape (nr*nz) x 3 of R, Z, N for each filament
- inductance.filaments.green_sum_over_filaments(gfunc, fil, r, z)¶
Calculate a greens function at position r, z, to an array of filaments.
- Parameters:
gfunc (function) – Green’s function to use
fil (array) – filament array N x (r, z, n)
r (float) – radial position of a point
z (float) – vertical position of a point
- Returns:
sum of the greens function at r, z, due to all filaments
- Return type:
float
- inductance.filaments.mutual_filaments_segmented(fils, pts)¶
Mutual inductance between a set of axisymmetric filaments and a segmented path.
- inductance.filaments.mutual_inductance_fil(rzn1, rzn2)¶
Mutual inductance of two filaments.
- Parameters:
rzn1 (array) – (r, z, n) of first filament
rzn2 (array) – (r, z, n) of second filament
- Returns:
mutual inductance in Henrys
- Return type:
float
- inductance.filaments.mutual_inductance_of_filaments(f1, f2)¶
Mutual inductance of two sets of filaments.
These should not be the same filament array, not setup to handle self inductance of filament.
- Parameters:
f1 (array) – first filament array
f2 (array) – second filament array
- Returns:
Mutual inductance of filament sets in Henries
- Return type:
float
- inductance.filaments.radial_force_fil(rzn1, rzn2)¶
Radial force between two filaments per conductor amp.
- Parameters:
rzn1 (array) – (r, z, n) of first filament
rzn2 (array) – (r, z, n) of second filament
- Returns:
force in Newtons/Amp^2
- Return type:
float
- inductance.filaments.radial_force_of_filaments(f1, f2)¶
Radial force between two sets of filaments.
These should not be the same filament array, not setup to handle self inductance of filament.
- Parameters:
f1 (array) – first filament array
f2 (array) – second filament array
- Returns:
Radial force in Newtons/Amp^2
- Return type:
float
- inductance.filaments.segment_path(pts, ds=0, close=False)¶
Segment a path into equal length segments.
- Parameters:
pts (array) – N x 3 array of points
ds (float, optional) – length between points. Defaults to minimun in pts.
close (bool, optional) – Should the path close the points. Defaults to False.
- Returns:
M x 3 array of points along the path array: M array of length along the path
- Return type:
array
- inductance.filaments.segmented_self_inductance(pts, s, a)¶
Self inductance of a segmented path by double integral.
- Parameters:
pts (N x 3 array) – points along the path
s (array) – length along the path
a (float) – radius of wire
- Returns:
self inductance of the path
- Return type:
float
Neumann’s formula.. double curve integral of mu_0/(4*pi) * int * int ( dx1 . dx2 / norm(x1-x2))) do all points except where x1-x2 blows up. instead follow Dengler https://doi.org/10.7716/aem.v5i1.331 which makes a approximation for that is good O(mu_0*a) lets just assume the thing is broken into small pieces and neglect end points
this code doesn’t work very well.. comparison test sorta fails phi = np.linspace(0,2*math.pi,100) test_xyz = np.array([[np.cos(p), np.sin(p), 0] for p in phi]) L_seg = L_approx_path_rect(test_xyz, .01, .01, 1, .1) L_maxwell(1, .01, .01, 1), L_lyle6(1, .01, .01, 1), L_seg (6.898558527897293e-06, 6.8985981243869525e-06, 6.907313505254537e-06) # Y=1/4 hmmm… (6.898558527897293e-06, 6.8985981243869525e-06, 7.064366776069971e-06) # Y=1/2
- inductance.filaments.sum_over_filaments(func, f1, f2)¶
Apply a function and sum over all combinations of two sets of filaments.
- Parameters:
func (function) – function to apply to each pair of filaments
f1 (array) – first filament array
f2 (array) – second filament array
- Returns:
result of summing over all combinations of func(f1[i], f2[j])
- Return type:
float
- inductance.filaments.vertical_force_fil(rzn1, rzn2)¶
Vertical force between two filaments per conductor amp.
- Parameters:
rzn1 (array) – (r, z, n) of first filament
rzn2 (array) – (r, z, n) of second filament
- Returns:
force in Newtons/Amp^2
- Return type:
float
- inductance.filaments.vertical_force_of_filaments(f1, f2)¶
Vertical force between two sets of filaments.
These should not be the same filament array, not setup to handle self inductance of filament.
- Parameters:
f1 (array) – first filament array
f2 (array) – second filament array
- Returns:
Vertical force in Newtons/Amp^2
- Return type:
float