helicon.fields
Magnetic field computation: Biot-Savart solver, field line tracing, FRC topology, and caching.
Core API
helicon.fields.biot_savart.Coil(z, r, I)
dataclass
A single circular current loop in SI units.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
z
|
float
|
Axial position of the coil centre [m]. |
required |
r
|
float
|
Radius of the coil [m]. |
required |
I
|
float
|
Current (ampere-turns) [A]. |
required |
helicon.fields.biot_savart.Grid(z_min, z_max, r_max, nz, nr)
dataclass
Axisymmetric (r, z) computation grid.
The grid spans z ∈ [z_min, z_max] and r ∈ [0, r_max].
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
z_min
|
float
|
Axial extent [m]. |
required |
z_max
|
float
|
Axial extent [m]. |
required |
r_max
|
float
|
Radial extent [m]. |
required |
nz
|
int
|
Number of grid points along each axis. |
required |
nr
|
int
|
Number of grid points along each axis. |
required |
helicon.fields.biot_savart.BField(Br, Bz, r, z, coils, backend)
dataclass
Result container — always stores NumPy arrays for interop.
Attributes:
| Name | Type | Description |
|---|---|---|
Br, Bz |
ndarray
|
Radial and axial magnetic field components [T], shape |
r, z |
ndarray
|
1-D coordinate arrays [m]. |
coils |
list[Coil]
|
Coil configuration that produced this field. |
backend |
str
|
Backend used ( |
Functions
load(path)
classmethod
Load a previously saved BField from HDF5.
Source code in src/helicon/fields/biot_savart.py
plot(*, component='Bz', field_lines=True, n_field_lines=8, ax=None, cmap='RdBu_r', figsize=(10, 5))
Plot the B-field on the (z, r) plane.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
component
|
``"Bz"`` | ``"Br"`` | ``"Bmag"``
|
Field component to show as a filled contour. |
'Bz'
|
field_lines
|
bool
|
Overlay field lines via contours of the flux function ψ. |
True
|
n_field_lines
|
int
|
Number of field line contours to draw. |
8
|
ax
|
matplotlib Axes
|
Axes to draw on. Creates new figure if None. |
None
|
cmap
|
str
|
Matplotlib colormap. |
'RdBu_r'
|
figsize
|
tuple
|
Figure size when creating a new figure. |
(10, 5)
|
Returns:
| Type | Description |
|---|---|
(fig, ax)
|
Matplotlib figure and axes. |
Source code in src/helicon/fields/biot_savart.py
114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 | |
save(path)
Write field data + coil metadata to an HDF5 file.
Source code in src/helicon/fields/biot_savart.py
helicon.fields.biot_savart.compute_bfield(coils, grid, *, backend='auto', n_phi=128, differentiable=False)
Compute the applied magnetic field from a set of circular coils.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
coils
|
sequence of Coil
|
Coil definitions in SI units. |
required |
grid
|
Grid
|
Computational domain. |
required |
backend
|
``"auto"`` | ``"mlx"`` | ``"numpy"``
|
|
'auto'
|
n_phi
|
int
|
Azimuthal quadrature points (MLX backend only). |
128
|
differentiable
|
bool
|
When True, use the loop-based MLX path so the computation graph is
preserved for |
False
|
Returns:
| Type | Description |
|---|---|
BField
|
Result container with NumPy arrays. |
Source code in src/helicon/fields/biot_savart.py
helicon.fields.biot_savart.compute_bfield_mlx_differentiable(coil_params, grid_r, grid_z, n_phi=128)
MLX-native differentiable Biot-Savart computation.
Designed for composition with mlx.core.grad(). Does not call
mx.eval so the computation graph is preserved.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
coil_params
|
(array, shape(N_coils, 3))
|
Each row is |
required |
grid_r
|
(array, shape(N_pts))
|
Flattened grid coordinates. |
required |
grid_z
|
(array, shape(N_pts))
|
Flattened grid coordinates. |
required |
n_phi
|
int
|
Number of azimuthal quadrature points. |
128
|
Returns:
| Type | Description |
|---|---|
Br, Bz : mx.array, shape (N_pts,)
|
|
Source code in src/helicon/fields/biot_savart.py
Field Cache
helicon.fields.cache.FieldCache(cache_dir=None)
On-disk HDF5 cache for BField results.
Source code in src/helicon/fields/cache.py
Functions
clear()
get(coils, grid)
Load a cached BField, or return None on miss / corrupt file.
Source code in src/helicon/fields/cache.py
key(coils, grid)
SHA-256 cache key (first 16 hex chars) from coil + grid params.
Source code in src/helicon/fields/cache.py
put(coils, grid, bfield)
helicon.fields.cache.compute_bfield_cached(coils, grid, *, cache=None, backend='auto', n_phi=128)
Like compute_bfield but with transparent HDF5 caching.
Source code in src/helicon/fields/cache.py
FRC Topology
helicon.fields.frc_topology.FRCTopologyResult(is_open, is_closed, psi_rz, o_point_z, o_point_r, x_point_z, psi_separatrix)
dataclass
Topology classification of an FRC magnetic field.
helicon.fields.frc_topology.find_frc_topology(bfield)
Classify every grid point as open or closed for an FRC field.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
bfield
|
BField
|
Magnetic field on an axisymmetric (r, z) grid. |
required |
Returns:
| Type | Description |
|---|---|
FRCTopologyResult
|
|
Source code in src/helicon/fields/frc_topology.py
52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 | |
Field Lines
helicon.fields.field_lines.trace_field_lines(bfield, *, n_lines=20, z_start=None, r_range=None, max_length=100.0, max_steps=10000, ds=0.001)
Trace multiple field lines from evenly spaced starting points.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
bfield
|
BField
|
Magnetic field data. |
required |
n_lines
|
int
|
Number of field lines to trace. |
20
|
z_start
|
float
|
Axial position to launch lines from. Defaults to domain midpoint. |
None
|
r_range
|
(r_min, r_max)
|
Radial range for starting points. Defaults to (0, r_max). |
None
|
max_length
|
float
|
Maximum arc length to trace [m]. |
100.0
|
max_steps
|
int
|
Maximum integration steps. |
10000
|
ds
|
float
|
Initial step size [m]. |
0.001
|
Source code in src/helicon/fields/field_lines.py
helicon.fields.field_lines.FieldLine(r, z, line_type, start_r, start_z, psi_start)
dataclass
A single traced field line.
helicon.fields.field_lines.FieldLineSet(lines, psi, r, z, separatrix_psi)
dataclass
Collection of traced field lines with topology info.
External Import
helicon.fields.import_external.load_csv_bfield(path, *, r_col='r', z_col='z', Br_col='Br', Bz_col='Bz', delimiter=',')
Load a B-field map from a CSV file.
The CSV must have a header row and at minimum four columns for
r, z, Br, Bz. Data must lie on a regular meshgrid —
all combinations of the unique r and z values must be present.
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
path
|
str or Path
|
|
required |
r_col
|
str
|
Column names in the header row. |
'r'
|
z_col
|
str
|
Column names in the header row. |
'r'
|
Br_col
|
str
|
Column names in the header row. |
'r'
|
Bz_col
|
str
|
Column names in the header row. |
'r'
|
delimiter
|
str
|
Field delimiter (default |
','
|
Returns:
| Type | Description |
|---|---|
BField
|
Arrays shaped |
Source code in src/helicon/fields/import_external.py
helicon.fields.import_external.load_femm_bfield(path, *, length_scale=0.001)
Load a B-field map exported from FEMM in text (.ans) format.
FEMM axisymmetric analyses can export field data as a tab-delimited
table. The expected columns are r z Br Bz |B| (coordinates
in mm by default).
Parameters:
| Name | Type | Description | Default |
|---|---|---|---|
path
|
str or Path
|
Path to the FEMM |
required |
length_scale
|
float
|
Multiply coordinates by this factor to convert to SI metres. Default 1e-3 converts FEMM's default millimetres → metres. |
0.001
|
Returns:
| Type | Description |
|---|---|
BField
|
|