ebisim.simulation package

This subpackage contains functions and classes provide an interface to run simulations and inspect their results.

class ebisim.simulation.AdvancedModel(device, targets, bg_gases, options, lb, ub, nq, q, a, eixs, rrxs, drxs, cxxs_bggas, cxxs_trgts)[source]

Bases: tuple

The advanced model class is the base for ebisim.simulation.advanced_simulation. It acts as a fast datacontainer for the underlying rhs function which represents the right hand side of the differential equation system. Since it is jitcompiled using numba, care is required during instatiation.

Parameters
  • device (ebisim.simulation.Device) – Container describing the EBIS/T and specifically the electron beam.

  • targets (numba.typed.List[ebisim.simulation.Target]) – List of ebisim.simulation.Target for which charge breeding is simulated.

  • bg_gases (numba.typed.List[ebisim.simulation.BackgroundGas], optional) – List of ebisim.simulation.BackgroundGas which act as CX partners, by default None.

  • options (ebisim.simulation.ModelOptions, optional) – Switches for effects considered in the simulation, see default values of ebisim.simulation.ModelOptions.

  • lb (ndarray) –

  • ub (ndarray) –

  • nq (int) –

  • q (ndarray) –

  • a (ndarray) –

  • eixs (ndarray) –

  • rrxs (ndarray) –

  • drxs (ndarray) –

  • cxxs_bggas (Any) –

  • cxxs_trgts (Any) –

Create new instance of AdvancedModel(device, targets, bg_gases, options, lb, ub, nq, q, a, eixs, rrxs, drxs, cxxs_bggas, cxxs_trgts)

count(value, /)

Return number of occurrences of value.

classmethod get(device, targets, bg_gases=None, options=ModelOptions(EI=True, RR=True, CX=True, DR=False, SPITZER_HEATING=True, COLLISIONAL_THERMALISATION=True, ESCAPE_AXIAL=True, ESCAPE_RADIAL=True, RECOMPUTE_CROSS_SECTIONS=False, RADIAL_DYNAMICS=False, IONISATION_HEATING=True, OVERRIDE_FWHM=False, RADIAL_SOLVER_MAX_STEPS=500, RADIAL_SOLVER_REL_DIFF=0.001))[source]
Parameters
Return type

AdvancedModel

index(value, start=0, stop=9223372036854775807, /)

Return first index of value.

Raises ValueError if the value is not present.

property a

Alias for field number 8

property bg_gases

Alias for field number 2

property cxxs_bggas

Alias for field number 12

property cxxs_trgts

Alias for field number 13

property device

Alias for field number 0

property drxs

Alias for field number 11

property eixs

Alias for field number 9

property lb

Alias for field number 4

property nq

Alias for field number 6

property options

Alias for field number 3

property q

Alias for field number 7

property rrxs

Alias for field number 10

property targets

Alias for field number 1

property ub

Alias for field number 5

class ebisim.simulation.AdvancedResult(*, t, N, device, target, kbT, model, id_, res=None, rates=None)[source]

Bases: BasicResult[Device]

Instances of this class are containers for the results of ebisim advanced_simulations and contain a variety of convenience methods for simple plot generation etc.

Parameters
  • t – An array holding the time step coordinates.

  • N – An array holding the occupancy of each charge state at a given time.

  • device – If coming from advanced_simulation, this is the machine / electron beam description.

  • target – If coming from advanced_simulation, this is the target represented by this Result.

  • kbT – An array holding the temperature of each charge state at a given time.

  • model – The AdvancedModel instance that was underlying the advanced_simulation

  • id – The position of this target in the list of all simulated targets.

  • res – The result object returned by scipy.integrate.solve_ivp. This can contain useful information about the solver performance etc. Refer to scipy documentation for details.

  • rates – A dictionary containing the different breeding rates in arrays shaped like N.

abundance_at_time(t)

Yields the abundance of each charge state at a given time

Parameters

t (float) – <s> Point of time to evaluate.

Returns

Abundance of each charge state, array index corresponds to charge state.

Return type

ndarray

plot(relative=False, **kwargs)

Plot the charge state evolution of this result object.

Parameters
  • relative (bool) – Flags whether the absolute numbers or a relative fraction should be plotted at each timestep, by default False.

  • kwargs (Any) – Keyword arguments are handed down to ebisim.plotting.plot_generic_evolution, cf. documentation thereof. If no arguments are provided, reasonable default values are injected.

Returns

Figure handle of the generated plot.

Raises

ValueError – If the required data (self.t, self.N) is not available, i.e. corresponding attributes of this Result instance have not been set correctly.

Return type

Figure

plot_charge_states(relative=False, **kwargs)

Plot the charge state evolution of this result object.

Parameters
  • relative (bool) – Flags whether the absolute numbers or a relative fraction should be plotted at each timestep, by default False.

  • kwargs (Any) – Keyword arguments are handed down to ebisim.plotting.plot_generic_evolution, cf. documentation thereof. If no arguments are provided, reasonable default values are injected.

Returns

Figure handle of the generated plot.

Raises

ValueError – If the required data (self.t, self.N) is not available, i.e. corresponding attributes of this Result instance have not been set correctly.

Return type

Figure

plot_energy_density(**kwargs)[source]

Plot the energy density evolution of this result object.

Parameters

kwargs (Any) – Keyword arguments are handed down to ebisim.plotting.plot_generic_evolution, cf. documentation thereof. If no arguments are provided, reasonable default values are injected.

Returns

Figure handle of the generated plot.

Raises

ValueError – If the required data (self.t, self.kbT) is not available, i.e. corresponding attributes of this Result instance have not been set correctly.

Return type

Figure

plot_radial_distribution_at_time(t, **kwargs)[source]

Plot the radial ion distribution at time t.

Parameters
  • t (float) – <s> Point of time to evaluate.

  • kwargs (Any) – Keyword arguments are handed down to ebisim.plotting.plot_radial_distribution and ebisim.plotting.plot_generic_evolution, cf. documentation thereof. If no arguments are provided, reasonable default values are injected.

Returns

Figure handle of the generated plot.

Return type

Figure

plot_rate(rate_key, dens_threshold=0.001, **kwargs)[source]

Plots the requested ionisation- or energy flow rates.

Parameters
  • rate_key (Rate) – The key identifying the rate to be plotted. See ebisim.simulation.Rate for valid values.

  • dens_threshold (float) – If given temperatures are only plotted where the particle denisty is larger than the threshold value.

  • kwargs – Keyword arguments are handed down to ebisim.plotting.plot_generic_evolution, cf. documentation thereof. If no arguments are provided, reasonable default values are injected.

Returns

Figure handle of the generated plot

Raises

ValueError – If the required data (self.rates) is not available, or an invalid key is requested.

Return type

Figure

plot_temperature(dens_threshold=0.001, **kwargs)[source]

Plot the temperature evolution of this result object.

Parameters
  • dens_threshold (float) – If given temperatures are only plotted where the particle denisty is larger than the threshold value.

  • kwargs (Any) – Keyword arguments are handed down to ebisim.plotting.plot_generic_evolution, cf. documentation thereof. If no arguments are provided, reasonable default values are injected.

Returns

Figure handle of the generated plot

Raises

ValueError – If the required data (self.t, self.kbT) is not available, i.e. corresponding attributes of this Result instance have not been set correctly.

Return type

Figure

radial_distribution_at_time(t)[source]

Yields the radial distribution information at time

Parameters

t (float) – <s> Point of time to evaluate.

Returns

  • phi – Radial potential

  • n3d – On axis 3D density for each charge state.

  • shapes – The Boltzmann shape factors for each charge state.

Return type

Tuple[ndarray, ndarray, ndarray]

temperature_at_time(t)[source]

Yields the temperature of each charge state at a given time

Parameters

t (float) – <s> Point of time to evaluate.

Returns

  • <eV>

  • Temperature of each charge state, array index corresponds to charge state.

Return type

ndarray

times_of_highest_abundance()

Yields the point of time with the highest abundance for each charge state

Returns

  • <s>

  • Array of times.

Return type

ndarray

class ebisim.simulation.BackgroundGas(name, ip, n0)[source]

Bases: tuple

Use the static get() factory methods to create instances of this class.

Simple datacontainer for a background gas for advanced simulations. A background gas only acts as a charge exchange partner to the Targets in the simulation.

Create new instance of BackgroundGas(name, ip, n0)

Parameters
  • name (str) –

  • ip (float) –

  • n0 (float) –

count(value, /)

Return number of occurrences of value.

classmethod get(element, p, T=300.0)[source]

Factory method for defining a background gas.

Parameters
  • element (Union[Element, str, int]) – An instance of the Element class, or an identifier for the element, i.e. either its name, symbol or proton number.

  • p (float) – <mbar> Gas pressure.

  • T (float) – <K> Gas temperature, by default 300 K (Room temperature)

Returns

ebisim.simulation.BackgroundGas – Ready to use BackgroundGas specification.

Return type

BackgroundGas

index(value, start=0, stop=9223372036854775807, /)

Return first index of value.

Raises ValueError if the value is not present.

property ip

<eV> Ionisation potential of this Gas.

property n0

<1/m^3> Gas number density.

property name

Name of the element.

class ebisim.simulation.BasicResult(*, t, N, device, target, res=None)[source]

Bases: Generic[GenericDevice]

Instances of this class are containers for the results of ebisim basic_simulations and contain a variety of convenience methods for simple plot generation etc.

Parameters
  • t – An array holding the time step coordinates.

  • N – An array holding the occupancy of each charge state at a given time.

  • device – If coming from advanced_simulation, this is the machine / electron beam description.

  • target – If coming from advanced_simulation, this is the target represented by this Result.

  • res – The result object returned by scipy.integrate.solve_ivp. This can contain useful information about the solver performance etc. Refer to scipy documentation for details.

abundance_at_time(t)[source]

Yields the abundance of each charge state at a given time

Parameters

t (float) – <s> Point of time to evaluate.

Returns

Abundance of each charge state, array index corresponds to charge state.

Return type

ndarray

plot(relative=False, **kwargs)

Alias for plot_charge_states

Parameters
  • relative (bool) –

  • kwargs (Any) –

Return type

Figure

plot_charge_states(relative=False, **kwargs)[source]

Plot the charge state evolution of this result object.

Parameters
  • relative (bool) – Flags whether the absolute numbers or a relative fraction should be plotted at each timestep, by default False.

  • kwargs (Any) – Keyword arguments are handed down to ebisim.plotting.plot_generic_evolution, cf. documentation thereof. If no arguments are provided, reasonable default values are injected.

Returns

Figure handle of the generated plot.

Raises

ValueError – If the required data (self.t, self.N) is not available, i.e. corresponding attributes of this Result instance have not been set correctly.

Return type

Figure

times_of_highest_abundance()[source]

Yields the point of time with the highest abundance for each charge state

Returns

  • <s>

  • Array of times.

Return type

ndarray

class ebisim.simulation.Device(current, e_kin, r_e, length, j, v_ax, v_ax_sc, v_ra, b_ax, r_dt, r_dt_bar, fwhm, rad_grid, rad_fd_l, rad_fd_d, rad_fd_u, rad_phi_uncomp, rad_phi_ax_barr, rad_re_idx)[source]

Bases: tuple

Use the static get() factory methods to create instances of this class.

Objects of this class are used to pass important EBIS/T parameters into the simulation.

Create new instance of Device(current, e_kin, r_e, length, j, v_ax, v_ax_sc, v_ra, b_ax, r_dt, r_dt_bar, fwhm, rad_grid, rad_fd_l, rad_fd_d, rad_fd_u, rad_phi_uncomp, rad_phi_ax_barr, rad_re_idx)

Parameters
  • current (float) –

  • e_kin (float) –

  • r_e (float) –

  • length (Optional[float]) –

  • j (float) –

  • v_ax (float) –

  • v_ax_sc (float) –

  • v_ra (float) –

  • b_ax (float) –

  • r_dt (float) –

  • r_dt_bar (float) –

  • fwhm (float) –

  • rad_grid (ndarray) –

  • rad_fd_l (ndarray) –

  • rad_fd_d (ndarray) –

  • rad_fd_u (ndarray) –

  • rad_phi_uncomp (ndarray) –

  • rad_phi_ax_barr (ndarray) –

  • rad_re_idx (int) –

count(value, /)

Return number of occurrences of value.

classmethod get(*, current, e_kin, r_e, v_ax, b_ax, r_dt, length=None, v_ra=None, j=None, fwhm=None, n_grid=400, r_dt_bar=None)[source]

Factory method for defining a device. All arguments are keyword only to reduce the chance of mixing them up.

Parameters
  • current (float) – <A> Electron beam current.

  • e_kin (float) – <eV> Uncorrected electron beam energy.

  • r_e (float) – <m> Electron beam radius.

  • v_ax (float) – <V> Axial barrier bias.

  • b_ax (float) – <T> Axial magnetic flux density in the trap.

  • r_dt (float) – <m> Drift tube radius.

  • length (Optional[float]) – <m> Trap length -> Currently not used in simulations.

  • v_ra (Optional[float]) – <V> Override for radial trap depth. Only effective if ModelOptions.RADIAL_DYNAMICS=False.

  • j (Optional[float]) – <A/cm^2> Override for current density.

  • fwhm (Optional[float]) – <eV> Override for the electron beam energy spread. Only effective if ModelOptions.RADIAL_DYNAMICS=False.

  • n_grid (int) – Approximate number of nodes for the radial mesh.

  • r_dt_bar (Optional[float]) – Radius of the barrier drift tubes. If not passed, assuming equal radius as trap drift tube.

Returns

ebisim.simulation.Device – The populated device object.

Return type

Device

index(value, start=0, stop=9223372036854775807, /)

Return first index of value.

Raises ValueError if the value is not present.

property b_ax

<T> Axial magnetic field density.

property current

<A> Beam current.

property e_kin

<eV> Uncorrected beam energy.

property fwhm

<eV> Electron beam energy spread.

property j

<A/cm^2> Current density.

property length

<m> Trap length.

property r_dt

<m> Drift tube radius.

property r_dt_bar

<m> Drift tube radius of the barrier drift tubes.

property r_e

<m> Beam radius.

property rad_fd_d

Diagonal vector of finite difference scheme.

property rad_fd_l

Lower diagonal vector of finite difference scheme.

property rad_fd_u

Upper diagonal vector of finite difference scheme.

property rad_grid

<m> Radial grid for finite difference computations.

property rad_phi_ax_barr

<V> Radial potential of the electron beam in the barrier tube

property rad_phi_uncomp

<V> Radial potential of the electron beam.

property rad_re_idx

Index of the radial grid point closest to r_e

property v_ax

<V> Axial barrier voltage.

property v_ax_sc

<V> Axial barrier space charge correction.

property v_ra

<V> Radial barrier voltage.

class ebisim.simulation.EnergyScanResult(sim_kwargs, energies, results)[source]

Bases: object

This class provides a convenient interface to access and evaluate the the results of the ebisim.simulation.energy_scan function. Abundances at arbitrary times are provided by performing linear interpolations of the solutions of the rate equation.

Parameters
  • sim_kwargs (dict) – The sim_kwargs dictionary as provided during the call to ebisim.simulation.energy_scan.

  • energies (numpy.array) – <eV> A sorted array containing the energies at which the energy scan has been evaluated.

  • results (list of ebisim.simulation.Result) – A list of Result objects holding the results of each individual simulation

abundance_at_time(t)[source]

Provides information about the charge state distribution at a given time for all energies.

Parameters

t (float) – <s> Point of time to evaluate.

Returns

  • energies (numpy.array) – The evaluated energies.

  • abundance (numpy.array) – Contains the abundance of each charge state (rows) for each energy (columns).

Raises

ValueError – If ‘t’ is not part of the simulated time domain

abundance_of_cs(cs)[source]

Provides information about the abundance of a single charge states at all simulated times and energies.

Parameters

cs (int) – The charge state to evaluate.

Returns

  • energies (numpy.array) – The evaluated energies.

  • times (numpy.array) – The evaluated timesteps.

  • abundance (numpy.array) – Abundance of charge state ‘cs’ at given times (rows) and energies (columns).

Raises

ValueError – If ‘cs’ is not a sensible charge state for the performed simulation.

get_result(e_kin)[source]

Returns the result object corresponding to the simulation at a given energy

Parameters

e_kin (float) – <eV> The energy of the simulation one wishes to retrieve.

Returns

ebisim.simulation.Result – The Result object for the polled scan step.

Raises

ValueError – If the polled energy is not available.

plot_abundance_at_time(t, cs=None, **kwargs)[source]

Produces a plot of the charge state abundance for different energies at a given time.

Parameters
  • t (float) – <s> Point of time to evaluate.

  • cs (list or None, optional) – If None, all charge states are plotted. By supplying a list of int it is possible to filter the charge states that should be plotted. By default None.

  • **kwargs – Keyword arguments are handed down to ebisim.plotting.plot_energy_scan, cf. documentation thereof. If no arguments are provided, reasonable default values are injected.

Returns

matplotlib.Figure – Figure handle of the generated plot.

plot_abundance_of_cs(cs, **kwargs)[source]

Produces a 2D contour plot of the charge state abundance for all simulated energies and times.

Parameters
  • cs (int) – The charge state to plot.

  • **kwargs – Keyword arguments are handed down to ebisim.plotting.plot_energy_scan, cf. documentation thereof. If no arguments are provided, reasonable default values are injected.

Returns

matplotlib.Figure – Figure handle of the generated plot.

class ebisim.simulation.ModelOptions(EI=True, RR=True, CX=True, DR=False, SPITZER_HEATING=True, COLLISIONAL_THERMALISATION=True, ESCAPE_AXIAL=True, ESCAPE_RADIAL=True, RECOMPUTE_CROSS_SECTIONS=False, RADIAL_DYNAMICS=False, IONISATION_HEATING=True, OVERRIDE_FWHM=False, RADIAL_SOLVER_MAX_STEPS=500, RADIAL_SOLVER_REL_DIFF=0.001)[source]

Bases: tuple

An instance of ModelOptions can be used to turn on or off certain effects in an advanced simulation.

Create new instance of ModelOptions(EI, RR, CX, DR, SPITZER_HEATING, COLLISIONAL_THERMALISATION, ESCAPE_AXIAL, ESCAPE_RADIAL, RECOMPUTE_CROSS_SECTIONS, RADIAL_DYNAMICS, IONISATION_HEATING, OVERRIDE_FWHM, RADIAL_SOLVER_MAX_STEPS, RADIAL_SOLVER_REL_DIFF)

Parameters
  • EI (bool) –

  • RR (bool) –

  • CX (bool) –

  • DR (bool) –

  • SPITZER_HEATING (bool) –

  • COLLISIONAL_THERMALISATION (bool) –

  • ESCAPE_AXIAL (bool) –

  • ESCAPE_RADIAL (bool) –

  • RECOMPUTE_CROSS_SECTIONS (bool) –

  • RADIAL_DYNAMICS (bool) –

  • IONISATION_HEATING (bool) –

  • OVERRIDE_FWHM (bool) –

  • RADIAL_SOLVER_MAX_STEPS (int) –

  • RADIAL_SOLVER_REL_DIFF (float) –

count(value, /)

Return number of occurrences of value.

index(value, start=0, stop=9223372036854775807, /)

Return first index of value.

Raises ValueError if the value is not present.

property COLLISIONAL_THERMALISATION

Switch for ion-ion thermalisation, default True.

property CX

Switch for charge exchange, default True.

property DR

Switch for dielectronic recombination, default False.

property EI

Switch for electron impact ionisation, default True.

property ESCAPE_AXIAL

Switch for axial escape from the trap, default True.

property ESCAPE_RADIAL

Switch for radial escape from the trap, default True.

property IONISATION_HEATING

Switch for ionisation heating/recombination cooling

property OVERRIDE_FWHM

If set use FWHM from device definition instead of computed value.

property RADIAL_DYNAMICS

Switch for effects of radial ion cloud extent. May be computationally very intensive

property RADIAL_SOLVER_MAX_STEPS

Alias for field number 12

property RADIAL_SOLVER_REL_DIFF

Alias for field number 13

property RECOMPUTE_CROSS_SECTIONS

Switch deciding whether EI, RR, and DR cross sections are recomputed on each call of the differential equation system. Advisable if electron beam energy changes over time and sharp transitions are expected, e.g. DR or ionisation thresholds for a given shell. Default False.

property RR

Switch for radiative recombination, default True.

property SPITZER_HEATING

Switch for Spitzer- or electron-heating, default True.

class ebisim.simulation.Rate(value)[source]

Bases: IntEnum

Enum for conveniently identifying rates produced in advanced simulations

AX_CO = 105
CHARGE_EXCHANGE = 104
COLLISIONAL_THERMALISATION = 301
COLLISION_RATE_SELF = 402
COLLISION_RATE_TOTAL = 401
CX = 104
DIELECTRONIC_RECOMBINATION = 103
DR = 103
EI = 101
ELECTRON_IONISATION = 101
E_KIN_FWHM = 532
E_KIN_MEAN = 531
F_EI = 501
IONISATION_HEAT = 303
LOSSES_AXIAL_COLLISIONAL = 105
LOSSES_RADIAL_COLLISIONAL = 107
OVERLAP_FACTORS_EBEAM = 501
RADIATIVE_RECOMBINATION = 102
RA_CO = 107
RR = 102
SPITZER_HEATING = 302
TRAPPING_PARAMETER_AXIAL = 521
TRAPPING_PARAMETER_RADIAL = 522
TRAP_DEPTH_AXIAL = 523
TRAP_DEPTH_RADIAL = 524
T_AX_CO = 205
T_COLLISIONAL_THERMALISATION = 301
T_CT = 301
T_LOSSES_AXIAL_COLLISIONAL = 205
T_LOSSES_RADIAL_COLLISIONAL = 207
T_RA_CO = 207
T_SH = 302
T_SPITZER_HEATING = 302
V_AX = 523
V_RA = 524
W_AX = 521
W_RA = 522
ebisim.simulation.advanced_simulation(device, targets, t_max, bg_gases=None, options=None, rates=False, solver_kwargs=None, verbose=True, n_threads=1)[source]

Interface for performing advanced charge breeding simulations.

For a list of effects refer to ebisim.simulation.ModelOptions.

Parameters
  • device (Device) – Container describing the EBIS/T and specifically the electron beam.

  • targets (Union[Element, List[Element]]) – Target(s) for which charge breeding is simulated.

  • t_max (float) – <s> Simulated breeding time

  • bg_gases (Optional[Union[BackgroundGas, List[BackgroundGas]]]) – Background gas(es) which act as CX partners.

  • rates (bool) – If true a ‘second run’ is performed to store the rates, this takes extra time and can create quite a bit of data.

  • options (Optional[ModelOptions]) – Switches for effects considered in the simulation, see default values of ebisim.simulation.ModelOptions.

  • solver_kwargs (Optional[Dict[str, Any]]) – If supplied these keyword arguments are unpacked in the solver call. Refer to the documentation of scipy.integrate.solve_ivp for more information. By default None.

  • verbose (bool) – Print a little progress indicator and some status messages, by default True.

  • n_threads (int) – How many threads to use (mostly for jacbion estimation which can evaluate the RHS in parallel with different inputs.)

Returns

  • An instance of the AdvancedResult class, holding the simulation parameters, timesteps and

  • charge state distribution including the species temperature.

Return type

Union[AdvancedResult, Tuple[AdvancedResult, …]]

ebisim.simulation.basic_simulation(element, j, e_kin, t_max, dr_fwhm=None, N_initial=None, CNI=False, solver_kwargs=None)[source]

Interface for performing basic charge breeding simulations.

These simulations only include the most important effects, i.e. electron ionisation, radiative recombination and optionally dielectronic recombination (for those transitions whose data is available in the resource directory). All other effects are ignored.

Continuous Neutral Injection (CNI) can be activated on demand.

The results only represent the proportions of different charge states, not actual densities.

Parameters
  • element (Union[Element, str, int]) – An instance of the Element class, or an identifier for the element, i.e. either its name, symbol or proton number.

  • j (float) – <A/cm^2> Current density

  • e_kin (float) – <eV> Electron beam energy

  • t_max (float) – <s> Simulated breeding time

  • dr_fwhm (Optional[float]) – <eV> If a value is given, determines the energy spread of the electron beam (in terms of Full Width Half Max) and hence the effective width of DR resonances. Otherwise DR is excluded from the simulation.

  • N_initial (Optional[ndarray]) – Determines the initial charge state distribution if given, must have Z + 1 entries, where the array index corresponds to the charge state. If no value is given the distribution defaults to 100% of 1+ ions at t = 0 s, or a small amount of neutral atoms in the case of CNI.

  • CNI (bool) – If Continuous Neutral Injection is activated, the neutrals are assumed to form an infinite reservoir. Their absolute number will not change over time and hence they act as a constant source of new singly charged ions. Therefore the absolute amount of ions increases over time.

  • solver_kwargs (Optional[Dict[str, Any]]) – If supplied these keyword arguments are unpacked in the solver call. Refer to the documentation of scipy.integrate.solve_ivp for more information.

Returns

  • An instance of the BasicResult class, holding the simulation parameters, timesteps and

  • charge state distribution.

Return type

BasicResult

@numba.jitebisim.simulation.boltzmann_radial_potential_linear_density(r, rho_0, nl, kT, q, first_guess=None, ldu=None)[source]

Solves the Boltzmann Poisson equation for a static background charge density rho_0 and particles with line number density n, Temperature kT and charge state q.

Below, nRS and nCS are the number of radial sampling points and charge states.

Solution is found through Newton Raphson iterations, cf. [PICNPSb].

Parameters
  • r (np.ndarray) – <m> Radial grid points, with r[0] = 0, r[-1] = r_max.

  • rho_0 (np.ndarray (nRS, )) – <C/m^3> Static charge density at r.

  • nl (np.ndarray (1, nCS)) – <1/m> Line number density of Boltzmann distributed particles.

  • kT (np.ndarray (1, nCS)) – <eV> Temperature of Boltzmann distributed particles.

  • q (np.ndarray (1, nCS)) – Charge state of Boltzmann distributed particles.

  • ldu ((np.ndarray, np.ndarray, np.ndarray)) – The lower diagonal, diagonal, and upper diagonal vector describing the finite difference scheme. Can be provided if they have been pre-computed.

Returns

  • phi (np.ndarray (nRS, )) – <V> Potential at r.

  • nax (np.ndarray (1, nCS)) – <1/m^3> On axis number densities.

  • shape (np.ndarray (nRS, nCS)) – Radial shape factor of the particle distributions.

References

PICNPSb

“Nonlinear Poisson Solver” https://www.particleincell.com/2012/nonlinear-poisson-solver/

@numba.jitebisim.simulation.boltzmann_radial_potential_linear_density_ebeam(r, current, r_e, e_kin, nl, kT, q, first_guess=None, ldu=None, max_step=500, rel_diff=0.001)[source]

Solves the Boltzmann Poisson equation for a static background charge density rho_0 and particles with line number density n, Temperature kT and charge state q. The electron beam charge density is computed from a uniform current density and the iteratively corrected velocity profile of the electron beam.

Below, nRS and nCS are the number of radial sampling points and charge states.

Solution is found through Newton Raphson iterations, cf. [PICNPS].

Parameters
  • r (np.ndarray) – <m> Radial grid points, with r[0] = 0, r[-1] = r_max.

  • current (float) – <A> Electron beam current (positive sign).

  • r_e (float) – <m> Electron beam radius.

  • e_kin (float) – <eV> Uncorrected electron beam energy.

  • nl (np.ndarray (1, nCS)) – <1/m> Line number density of Boltzmann distributed particles.

  • kT (np.ndarray (1, nCS)) – <eV> Temperature of Boltzmann distributed particles.

  • q (np.ndarray (1, nCS)) – Charge state of Boltzmann distributed particles.

  • ldu ((np.ndarray, np.ndarray, np.ndarray)) – The lower diagonal, diagonal, and upper diagonal vector describing the finite difference scheme. Can be provided if they have been pre-computed.

Returns

  • phi (np.ndarray (nRS, )) – <V> Potential at r.

  • nax (np.ndarray (1, nCS)) – <1/m^3> On axis number densities.

  • shape (np.ndarray (nRS, nCS)) – Radial shape factor of the particle distributions.

References

PICNPS

“Nonlinear Poisson Solver” https://www.particleincell.com/2012/nonlinear-poisson-solver/

@numba.jitebisim.simulation.boltzmann_radial_potential_onaxis_density(r, rho_0, n, kT, q, first_guess=None, ldu=None)[source]

Solves the Boltzmann Poisson equation for a static background charge density rho_0 and particles with a fixed on axis density n, Temperature kT and charge state q.

Below, nRS and nCS are the number of radial sampling points and charge states.

Solution is found through Newton iterations, cf. [PICNPSa].

Parameters
  • r (np.ndarray) – <m> Radial grid points, with r[0] = 0, r[-1] = r_max.

  • rho_0 (np.ndarray (nRS, )) – <C/m^3> Static charge density at r.

  • n (np.ndarray (1, nCS)) – <1/m> On axis density of Boltzmann distributed particles.

  • kT (np.ndarray (1, nCS)) – <eV> Temperature of Boltzmann distributed particles.

  • q (np.ndarray (1, nCS)) – Charge state of Boltzmann distributed particles.

  • ldu ((np.ndarray, np.ndarray, np.ndarray)) – The lower diagonal, diagonal, and upper diagonal vector describing the finite difference scheme. Can be provided if they have been pre-computed.

Returns

  • phi (np.ndarray (nRS, )) – <V> Potential at r.

  • nax (np.ndarray (1, nCS)) – <1/m^3> On axis number densities.

  • shape (np.ndarray (nRS, nCS)) – Radial shape factor of the particle distributions.

References

PICNPSa

“Nonlinear Poisson Solver” https://www.particleincell.com/2012/nonlinear-poisson-solver/

ebisim.simulation.energy_scan(sim_func, sim_kwargs, energies, parallel=False)[source]

This function provides a convenient way to repeat the same simulation for a number of different electron beam energies. This can reveal variations in the charge state balance due to weakly energy dependent ionisation cross sections or even resonant phenomena like dielectronic recombination.

Parameters
  • sim_func (callable) – The function handle for the simulation e.g. ebisim.simulation.basic_simulation.

  • sim_kwargs (dict) – A dictionary containing all the required and optional parameters of the simulations (except for the kinetic electron energy) as key value pairs. This is unpacked in the function call.

  • energies (list or numpy.array) – A list or array of the energies at which the simulation should be performed.

  • parallel (bool, optional) – Determine whether multiple simulations should be run in parallel using pythons multiprocessing.pool. This may accelerate the scan when performing a large number of simulations. By default False.

Returns

ebisim.simulation.EnergyScanResult – An object providing convenient access to the generated scan data.

@numba.jitebisim.simulation.fd_system_nonuniform_grid(r)[source]

Sets up the three diagonal vectors for a finite Poisson problem with radial symmetry on a nonuniform grid. d phi/dr = 0 at r = 0, and phi = phi0 at r = (n-1) * dr = r_max The finite differences are developed according to [Sundqvist1970].

Parameters

r (np.ndarray) – <m> Radial grid points, with r[0] = 0, r[-1] = r_max.

Returns

  • l (np.ndarray) – Lower diagonal vector.

  • d (np.ndarray) – Diagonal vector.

  • u (np.ndarray) – Upper diagonal vector.

References

Sundqvist1970

“A simple finite-difference grid with non-constant intervals”, Sundqvist, H., & Veronis, G., Tellus, 22(1), 26–31 (1970), https://doi.org/10.3402/tellusa.v22i1.10155

See also

ebisim.simulation._radial_dist.fd_system_uniform_grid

@numba.jitebisim.simulation.fd_system_uniform_grid(r)[source]

Sets up the three diagonal vectors for a finite Poisson problem with radial symmetry on a uniform grid. d phi/dr = 0 at r = 0, and phi = phi0 at r = (n-1) * dr = r_max

Parameters

r (np.ndarray) – <m> Radial grid points, must be evenly spaced, with r[0] = 0, r[-1] = r_max.

Returns

  • l (np.ndarray) – Lower diagonal vector.

  • d (np.ndarray) – Diagonal vector.

  • u (np.ndarray) – Upper diagonal vector.

See also

ebisim.simulation._radial_dist.fd_system_nonuniform_grid

@numba.jitebisim.simulation.heat_capacity(r, phi, q, kT)[source]

Computes the heat capacity of an ion cloud with a given charge state and temperature inside the external potenital phi(r). According to Lu, Currell cloud expansion.

Parameters
  • r (np.ndarray) – <m> Radial grid points, with r[0] = 0, r[-1] = r_max.

  • phi (np.ndarray) – <V> Potential at r.

  • q (int) – Ion charge state.

  • kT (float) – <eV> Ion temperature

Returns

float – <eV/eV> Constant volume heat capacity

@numba.jitebisim.simulation.radial_potential_nonuniform_grid(r, rho)[source]

Solves the radial Poisson equation on a nonuniform grid. Boundary conditions are d phi/dr = 0 at r = 0 and phi(rmax) = 0

Parameters
  • r (np.ndarray) – <m> Radial grid points, with r[0] = 0, r[-1] = r_max.

  • rho (np.ndarray) – <C/m^3> Charge density at r.

Returns

np.ndarray – Potential at r.

@numba.jitebisim.simulation.radial_potential_uniform_grid(r, rho)[source]

Solves the radial Poisson equation on a uniform grid. Boundary conditions are d phi/dr = 0 at r = 0 and phi(rmax) = 0

Parameters
  • r (np.ndarray) – <m> Radial grid points, must be evenly spaced, with r[0] = 0, r[-1] = r_max.

  • rho (np.ndarray) – <C/m^3> Charge density at r.

Returns

np.ndarray – <V> Potential at r.

@numba.jitebisim.simulation.tridiagonal_matrix_algorithm(l, d, u, b)[source]

Tridiagonal Matrix Algorithm [TDMA]. Solves a system of equations M x = b for x, where M is a tridiagonal matrix.

M = np.diag(d) + np.diag(u[:-1], 1) + np.diag(l[1:], -1)

Parameters
  • l (np.ndarray) – Lower diagonal vector l[i] = M[i, i-1].

  • d (np.ndarray) – Diagonal vector d[i] = M[i, i].

  • u (np.ndarray) – Upper diagonal vector u[i] = M[i, i+1].

  • b (np.ndarray) – Inhomogenety term.

Returns

x (np.ndarray) – Solution of the linear system.

References

TDMA

“Tridiagonal matrix algorithm” https://en.wikipedia.org/wiki/Tridiagonal_matrix_algorithm