Inductance¶
Closed-form partial-inductance formulas of F. W. Grover.
These are the same expressions that the binary evaluates in
grover_segment_self_inductance, mutual_inductance_4corner_grover,
and coupled_wire_self_inductance_grover (decompiled output:
decomp/output/asitic_kernel.c lines 5650, 4010, 19). The
binary computes everything in extended precision using x87 f2xm1
/ fscale / fpatan intrinsics. We just call the libm
functions; float is IEEE-754 binary64 which is sufficient for
all but the most aggressive cancellation paths (Grover’s are not
those).
Conventions:
Lengths are passed in microns.
The closed forms internally convert to cm (multiplying by
1e-4) because Grover’s tables are tabulated in those units.Returns inductance in nH.
- Reference:
F. W. Grover, Inductance Calculations: Working Formulas and Tables, Dover, 1946. Sections on parallel filaments and rectangular bars (Tables 24, 25, formulas 25.4, 7.16).
The Greenhouse formulation used by ASITIC computes the total
inductance of a planar coil as the sum of self-inductances of each
straight segment plus all pairwise mutual inductances. This module
provides the per-segment primitives; the summation lives in
inductance.partial.
- reasitic.inductance.grover.rectangular_bar_self_inductance(length_um, width_um, thickness_um)[source]¶
Self-inductance of a rectangular conductor bar.
\[L = 2\ell\,\Bigl[ \ln\!\bigl(\tfrac{2\ell}{W+T}\bigr) + 0.50049 + \frac{W+T}{3\ell} \Bigr]\]where
ℓ,W,Tare in cm and the result is in nH. The constant0.50049is Grover’s value for the GMD correction of a thin rectangular cross-section. This is the exact formula used by the original ASITIC binary incmd_inductance_compute(decompasitic_repl.c:1417).
- reasitic.inductance.grover.segment_self_inductance(length_um, radius_um)[source]¶
Self-inductance of a single straight round wire (Grover §7.4).
\[L = 2 \ell\, \Bigl[ \sinh^{-1}(\ell/r) + r/\ell - \sqrt{1 + (r/\ell)^2} \Bigr]\]where
ℓis length in cm,ris the equivalent round-wire radius in cm, and the result is in nH.Mirrors the decompiled
grover_segment_self_inductanceat0x08064308.Example
>>> from reasitic.inductance import segment_self_inductance >>> round(segment_self_inductance(100.0, 0.5), 4) 0.0999
- reasitic.inductance.grover.coupled_wire_self_inductance(width_um, thickness_um, separation_um)[source]¶
Self-inductance of a rectangular bar of width
width_um, thicknessthickness_um, with a parallel-bar separation ofseparation_um(the latter affects the GMD via the proximity correction).Mirrors the decompiled
coupled_wire_self_inductance_groverat0x0804cb90(Grover Table 24).Parameters and result are in microns / nH.
- reasitic.inductance.grover.perpendicular_segment_mutual(length1_um, length2_um, *, common_distance_um=0.0)[source]¶
Mutual inductance between two perpendicular straight filaments.
Two filaments meeting at right angles (e.g. adjacent legs of a square spiral): segment 1 along the x-axis on
[0, L1], segment 2 along the y-axis on[0, L2], sharing the origin ifcommon_distance_umis zero, otherwise displaced along the common axis by that distance.For filamentary perpendicular wires Maxwell’s formula reduces to 0 because
ds₁ · ds₂ = 0. The binary’smutual_inductance_orthogonal_segments(asitic_kernel.c:0x08061b84) is a closed-form for finite- width bars — it captures the fringe / GMD contribution that appears at corners. For axis-aligned 2-D spirals this term is typically <1 % of the total inductance; we return zero here as the standard Greenhouse approximation. Callers that need the corner correction can substitute their own kernel.
- reasitic.inductance.grover.mohan_modified_wheeler(*, n_turns, d_outer_um, d_inner_um, shape='square')[source]¶
Mohan 1999 modified-Wheeler closed-form L estimate.
Example
>>> from reasitic.inductance import mohan_modified_wheeler >>> round(mohan_modified_wheeler( ... n_turns=5, d_outer_um=200, d_inner_um=100, ... ), 2) 5.75 >>> mohan_modified_wheeler( ... n_turns=5, d_outer_um=50, d_inner_um=100, ... ) 0.0
\[L_\text{mw} = K_1 \mu_0 n^2 d_\text{avg} / (1 + K_2 \rho)\]where
ρ = (d_out − d_in) / (d_out + d_in),d_avg = (d_out + d_in) / 2, and(K_1, K_2)come from Mohan 1999 Table 1:square: K_1 = 2.34, K_2 = 2.75hexagonal: K_1 = 2.33, K_2 = 3.82octagonal: K_1 = 2.25, K_2 = 3.55circular: K_1 = 2.40, K_2 = 1.75
Returns L in nH. Fast first-order estimate, useful for sanity-checking the Greenhouse summation.
- reasitic.inductance.grover.hoer_love_perpendicular_mutual(*, L1_um, L2_um, a_um, b_um, c_um)[source]¶
Hoer-Love mutual-inductance integral for perpendicular bars.
Two filaments meeting at a corner: filament 1 from
(0,0,0)to(L1,0,0), filament 2 from(a, b, c)to(a, b+L2, c). For perpendicular filaments the dot product is zero so M=0; for finite-width bars the proximity integral\[M_\perp = \frac{\mu_0}{4\pi}\int\int \frac{(\hat{x}\cdot\hat{y}) \,dx\,dy} {\sqrt{(x-a)^2 + (y-b)^2 + c^2}}\]is identically zero (the dot-product vanishes element-wise). Hoer & Love’s 1965 result captures the same vanishing for arbitrary 3-D orientations. We provide this function so the dispatch surface in
reasitic.inductance.partial._segment_pair_mutual()can call it without conditional logic; it always returns 0.The non-zero “corner” contribution that the binary’s
mutual_inductance_orthogonal_segmentsreports is actually a self-inductance artifact — the L_self of the bend itself, folded into a per-pair lookup. Callers wanting that should add a corner-correction term to the diagonal of the partial-L matrix, not the off-diagonal.
- reasitic.inductance.grover.parallel_segment_mutual(length1_um, length2_um, sep_um, offset_um=0.0)[source]¶
Mutual inductance between two parallel filamentary segments.
Both segments lie along the same axis. With segment 1 occupying
[0, L1]along the axis and segment 2 occupying[offset, offset + L2], separated by perpendicular distanced, the closed form is\[M = \tfrac{\mu_0}{4\pi}\, \bigl[\phi(L_1-o) - \phi(-o) - \phi(L_1-o-L_2) + \phi(-o-L_2)\bigr]\]where
φ(t) = t·asinh(t/d) − √(t² + d²)is the antiderivative of the double-integral kernel.φis even, so absolute values are used.The prefactor
μ₀/(4π)evaluates to 1 nH/cm, hence with all lengths in cm the result is in nH directly.For two equal-length parallel filaments with no axial offset this reduces to the canonical Greenhouse Eq. 8:
M = 2L·[asinh(L/d) − √(1 + (d/L)²) + d/L] (nH)
Greenhouse partial-inductance summation.
Total inductance of a planar coil is the sum of the per-segment self-inductances plus all signed pairwise mutual inductances:
The sign σ_ij is +1 when the two segments carry current in the
same direction along their common axis, -1 when they oppose, and 0
when they are perpendicular (parallel-segment Grover formula yields
zero in that case anyway). Reference: H. M. Greenhouse, “Design of
Planar Rectangular Microelectronic Inductors,” IEEE Trans. PHP-10
(1974).
The per-pair dispatch mirrors the binary’s check_segments_intersect
(decomp 0x08061110) and mutual_inductance_3d_segments
(0x08062ebc):
perpendicular pairs (
|û_a · û_b| < 1e-10) → 0axis-aligned parallel pairs → closed-form Grover 4-corner via
reasitic.inductance.grover.parallel_segment_mutual()general parallel-not-axis-aligned and skew pairs → Maxwell double integral via
reasitic.inductance.skew.mutual_inductance_skew_segments()
The skew kernel internally short-circuits exactly-parallel pairs to its own closed form, so the only pairs that hit the numerical integrator are genuinely non-parallel ones (typical only for polygon spirals and 3D / multi-metal structures).
- reasitic.inductance.partial.compute_self_inductance(shape)[source]¶
Total self-inductance (in nH) of a shape.
Uses Greenhouse partial-inductance summation. Mutual terms are routed through
_segment_pair_mutual(), which dispatches to the axis-aligned-parallel closed form for Manhattan spirals and to the general skew kernel for polygon / 3D geometry.
- reasitic.inductance.partial.compute_mutual_inductance(shape_a, shape_b)[source]¶
Mutual inductance between two distinct shapes, in nH.
Sums the signed Greenhouse mutual contribution over every cross-pair of segments. The double-counting factor of 2 used in the self case does not apply here because each (i, j) pair appears exactly once.
Mirrors the binary’s
cmd_coupling_compute(asitic_repl.c:1478) — parallel, perpendicular, and general skew geometries are all handled via_segment_pair_mutual().
- reasitic.inductance.partial.coupling_coefficient(shape_a, shape_b)[source]¶
Magnetic coupling coefficient
k = M / sqrt(L₁·L₂).Returns 0 if either self-inductance is non-positive. The result is dimensionless and lies in (-1, +1) for physically realisable geometries.
Filament discretisation and impedance-matrix-based L extraction.
The closed-form Greenhouse summation in reasitic.inductance.partial
treats every conductor as a single filament. For finite-width / finite-
thickness conductors at high frequency this misses current crowding
(skin and proximity effects). The original ASITIC binary handles this
by splitting each segment into N_w × N_t parallel filaments,
filling a complex impedance matrix Z = R + jωM, and solving for the
per-filament currents under a unit-voltage excitation. The net
inductance follows from the parallel sum.
This module mirrors that pipeline at a Python level:
Subdivide each
Segmentinton_w × n_taxial filaments (filament_grid).- Build the per-filament Z matrix:
Diagonal:
R_i + jωL_self_ifrom the rectangular-bar formula.Off-diagonal:
jωM_ijfrom the parallel-segment Grover form for parallel filaments, 0 for perpendicular pairs.
Solve
Z·I = V·1(uniform unit voltage), thenL_eff = Im(V/I_total) / ω.
Mirrors the binary’s solve_inductance_matrix
(asitic_kernel.c:5687, address 0x08064360) and
set_cell_size_normal (asitic_kernel.c:7602,
0x0807043c).
The sub-filaments are not used for the closed-form Greenhouse path
in reasitic.inductance.compute_self_inductance(); that path
still treats each segment as one filament and is faster but less
accurate at high frequency.
- class reasitic.inductance.filament.Filament[source]¶
Bases:
objectOne sub-filament after discretising a
Segment.
- reasitic.inductance.filament.auto_filament_subdivisions(segment, rsh_ohm_per_sq, freq_ghz, *, cells_per_skin_depth=2.0, n_max=8)[source]¶
Pick
(n_w, n_t)so each filament is ≤δ_skin / cells_per_skin_depth.At the operating frequency
freq_ghzand metal sheet resistancersh, computes the skin depthδand returns subdivisions that keep every filament ≤ that fraction of δ. Capped atn_maxto avoid runaway grid sizes.Returns
(1, 1)at DC (freq_ghz ≤ 0).
- reasitic.inductance.filament.auto_filament_subdivisions_critical(segment, rsh_ohm_per_sq, freq_ghz, *, cells_per_skin_depth=4.0, n_max=16)[source]¶
Critical-mode filament subdivisions (finer than the normal mode).
Mirrors
set_cell_size_critical(decomp0x08071cec, 3454 B): the binary’s stricter cell-sizer for high-accuracy regimes. Same formula asauto_filament_subdivisions()but with twice the cells-per-skin-depth target (default 4× rather than 2×) and a higher cap (n_max=16rather than 8). The binary logs the result with a leading!#marker that routes the line to the per-spiral log file.Returns the
(n_w, n_t)subdivision counts.
- reasitic.inductance.filament.filament_grid(segment, n_w=1, n_t=1)[source]¶
Split
segmentinto an_w × n_tgrid of filaments.Each filament inherits the parent segment’s length and direction; they’re translated perpendicular to the segment axis to tile the cross-section. Width and thickness divide evenly.
Returns at least one filament. Raises
ValueErrorif either division count is non-positive.
- reasitic.inductance.filament.build_inductance_matrix(filaments)[source]¶
Build the symmetric N×N partial-inductance matrix in nH.
Diagonal entries are per-filament self-inductances; off-diagonal entries are the signed Grover parallel-segment mutuals.
- reasitic.inductance.filament.build_resistance_vector(filaments, tech, freq_ghz)[source]¶
Per-filament AC resistance, ready to fill the Z diagonal.
- reasitic.inductance.filament.solve_inductance_mna(shape, tech, freq_ghz, *, n_w=1, n_t=1)[source]¶
Modified-nodal-analysis solve for net (L, R) of a series spiral.
Topology (rigorous version):
The spiral is a single mesh: every parent segment carries the same total current
I_spiral(chosen as 1 here).Within parent
pthen_par = n_w·n_tsub-filaments are in parallel: each carries some fractioni_kofI_spiral, and they all share the parent’s terminal voltage dropV_p.
For each parent we get
n_parequations:Σ_{k in p} i_k = 1 (current-conservation) (Z i)_k1 - (Z i)_k = 0 for k = k2 ... kn_par (voltages must match at the parent's end nodes)
where
Z = diag(R) + jω Mis the complex per-filament impedance matrix. Stacking these gives anN × Ncomplex linear system in the filament currents. The total spiral impedance is thenZ_eff = Σ_p V_p(withI_spiral = 1implicit).Mirrors
solve_inductance_matrix(0x08064360) more accurately than the oldersolve_inductance_matrixSchur- complement method.
- reasitic.inductance.filament.solve_inductance_matrix(shape, tech, freq_ghz, *, n_w=1, n_t=1)[source]¶
Filament-based solve for net (L, R) of a shape at one frequency.
Returns
(L_nH, R_ohm).Topology: the spiral is a single series chain of conductor segments. Within one parent segment the
n_w × n_tsub-filaments are in parallel (sharing the segment’s voltage drop); across segments the parents are in series (sharing the spiral’s terminal current). The net impedance is built by folding the per-filament Z matrix back through this topology.For
n_w = n_t = 1this reduces to the closed-form Greenhouse summation and matchescompute_self_inductance()exactly (within floating-point round-off).For higher subdivisions, the per-segment effective L includes proximity / current-crowding corrections via the parallel sharing inside each parent. Mirrors the binary’s
solve_inductance_matrix(0x08064360).
Eddy-current matrix and inductance fold-in.
When an inductor sits over a conductive substrate, image currents flow in the substrate that reduce the effective inductance and add loss. The original ASITIC binary models this via a ground- image method: each conductor’s current induces a phantom anti- parallel current at the substrate-image position; the resulting eddy-mutual-inductance matrix is folded into the diagonal of the spiral’s impedance matrix.
This is a rough first-cut implementation that places one image
filament per source filament, mirrored vertically about the
substrate’s top surface (z = 0 in our convention). The
substrate-conductivity-dependent attenuation factor is the standard
exp(-2|z|/δ) skin-depth roll-off applied to the image current.
Mirrors the simpler half of gen_eddy_current_matrix
(asitic_kernel.c:0x080b0e50) and inductance_eddy_fold
(asitic_kernel.c:12578).
- reasitic.inductance.eddy.assemble_eddy_matrix(shape, tech, freq_ghz, *, n_w=1, n_t=1)[source]¶
Build the per-filament substrate-image eddy mutual-L matrix.
Mirrors the binary’s
eddy_matrix_assemble(decomp0x080930a4, 1501 bytes): produces anN × Ncomplex matrix whose off-diagonals are the source-to-image mutual inductances between every filament pair, attenuated by the skin-depth factor through the substrate.- Parameters:
- Returns:
(N, N)real matrix of source-image mutual-L coupling in nH. Diagonal entries are the self-image term; off- diagonals are the cross-image. Multiply byjωand add to the partial-L matrix to get the full eddy-correctedZ.- Return type:
- reasitic.inductance.eddy.eddy_packed_index(i, j)[source]¶
Index calculator for the binary’s packed eddy matrix layout.
Mirrors
eddy_packed_index(decomp0x080941ec):if j < i: index = 8 * j − 4 * i + 3 else: index = 4 * i − 3
The eddy matrix is symmetric, so the off-diagonal block is stored in the lower-triangular half with a stride of 8 (two doubles per complex entry); the diagonal block has a stride of 4.
Returns the byte-stride-aware offset into the packed buffer.
- reasitic.inductance.eddy.eddy_correction(shape, tech, freq_ghz, *, n_w=1, n_t=1, finite_thickness=True)[source]¶
Eddy-current correction to (L, R) at
freq_ghz.Returns
(ΔL_nH, ΔR_ohm)— the change in inductance and resistance due to substrate eddy currents. Apply by adding to the un-foldedsolve_inductance_matrix()result.When
finite_thicknessis True the substrate is modelled as a finite-thickness ground plane: the image-current attenuation depends not only on the source-to-surface depth but also on the substrate’s finite thicknesst_sub, which limits the eddy current path. The effective attenuation factor becomes:attn = exp(-2·z/δ) · (1 - exp(-2·t_sub/δ))
where
δis the substrate skin depth andt_subis the bulk-layer thickness from the tech file. For thin substrates (t_sub ≪ δ) this reduces to a small correction; for thick substrates (t_sub ≫ δ) it recovers the half-space half-space limit.A negative
ΔL(inductance reduction) and positiveΔR(loss increase) are the typical signs.
- reasitic.inductance.eddy.solve_inductance_with_eddy(shape, tech, freq_ghz, *, n_w=1, n_t=1, include_eddy=True)[source]¶
Like
solve_inductance_matrix()but also folds in the eddy correction fromeddy_correction(). Returns(L_nH, R_ohm).