API Documentation

Classes

ElectrodeConfig

class ElectrodeConfig

Dataclass for configuring a cylindrically symmetric electrostatic lens element.

Parameters

Parameters:
  • start (float) – Z-axis starting position in grid units

  • width (float) – Z-axis width in grid units

  • ap_start (float) – R-axis aperture starting position in grid units

  • ap_width (float) – R-axis aperture width in grid units

  • outer_diameter (float) – Full electrode diameter in grid units

  • voltage (float) – Applied voltage in volts

MagneticLensConfig

class MagneticLensConfig

Dataclass for configuring a cylindrically symmetric magnetic lens element.

Parameters

Parameters:
  • start (float) – Position where the magnetic lens begins on the z-axis in grid units

  • length (float) – Length of the magnetic lens on the z-axis in grid units

  • ap_start (float) – Starting position of the aperture on the r-axis in grid units

  • ap_width (float) – Width of the aperture in the lens on the r-axis in grid units

  • outer_diameter (float) – Diameter of the magnetic lens on the r-axis in grid units

  • mu_r (float) – Relative magnetic permeability of the lens material in dimensionless units

  • mmf (float) – Magnetomotive force of the lens in ampere-turns

ElectricField

class ElectricField(nz, nr, axial_size, radial_size)

Manages electric potential fields, electrode configurations, and field solving operations.

Parameters

Parameters:
  • nz (float) – Number of z-axis grid points

  • nr (float) – Number of r-axis grid points

  • axial_size (float) – Physical z-axis length in meters

  • radial_size (float) – Physical r-axis length in meters

Attributes

nz: int

Number of z-axis grid points

nr: int

Number of r-axis grid points

axial_size: float

Z-axis size in meters

radial_size: float

R-axis size in meters

dz: float

Z-axis grid spacing (meters), calculated as axial_size/nz

dr: float

R-axis grid spacing (meters), calculated as radial_size/nr

potential: numpy.ndarray

2D array of electric potential values in volts

electrode_mask: numpy.ndarray

Boolean mask indicating electrode positions

Ez: numpy.ndarray

Z-component of electric field in V/m

Er: numpy.ndarray

R-component of electric field in V/m

Methods

add_electrode(config)

Adds an electrode to the field configuration.

Parameters:

config (ElectrodeConfig) – Electrode configuration parameters

build_laplacian_matrix(mask, dirichlet_values=None)

Constructs sparse matrix representation of the Laplacian operator for solving Laplace’s equation. Implements Dirichlet boundary conditions (0V) at all boundaries to simulate grounded metal boundaries.

Parameters:
  • mask (numpy.ndarray) – Boolean array marking electrode positions

  • dirichlet_values (numpy.ndarray) – Optional potential values at masked positions

Returns:

Sparse matrix A and right-hand side vector b

Return type:

tuple(scipy.sparse.csr_matrix, numpy.ndarray)

solve_potential(max_iterations=500, convergence_threshold=1e-6)

Solves the electrostatic potential field using PyAMG multigrid methods. Complexity is O(N) with grid size.

Parameters:
  • max_iterations (float) – Maximum solver iterations (default: 500)

  • convergence_threshold (float) – Convergence criterion (default: 1e-6)

Returns:

Solved potential field array

Return type:

numpy.ndarray

get_field_at_position(z, r)

Returns electric field components at a specific position.

Parameters:
  • z (float) – Z-axis position in meters

  • r (float) – R-axis position in meters

Returns:

Electric field components (Ez, Er)

Return type:

tuple(float, float)

MagneticField

class MagneticField(electron_optics)

Manages magnetic vector potential fields, magnetic lens configurations, and field solving operations.

Parameters

Parameters:

electron_optics (ElectronOptics) – The parent electron optics system containing field dimensions and grid parameters

Attributes

nz: int

Number of z-axis grid points

nr: int

Number of r-axis grid points

axial_size: float

Z-axis size in meters

radial_size: float

R-axis size in meters

dz: float

Z-axis grid spacing (meters)

dr: float

R-axis grid spacing (meters)

vector_potential: numpy.ndarray

2D array of magnetic vector potential values in Tesla-meters

magnetic_mask: numpy.ndarray

Boolean mask indicating magnetic material positions

mu_r: numpy.ndarray

2D array of relative permeability values

current_density: numpy.ndarray

2D array of current density values in amperes per square meter

Bz: numpy.ndarray

Z-component of magnetic field in Tesla

Br: numpy.ndarray

R-component of magnetic field in Tesla

lens_config: MagneticLensConfig

Configuration of the magnetic lens

Methods

add_magnetic_lens(config)

Adds a magnetic lens to the field and handles all necessary calculations.

Parameters:

config (MagneticLensConfig) – Configuration parameters for the magnetic lens

build_laplacian_matrix(mask, dirichlet_values=None)

Builds a sparse matrix for the Laplacian ∇²A = -μ₀μᵣJ, and implements Neumann boundary conditions at all boundaries.

Parameters:
  • mask (numpy.ndarray) – Boolean array, true where magnetic materials exist

  • dirichlet_values (numpy.ndarray) – Vector potential values where mask is True

Returns:

Sparse matrix A and right-hand side vector b

Return type:

tuple(scipy.sparse.csr_matrix, numpy.ndarray)

solve_vector_potential(max_iterations=500, convergence_threshold=1e-6)

Solves the magnetic vector potential field using Multigrid methods with PyAMG.

Parameters:
  • max_iterations (float) – Maximum number of iterations for the solver (default: 500)

  • convergence_threshold (float) – Convergence criterion for the solution (default: 1e-6)

Returns:

The solved vector potential field

Return type:

numpy.ndarray

calculate_b_from_a()

Calculates the magnetic field components from the vector potential using B = ∇ × A, with special handling of differentiation at boundaries and at r = 0.

get_field_at_position(z, r)

Returns the magnetic field components at a specific position.

Parameters:
  • z (float) – Position along the z-axis in meters

  • r (float) – Position along the r-axis in meters

Returns:

The magnetic field components (Bz, Br) at the specified position

Return type:

tuple(float, float)

ParticleTracer

class ParticleTracer(electric_field)

Handles charged particle trajectory calculations and dynamics simulations.

Parameters

Parameters:

electric_field (ElectricField) – Electric field for particle simulation

Class Constants

ELECTRON_CHARGE = -1.60217663e-19

Elementary charge in Coulombs

ELECTRON_MASS = 9.1093837e-31

Electron rest mass in kilograms

SPEED_OF_LIGHT = 299792458.0

Speed of light in vacuum (m/s)

Attributes

field: ElectricField

Associated potential field instance

current_ion: dict

Current particle properties including symbol, atomic_number, mass, charge, and charge_mass_ratio

q: float

Particle charge in Coulombs

m: float

Particle mass in kilograms

Methods

set_ion(symbol='e-', charge_state=1)

Configures the tracer for a specific ion or electron. Integrates with Mendeleev library for atomic data.

Parameters:
  • symbol (str) – Chemical symbol or ‘e-’ for electrons

  • charge_state (float) – Ion charge state

Returns:

Self for method chaining

Return type:

ParticleTracer

get_velocity_from_energy(energy_eV)

Converts kinetic energy to velocity with relativistic corrections. Accurate for all energy scales from single-digit eV to GeV.

Parameters:

energy_eV (float) – Kinetic energy in electronvolts

Returns:

Particle velocity in m/s

Return type:

float

particle_dynamics(t, state)

Calculates particle dynamics for ODE solver. Uses relativistic equations of motion.

Parameters:
  • t (float) – Current time

  • state (list) – State vector [z, r, pz, pr]

Returns:

State derivatives [vz, vr, dpz_dt, dpr_dt]

Return type:

list

trace_trajectory(initial_position, initial_velocity, simulation_time, method='BDF', rtol=1e-9, atol=1e-12)

Solves particle equations of motion in the electric field.

Parameters:
  • initial_position (tuple) – Initial (z, r) position in meters

  • initial_velocity (tuple) – Initial (vz, vr) velocity in m/s

  • simulation_time (float) – Total simulation time in seconds

  • method (str) – ODE solver method (default: ‘BDF’)

  • rtol (float) – Relative tolerance (default: 1e-9)

  • atol (float) – Absolute tolerance (default: 1e-12)

Returns:

ODE solution object

Return type:

scipy.integrate.OdeResult

EinzelLens

class EinzelLens(position, width, aperture_center, aperture_width, outer_diameter, focus_voltage, gap_size=1)

Three-electrode einzel (unipotential) lens system for particle focusing.

Parameters

Parameters:
  • position (float) – Starting z-axis position in grid units

  • width (float) – Total lens width in grid units

  • aperture_center (float) – Aperture center on r-axis in grid units

  • aperture_width (float) – Aperture size in grid units

  • outer_diameter (float) – Electrode diameter in grid units

  • focus_voltage (float) – Central electrode voltage in volts

  • gap_size (int) – Inter-electrode gap size in grid units (default: 1)

Attributes

electrode1: ElectrodeConfig

First electrode configuration (0V)

electrode2: ElectrodeConfig

Central electrode configuration (focus_voltage)

electrode3: ElectrodeConfig

Third electrode configuration (0V)

Methods

add_to_field(field)

Adds all three electrodes to the electric field.

Parameters:

field (ElectricField) – Target electric field

ElectronOptics

class ElectronOptics(nr, nz, axial_size=0.1, radial_size=0.1)

Complete electron optics simulation environment integrating field calculations, particle tracing, and visualization capabilities.

Parameters

Parameters:
  • nr (float) – Number of r-axis grid points

  • nz (float) – Number of z-axis grid points

  • axial_size (float) – Z-axis size in meters (default: 0.1)

  • radial_size (float) – R-axis size in meters (default: 0.1)

Attributes

field: ElectricField

Electric field instance

tracer: ParticleTracer

Particle trajectory calculator

elements: list

List of all system electrodes and lenses

magnetic_lenses: MagneticField

Magnetic field instance

Methods

add_magnetic_lens(config)

Adds a magnetic lens to the system.

Parameters:

config (MagneticLensConfig) – Magnetic lens configuration

add_electrode(config)

Adds an electrode to the system.

Parameters:

config (ElectrodeConfig) – Electrode configuration

add_einzel_lens(position, width, aperture_center, aperture_width, outer_diameter, focus_voltage, gap_size=1)

Adds an einzel lens to the system.

Parameters:
  • position (float) – Starting z-axis position in grid units

  • width (float) – Total lens width in grid units

  • aperture_center (float) – Aperture center on r-axis in grid units

  • aperture_width (float) – Aperture size in grid units

  • outer_diameter (float) – Electrode diameter in grid units

  • focus_voltage (float) – Central electrode voltage in volts

  • gap_size (int) – Inter-electrode gap size (default: 1)

solve_fields()

Solves the potential and/or vector potential fields using PyAMG multigrid solver.

Returns:

Solved potential field, vector potential field, or a dictionary containing both

Return type:

numpy.ndarray or dict

simulate_beam(energy_eV, start_z, r_range, angle_range, num_particles, simulation_time, n_jobs=-1)

Simulates a particle beam with specified parameters. Uses parallel processing for multiple trajectories.

Parameters:
  • energy_eV (float) – Particle kinetic energy in electronvolts

  • start_z (float) – Starting z-position in meters

  • r_range (tuple) – Range of initial r-positions in meters

  • angle_range (tuple) – Range of initial angles in degrees

  • num_particles (float) – Number of particles to simulate

  • simulation_time (float) – Total simulation time in seconds

  • n_jobs (int) – Number of parallel jobs (default: -1 for all cores)

Returns:

List of trajectory solutions

Return type:

list

visualize_system(trajectories=None, r_limits=None, figsize=(16, 10), title='Picht')

Creates interactive visualization of the system and particle trajectories.

Parameters:
  • trajectories (list) – Optional trajectory solutions to display

  • r_limits (tuple) – Optional y-axis limits in meters

  • figsize (tuple) – Figure dimensions in inches

  • title (str) – Plot title

Returns:

Matplotlib figure object

Return type:

matplotlib.figure.Figure

Export

class Export(electron_optics)

Handles data export to HDF5 and CAD formats.

Parameters

Parameters:

electron_optics (ElectronOptics) – System instance to export

Attributes

system

Reference to the electron optics system

field

Shortcut to system’s potential field

magnetic_lenses

Shortcut to system’s magnetic field

Methods

export_traj(trajectories)

Exports simulation data to HDF5 format including fields, trajectories, system configuration, and particle properties.

Parameters:

trajectories (list) – List of trajectory solutions

Output:

Creates ‘simulation_results.h5’ file

cad_export()

Exports electrode geometry to STEP format for CAD integration.

Output:

Creates ‘save.step’ file

_create_electrode_shape(electrode_config)

Converts 2D electrode config data to 3D CAD geometric data.

Parameters:

electrode_config (ElectrodeConfig) – ElectrodeConfig data in 2D axisymmetric data

Returns:

3D CAD solid representation of the electrode

_create_magnetic_lens_shape(magnetic_config)

Converts 2D magnetic lens config data to 3D CAD geometric data.

Parameters:

magnetic_config (MagneticLensConfig) – MagneticLensConfig data in 2D axisymmetric format

Returns:

3D CAD solid representation of the magnetic lens

Misc. Functions

get_field(z, r, Ez, Er, axial_size, radial_size, dz, dr, nz, nr)

Retrieves electric field values at fractional grid positions using nearest-neighbor interpolation.

Parameters:
  • z (float) – Z-axis position in meters

  • r (float) – R-axis position in meters

  • Ez (numpy.ndarray) – Z-component field array

  • Er (numpy.ndarray) – R-component field array

  • axial_size (float) – Total z-axis size in meters

  • radial_size (float) – Total r-axis size in meters

  • dz (float) – Z-axis grid spacing

  • dr (float) – R-axis grid spacing

  • nz (int) – Number of z-axis grid points

  • nr (int) – Number of r-axis grid points

Returns:

Electric field components (Ez, Er) in V/m

Return type:

tuple(float, float)

calc_dynamics(z, r, pz, pr, Ez, Er, Bz, Br, q, m, c, r_axis=0.0)

Calculates charged particle acceleration by applying the Lorentz force with special-relativistic corrections. Uses energy momentum formalism for full eV to TeV support.

Parameters:
  • z (float) – Z-axis position in meters

  • r (float) – R-axis position in meters

  • pz (float) – Z-axis momentum in kg⋅m/s

  • pr (float) – R-axis momentum in kg⋅m/s

  • Ez (float) – Z-axis electric field in V/m

  • Er (float) – R-axis electric field in V/m

  • Bz (float) – Z-axis magnetic field in Tesla

  • Br (float) – R-axis magnetic field in Tesla

  • q (float) – Particle charge in Coulombs

  • m (float) – Particle mass in kilograms

  • c (float) – Speed of light in m/s

  • r_axis (float) – Reference r-axis position (default: 0.0)

Returns:

Array [vz, vr, dpz_dt, dpr_dt] with velocities and forces

Return type:

numpy.ndarray