Using XrayDB from Python

The python directory contains the source code for a Python module for XrayDB. This module gives a user-friendly wrapping of the XrayDB, and automates the the conversion of data from sqlite database into Python and numpy arrays. The module requires the numpy, scipy and sqlalchemy modules, all of which are readily available and can be installed with:

pip install xraydb

The current version of the Python module is 4.5.4, corresponding to version 6 of xraydb.sqlite.

The Python xraydb module

To use the XrayDB from Python, you can import the xraydb module and start using it:

>>> import xraydb
>>> xraydb.atomic_number('Ag')
47
#
# X-ray elastic (Thomson) scattering factors:
>>> import numpy as np
>>> q =np.linspace(0, 5, 11)
>>> xraydb.f0('Fe', q)
array([25.994603  , 11.50848469,  6.55945765,  4.71039413,  3.21048827,
       2.20939146,  1.65112769,  1.36705887,  1.21133507,  1.10155689,
       1.0035555 ])
#
# X-ray emission lines:
>>> for name, line in xraydb.xray_lines('Zn', 'K').items():
...     print(name, ' = ', line)
...
Ka3  =  XrayLine(energy=8462.8, intensity=0.000316256, initial_level='K', final_level='L1')
Ka2  =  XrayLine(energy=8614.1, intensity=0.294353, initial_level='K', final_level='L2')
Ka1  =  XrayLine(energy=8637.2, intensity=0.576058, initial_level='K', final_level='L3')
Kb3  =  XrayLine(energy=9567.6, intensity=0.0438347, initial_level='K', final_level='M2')
Kb1  =  XrayLine(energy=9570.4, intensity=0.0846229, initial_level='K', final_level='M3')
Kb5  =  XrayLine(energy=9648.8, intensity=0.000815698, initial_level='K', final_level='M4,5')
#
# X-ray absorption edges:
>>> xraydb.xray_edge('As', 'K')
XrayEdge(energy=11867.0, fyield=0.548989, jump_ratio=7.314)
#
# X-ray attenuation factors:
>>> as_kedge = xraydb.xray_edge('As', 'K').energy
>>> energies = np.linspace(-50, 50, 5) + as_kedge
>>> muvals = xraydb.mu_elam('As', energies)
>>> for en, mu in zip(energies, muvals):
...     print("{:.0f}   {:8.2f}".format(en, mu))
...
11817      26.07
11842      25.92
11867      25.77
11892     178.32
11917     177.38

Table of XrayDB function for Atomic and X-ray data for the elements

Most of these function return some element-specific property from the element symbol or atomic number. Some of the data extends to Z=98 (Cf), but some data may not be available for Z > 92 (U). Except where noted, the data comes from [Elam, Ravel, and Sieber (2002)].

xraydb functions

description

atomic_number()

atomic number from symbol

atomic_symbol()

atomic symbol from number

atomic_mass()

atomic mass

atomic_name()

atomic name (English)

atomic_density()

density of pure element

f0()

elastic scattering factor ([Waasmaier and Kirfel (1995)])

f0_ions()

list of valid “ions” for f0() ([Waasmaier and Kirfel (1995)])

xray_edge()

xray edge data for a particular element and edge

xray_edges()

dictionary of all X-ray edges data for an element

xray_lines()

dictionary of all X-ray emission line data for an element

fluor_yield()

fluorescent yield for an X-ray emission line

ck_probability()

Coster-Kronig transition probability between two atomic levels

mu_elam()

absorption cross-section, photo-electric or total for an element

coherent_cross_section_elam()

coherent scattering cross-section for an element

incoherent_cross_section_elam()

incoherent scattering cross-section for an element

chantler_energies()

energies of tabulation for Chantler data ([Chantler (2000)])

f1_chantler()

\(f'(E)\) anomalous scattering factor ([Chantler (2000)])

f2_chantler()

\(f"(E)\) anomalous scattering factor ([Chantler (2000)])

mu_chantler()

absorption cross-section ([Chantler (2000)])

guess_edge()

guess element and edge from energy of absorption edge

chemparse()

parse a chemical formula to atomic abundances

validate_formula()

test whether a chemical formula can be parsed.

get_materials()

get a dictionary of known materials {name:(formula, density)}

get_material()

get a (formula, density) tuple for a material in the materials database

find_material()

get a material instance for a material in the materials database

add_material()

add a material to local materials database

material_mu()

absorption cross-section for a material at X-ray energies

material_mu_components()

dictionary of elemental components of mu for material

xray_delta_beta()

anomalous index of refraction for material and energy

darwin_width()

Darwin widths for monochromator crystals

mirror_reflectivity()

X-ray reflectivities for mirror materials (thick slab limit)

ionization_potential()

effective ionization potential for a gas, as for ion chambers

ionchamber_fluxes()

calculate fluxes from ion chamber voltages, gases, and sensitivities

get_xraydb()

return instance of the XrayDB

Returns:

XrayDB

Example

>>> import xraydb
>>> xdb = xraydb.get_xraydb()

Atomic Properties

atomic_number(element)

z for element name

Parameters:

element (str) – atomic symbol

Returns:

atomic number

atomic_symbol(z)

atomic symbol for atomic number

Parameters:

z (int) – atomic number

Returns:

atomic symbol

atomic_mass(element)

molar mass for an element

Parameters:

element (int, str) – atomic number, atomic symbol for element

Returns:

atomic mass, in AMU

atomic_name(z)

atomic name for atomic number

Parameters:

z (int) – atomic number

Returns:

atomic name (English)

atomic_density(element)

density (gr/cm^3) for common for of an element

Parameters:

element (int, str) – atomic number, atomic symbol for element

Returns:

density in gm/cm^3

Elastic Scattering Factors

f0(ion, q)

elastic X-ray scattering factor, f0(q), for an ion.

Parameters:
  • ion (int or str) – atomic number, atomic symbol or ionic symbol of scatterer

  • q (float, ndarray) – Q value(s) for scattering

Returns

scattering factor for each Q value

Notes

  1. from D. Waasmaier and A. Kirfel, Acta Cryst. A51 p416 (1995) and International Tables for Crystallography, Vol. C.

  2. ion can be of the form: 26, Fe, Fe2+. For a full list of ions use f0_ions()

  3. elements supported are from Z = 1 to 98 (‘H’ to ‘Cf’)

  4. q = sin(theta) / lambda, where theta=incident angle, lambda=X-ray wavelength

f0_ions(element=None)

list ion names supported in the f0() calculation from Waasmaier and Kirfel.

Parameters:

element (None, int, str) – scatterer

Returns:

list of strings for matching ion names

Notes

if element is None, all 211 ions are returned.

X-ray Edges

xray_edge(element, edge, energy_only=False)

get x-ray absorption edge data for an element: (energy(in eV), fluorescence yield, jump ratio)

Parameters:
  • element (int, str) – atomic number, atomic symbol for element

  • edge (str) – iupac symbol of X-ray edge

  • energy_only (bool) – whether to return only the energy [False]

Returns:

XrayEdge namedtuple containing (energy, fluorescence_yield, edge_jump) or float of energy

xray_edges(element)
get dictionary of x-ray absorption edges:

energy(in eV), fluorescence yield, and jump ratio for an element.

Parameters:

element (int, str) – atomic number, atomic symbol for element

Returns:

dictionary of XrayEdge named tuples.

Notes

  1. The dictionary will have keys of edge (iupac symbol) and values containing an XrayEdge namedtuple containing (energy, fluorescence_yield, edge_jump)

core_width(element, edge=None)

returns core hole width for an element and edge

Parameters:
  • element (int or str) – element

  • edge (None or str) – edge to consider

Returns:

core width or list of core widths

Notes

  1. if edge is None, values are return for all edges

  2. Data from Krause and Oliver (1979) and Keski-Rahkonen and Krause (1974)

guess_edge(energy, edges=['K', 'L3', 'L2', 'L1', 'M5'])

guess an element and edge based on energy (in eV)

Parameters:
  • energy (float) – approximate edge energy (in eV)

  • edges (None or list of strings) – edges to consider

Returns:

a tuple of (atomic symbol, edge) for best guess

Notes

by default, the list of edges is [‘K’, ‘L3’, ‘L2’, ‘L1’, ‘M5’]

X-ray Emission Lines

xray_lines(element, initial_level=None, excitation_energy=None)

get dictionary of X-ray emission lines of an element

Parameters:
  • element (int, str) – atomic number, atomic symbol for element

  • initial_level (None or str) – iupac symbol of initial level

  • excitation_energy (None or float) – exciation energy

Returns:

dict of X-ray lines with keys of siegbahn notation and values of XrayLine tuples of (energy, intensity, initial level, final level)

Notes

  1. excitation energy will supercede initial_level, as it means ‘all intial levels with below this energy

Exaample:
>>> for name, line in xraydb.xray_lines('Mn', 'K').items():
...     print(name, line)
...
Ka3 XrayLine(energy=5769.9, intensity=0.000265963, initial_level='K', final_level='L1')
Ka2 XrayLine(energy=5889.1, intensity=0.293941, initial_level='K', final_level='L2')
Ka1 XrayLine(energy=5900.3, intensity=0.58134, initial_level='K', final_level='L3')
Kb3 XrayLine(energy=6491.8, intensity=0.042234, initial_level='K', final_level='M2')
Kb1 XrayLine(energy=6491.8, intensity=0.0815329, initial_level='K', final_level='M3')
Kb5 XrayLine(energy=6537.0, intensity=0.000685981, initial_level='K', final_level='M4,5')
fluor_yield(element, edge, line, energy)

fluorescence yield for an X-ray emission line or family of lines.

Parameters:
  • element (int, str) – atomic number, atomic symbol for element

  • edge (str) – iupac symbol of X-ray edge

  • line (str) – siegbahn notation for emission line

  • energy (float) – incident X-ray energy

Returns:

fluorescence yield, weighted average fluorescence energy, net_probability

Examples

>>> xraydb.fluor_yield('Fe', 'K', 'Ka', 8000)
0.350985, 6400.752419799043, 0.874576096
>>> xraydb.fluor_yield('Fe', 'K', 'Ka', 6800)
0.0, 6400.752419799043, 0.874576096
>>> xraydb.fluor_yield('Ag', 'L3', 'La', 6000)
0.052, 2982.129655446868, 0.861899000000000

See also

xray_lines which gives the full set of emission lines (‘Ka1’, ‘Kb3’, etc) and probabilities for each of these.

ck_probability(element, initial, final, total=True)

transition probability for an element, initial, and final levels.

Parameters:
  • element (int, str) – atomic number, atomic symbol for element

  • initial (str) – iupac symbol for initial level

  • final (str) – iupac symbol for final level

  • total (bool) – whether to include transitions via possible intermediate levels [True]

Returns:

transition probability, or 0 if transition is not allowed.

Absorption and Scattering Cross-sections

mu_elam(element, energy, kind='total')

X-ray mass attenuation coefficient, mu/rho, for an element and energy or array of energies. Data is from the Elam tables.

Parameters:
  • element (int, str) – atomic number, atomic symbol for element

  • energy (float or ndarray) – energy or array of energies

  • kind (str) – type of cross-section to use, one of (‘total’, ‘photo’, ‘coh’, ‘incoh’) [‘total’]

Returns:

float value or ndarray

Notes

  1. Values returned are in units of cm^2/gr

  2. The default is to return total attenuation coefficient.

coherent_cross_section_elam(element, energy)

coherent scaattering cross-section for an element and energy or array of energies. Data is from the Elam tables.

Parameters:
  • element (int, str) – atomic number, atomic symbol for element

  • energy (float or ndarray) – energy or array of energies

Returns:

float value or ndarray

Notes

  1. Values returned are in units of cm^2/gr

incoherent_cross_section_elam(element, energy)

incoherent scaattering cross-section for an element and energy or array of energies. Data is from the Elam tables.

Parameters:
  • element (int, str) – atomic number, atomic symbol for element

  • energy (float or ndarray) – energy or array of energies

Returns:

float value or ndarray

Notes

  1. Values returned are in units of cm^2/gr

chantler_energies(element, emin=0, emax=1000000000.0)

energies at which Chantler data is tabulated for a particular element.

Parameters:
  • element (int, str) – atomic number, atomic symbol for element

  • emin (float) – lower bound of energies (default=0)

  • emax (float) – upper bound of energies (default=1.e9)

Returns:

ndarray of energies

Notes

energies are in eV

f1_chantler(element, energy, **kws)

real part of anomalous x-ray scattering factor for an element and energy or array of energies. Data is from the Chantler tables.

Parameters:
  • element (int, str) – atomic number, atomic symbol for element

  • energy (float or ndarray) – energy or array of energies

Returns:

float value or ndarray

Notes

  1. Values returned are in units of electrons

f2_chantler(element, energy)

imaginary part of anomalous x-ray scattering factor for an element and energy or array of energies. Data is from the Chantler tables.

Parameters:
  • element (int, str) – atomic number, atomic symbol for element

  • energy (float or ndarray) – energy or array of energies

Returns:

float value or ndarray

Notes

  1. Values returned are in units of electrons

mu_chantler(element, energy, incoh=False, photo=False)

X-ray mass attenuation coeficient, mu/rho, for an element and energy or array of energies. Data is from the Chantler tables.

Parameters:
  • element (int, str) – atomic number, atomic symbol for element

  • energy (float or ndarray) – energy or array of energies

  • incoh (bool) – whether to return only the incoherent contribution [False]

  • photo (bool) – whether to return only the photo-electric contribution [False]

Returns:

float value or ndarray

Notes

  1. Values returned are in units of cm^2/gr

  2. The default is to return total attenuation coefficient.

Chemical and Materials database

chemparse(formula)

parse a chemical formula to a dictionary of elemental abundances

Parameters:

formula (str) – chemical formula

Returns:

dict of element symbol and abundance.

Examples

>>> from xraydb import chemparse
>>> chemparse('Mn(SO4)2(H2O)7)')
{'H': 14.0, 'S': 2.0, 'Mn': 1, 'O': 15.0}
>>> chemparse('Zn1.e-5Fe3O4')
{'Zn': 1e-05, 'Fe': 3.0, 'O': 4.0}
>>> chemparse('CO')
{'C': 1, 'O': 1}
>>> chemparse('Co')
{'Co': 1}
>>> chemparse('co')
ValueError: unrecognized element or number:
co
validate_formula(formula)

return whether a chemical formula is valid and can be parsed to a dictionary with chemparse()

Parameters:

formula (str) – chemical formula

Returns:

bool (True or False) for whether chemparse() will succeed

Examples

>>> from xraydb import validate_formula
>>> validate_formula('Mn(SO4)2(H2O)7)')
True
>>> validate_formula('Mn(SO4)2(H2O7')
False
>>> validate_formula('Z')
False
get_materials(force_read=False, categories=None)

get dictionary of all available materials

Parameters:
  • force_read (bool) – whether to force a re-reading of the materials database [False]

  • categories (list of strings or None) – restrict results to those that match category names

Returns:

dict with keys of material name and values of Materials instances

Examples

>>> for name, m in xraydb.get_materials().items():
...      print(name, m)
...
water H2O 1.0
lead Pb 11.34
aluminum Al 2.7
kapton C22H10N2O5 1.42
polyimide C22H10N2O5 1.42
nitrogen N 0.00125
argon Ar 0.001784
...
find_material(name)

look up material name, return material instance

Parameters:

name (str) – name of material or chemical formula

Returns:

material instance

Examples

>>> xraydb.find_material('kapton')
Material(formula='C22H10N2O5', density=1.42, name='kapton', categories=['polymer'])

See also

get_material()

get_material(name)

look up material name, return formula and density

Parameters:

name (str) – name of material or chemical formula

Returns:

chemical formula, density of material

Examples

>>> xraydb.get_material('kapton')
('C22H10N2O5', 1.43)

See also

find_material()

add_material(name, formula, density, categories=None)

add a material to the users local material database

Parameters:
  • name (str) – name of material

  • formula (str) – chemical formula

  • density (float) – density

  • categories (list of strings or None) – list of category names

Returns:

None

Notes

the data will be saved to the file ‘xraydb/materials.dat’ in the users configuration folder, and will be useful in subsequent sessions.

Examples

>>> xraydb.add_material('becopper', 'Cu0.98e0.02', 8.3, categories=['metal'])

X-ray properties of materials

For some further examples, see Example Calculations of X-ray properties of materials.

material_mu(name, energy, density=None, kind='total')

X-ray attenuation length (in 1/cm) for a material by name or formula

Parameters:
  • name (str) – chemical formul or name of material from materials list.

  • energy (float or ndarray) – energy or array of energies in eV

  • density (None or float) – material density (gr/cm^3).

  • kind (str) – ‘photo’ or ‘total’ for whether to return the photo-absorption or total cross-section [‘total’]

Returns:

absorption length in 1/cm

Notes

  1. material names are not case sensitive, chemical compounds are case sensitive.

  2. mu_elam() is used for mu calculation.

  3. if density is None and material is known, that density will be used.

Examples

>>> material_mu('H2O', 10000.0)
5.32986401658495
material_mu_components(name, energy, density=None, kind='total')

material_mu_components: absorption coefficient (in 1/cm) for a compound

Parameters:
  • name (str) – chemical formul or name of material from materials list.

  • energy (float or ndarray) – energy or array of energies in eV

  • density (None or float) – material density (gr/cm^3).

  • kind (str) – ‘photo’ or ‘total’for whether to return photo-absorption or total cross-section [‘total’]

Returns:

dict for constructing mu per element,

with elements ‘mass’ (total mass), ‘density’, and

’elements’ (list of atomic symbols for elements in material).

For each element, there will be an item (atomic symbol as key) with tuple of (stoichiometric fraction, atomic mass, mu)

Examples

>>> xraydb.material_mu('quartz', 10000)
50.36774553547068
>>> xraydb.material_mu_components('quartz', 10000)
{'mass': 60.0843, 'density': 2.65, 'elements': ['Si', 'O'],
'Si': (1, 28.0855, 33.87943243018506), 'O': (2.0, 15.9994, 5.952824815297084)}
xray_delta_beta(material, density, energy)

anomalous components of the index of refraction for a material, using the tabulated scattering components from Chantler.

Parameters:
  • material – chemical formula (‘Fe2O3’, ‘CaMg(CO3)2’, ‘La1.9Sr0.1CuO4’)

  • density – material density in g/cm^3

  • energy – x-ray energy in eV

Returns:

(delta, beta, atlen)

where

delta : real part of index of refraction beta : imag part of index of refraction atlen : attenuation length in cm

These are the anomalous scattering components of the index of refraction:

n = 1 - delta - i*beta = 1 - lambda**2 * r0/(2*pi) Sum_j (n_j * fj)

Adapted from code by Yong Choi

darwin_width(energy, crystal='Si', hkl=(1, 1, 1), a=None, polarization='s', ignore_f2=False, ignore_f1=False, m=1)

darwin width for a crystal reflection and energy

Parameters:
  • energy (float) – X-ray energy in eV

  • crystal (string) – name of crystal (one of ‘Si’, ‘Ge’, or ‘C’) [‘Si’]

  • hkl (tuple) – h, k, l for reflection [(1, 1, 1)]

  • a (float or None) – lattice constant [None - use built-in value]

  • polarization ('s','p', 'u') – mono orientation relative to X-ray polarization [‘s’]

  • ignore_f1 (bool) – ignore contribution from f1 - dispersion (False)

  • ignore_f2 (bool) – ignore contribution from f2 - absorption (False)

  • m (int) – order of reflection [1]

Returns:

A named tuple ‘DarwinWidth’ with the following fields

theta: float, nominal Bragg angle, in rad,

theta_offset: float, angular offset from Bragg angle, in rad,

theta_width: float, estimated angular Darwin width, in rad,

theta_fwhm: float, estimated FWHM of angular intensity, in rad,

energy_width: float, estimated angular Darwin width, in rad,

energy_fwhm: float, estimated FWHM of energy intensity, in eV,

zeta: nd-array of Zeta parameter (delta_Lambda / Lambda),

dtheta: nd-array of angles away from Bragg angle, theta in rad,

denergy: nd-array of energies away from Bragg energy, in eV,

intensity: nd-array of reflected intensity

Notes

  1. This follows the calculation from section 6.4 of Elements of Modern X-ray Physics, 2nd Edition J Als-Nielsen, and D. McMorrow.

  2. Only diamond structures (Si, Ge, diamond) are currently supported. Default values of lattice constant a are in Angstroms: 5.4309 for Si, 5.6578, for ‘Ge’, and 3.567 for ‘C’.

  3. The theta_width and energy_width values will closely match the width of the intensity profile that would = 1 when ignoring the effect of absorption. These are the values commonly reported as ‘Darwin Width’. The value reported for theta_fwhm’ and `energy_fwhm are larger than this by sqrt(9/8) ~= 1.06.

  4. Polarization can be ‘s’, ‘p’, ‘u’, or None. ‘s’ means vertically deflecting crystal and a horizontally-polarized source, as for most synchrotron beamlines. ‘p’ is for a horizontally-deflecting crystal. ‘u’ or None is for unpolarized light, as for most fluorescence/emission.

Examples

>>> dw = darwin_width(10000, crystal='Si', hkl=(1, 1, 1))
>>> dw.theta_width, dw.energy_width
(2.8593683930207114e-05, 1.4177346002236872)
mirror_reflectivity(formula, theta, energy, density=None, roughness=0.0, polarization='s')

mirror reflectivity for a thick, singl-layer mirror.

Parameters:
  • formula (string) – material name or formula (‘Si’, ‘Rh’, ‘silicon’)

  • theta (float or nd-array) – mirror angle in radians

  • energy (float or nd-array) – X-ray energy in eV

  • density (float or None) – material density in g/cm^3

  • roughness (float) – mirror roughness in Angstroms

  • polarization ('s' or 'p') – mirror orientation relative to X-ray polarization

Returns:

mirror reflectivity values

Notes

  1. only one of theta or energy can be an nd-array

  2. density can be None for known materials

  3. polarization of ‘s’ puts the X-ray polarization along the mirror surface, ‘p’ puts it normal to the mirror surface. For horizontally polarized X-ray beams from storage rings, ‘s’ will usually mean ‘vertically deflecting’ and ‘p’ will usually mean ‘horizontally deflecting’.

ionization_potential(gas)

return effective ionization potential for a gas or diode semiconductor, as appropriate for ionization chambers in the linear regime (not in the ‘proportional counter’ regime) or for PIN photodiodes (not in ‘avalanche’ mode).

Parameters:

gas (string) – name of gas or ‘Si’ or ‘Ge’

Returns:

ionization potential in eV

Notes

Data from G. F. Knoll, Radiation Detection and Measurement, Table 5-1, and from ICRU Report 31, 1979. Supported gas names and effective potentials:

gas names

potential (eV)

hydrogen, H

36.5

helium, He

41.3

nitrogen, N, N2

34.8

oxygen, O, O2

30.8

neon, Ne

35.4

argon, Ar

26.4

krypton, Kr

24.4

xenon, Xe

22.1

air

33.8

methane, CH4

27.3

carbondioxide, CO2

33.0

silicon, Si

3.68

germanium, Ge

2.97

If the gas is not recognized the default value of 32.0 eV will be returned.

ionchamber_fluxes(gas='nitrogen', volts=1.0, length=100.0, energy=10000.0, sensitivity=1e-06, sensitivity_units='A/V', with_compton=True, both_carriers=True)

return ion chamber and PIN diode fluxes for a gas, mixture of gases, or semiconductor material, ion chamber length (or diode thickness), X-ray energy, recorded voltage and current amplifier sensitivity. See note for details.

Parameters:
  • gas (string or dict) – name or formula of fill gas (see note 1) [‘nitrogen’]

  • volts (float) – measured voltage output of current amplifier [1.0]

  • length (float) – active length of ion chamber in cm [100]

  • energy (float) – X-ray energy in eV [10000]

  • sensitivity (float) – current amplifier sensitivity [1.e-6]

  • sensitivity_units (string) – units of current amplifier sensitivity (see note 2 for options) [‘A/V’]

  • with_compton (bool) – switch to control the contribution of Compton scattering (see note 3) [True]

  • both_carriers (bool) – switch to control whether to count both electron and ion current (see note 4) [True]

Returns:

named tuple IonchamberFluxes with fields

incident flux of beam incident on ion chamber in Hz

transmitted flux of beam output of ion chamber in Hz

photo flux absorbed by photo-electric effect in Hz

incoherent flux attenuated by incoherent scattering in Hz

Examples

>>> from xraydb import ionchamber_fluxes
>>> fl = ionchamber_fluxes(gas='nitrogen', volts=1.25,
                           length=20.0, energy=10e3, sensitivity=1.e-6)
>>> print(f"Fluxes: In={fl.incident:g}, Out={fl.transmitted:g}, Transmitted={100*fl.transmitted/fl.incident:.2f}%")
Fluxes: In=3.20045e+11, Out=2.90464e+11, Transmitted=90.76%
>>> fl = ionchamber_fluxes(gas={'nitrogen':0.5, 'helium': 0.5},
                           volts=1.25, length=20.0, energy=10000.0,
                           sensitivity=1.e-6)
>>> print(f"Fluxes: In={fl.incident:g}, Out={fl.transmitted:g}, Transmitted={100*fl.transmitted/fl.incident:.2f}%")
Fluxes: In=6.83845e+11, Out=6.51188e+11, Transmitted=95.22%

Notes

  1. The gas value can either be a string for the name of chemical formula for the gas or diode material, or dictionary with keys that are gas names or formulas and values that are the relative fraction for mixes gases. For diode materials, mixtures are not supported.

    The gas formula is used both the contributions for mu and to get the weighted effective ionization potential for the material.

    The effective ionization potentials are known for a handful of gases and diodes (see ionization_potential function), and range between 20 and 45 eV for gases, and around 3 eV for semiconductors. For unknown gases the value of 32.0 eV will be used.

  2. The sensitivity and sensitivity_units arguments have some overlap to specify the sensitivity of the current amplifier. Generally, the units are in A/V, but you can add a common SI prefix of ‘p’, ‘pico’, ‘n’, ‘nano’, (unicode ‘u03bc’), ‘u’, ‘micro’, ‘m’, ‘milli’ so that, ionchamber_fluxes(…, sensitivity=1.e-6) and ionchamber_fluxes(…, sensitivity=1, sensitivity_units=’uA/V’) will both give a sensitivity of 1 microAmp / Volt.

  3. The effect of Compton scattering on the ion chamber current can be approximated using the mean energy of the Compton-scattered electron. See the documentation for more details. Set with_compton=False to turn off this correction.

  4. The effective ionization potential generates an electron and ions pair, and normally both carriers will contribute to the current. Thus, the number of carries below, N_carriers is 2. To consider the current from 1 carrier, for example if using a Frisch grid, use both_carries=False, which will set N_carriers to 1.