Charge Generation models#

Charge generation models are used to add to and manipulate data in Charge array inside the Detector object. The values in the Charge array represent charge in electron. If the photon collection model group is used, a model like Simple photoconversion needs to be enabled in the pipeline to make the conversion from Photon to Charge. Otherwise, a model like Load charge needs to be enabled to initialize the Charge array.

Create and Store a detector#

The models Save detector and Load detector can be used respectively to create and to store a Detector to/from a file.

These models can be used when you want to store or to inject a Detector into the current Pipeline.

Save detector#

This model saves the current Detector into a file. Accepted file formats are .h5, .hdf5, .hdf and .asdf.

- name: save_detector
  func: pyxel.models.save_detector
  enabled: true
  arguments:
    filename: my_detector.h5
pyxel.models.save_detector(detector, filename)[source]

Save the current detector into a file.

Load detector#

This model loads a Detector from a file and injects it in the current pipeline. Accepted file formats are .h5, .hdf5, .hdf and .asdf.

- name: load_detector
  func: pyxel.models.load_detector
  enabled: true
  arguments:
    filename: my_detector.h5
pyxel.models.load_detector(detector, filename)[source]

Load a new detector from a file.

Raises:

TypeError – If the loaded detector has not the same type of the current detector.

Simple photoconversion#

PhotonCharge

With this model you can create and add charge to Detector via photoelectric effect by converting photons to charge. This model supports both monochromatic and multiwavelength photons, converting either a 2D photon array or 3D photon array to the 2D charge array. If the previous model group photon collection returns a 3D photon array, the photon array will be integrated along the wavelength dimension before applying the quantum efficiency (QE).

Binomial sampling of incoming Poisson distributed photons is used in the conversion by default, with probability QE. It can be turned off by setting the argument binomial_sampling to False. User can provide an optional quantum efficiency (quantum_efficiency) parameter. If not provided, quantum efficiency from detector Characteristics is used. It is also possible to set the seed of the random generator with the argument seed.

Basic example of YAML configuration model:

- name: simple_conversion
  func: pyxel.models.charge_generation.simple_conversion
  enabled: true
  arguments:
    quantum_efficiency: 0.8  # optional
pyxel.models.charge_generation.simple_conversion(detector, quantum_efficiency=None, seed=None, binomial_sampling=True)[source]#

Generate charge from incident photon via photoelectric effect, simple model.

Parameters:
  • detector (Detector) – Pyxel Detector object.

  • quantum_efficiency (float, optional) – Quantum efficiency.

  • seed (int, optional)

  • binomial_sampling (bool) – Binomial sampling. Default is True.

Warning

Model assumes shot noise model was applied to photon array when using binomial sampling.

Conversion with custom QE map#

PhotonCharge

With this model you can create and add charge to Detector via photoelectric effect by converting photons in charge. Binomial sampling of incoming Poisson distributed photons is used in the conversion by default, with probability QE. It can be turned off by setting the argument binomial_sampling to False. Besides that, user can input a custom quantum efficiency map by providing a filename of the QE map. Accepted file formats for QE map are .npy, .fits, .txt, .data, .jpg, .jpeg, .bmp, .png and .tiff. Use argument position to set the offset from (0,0) pixel and set where the input QE map is placed onto detector. You can set preset positions with argument align. Values outside of detector shape will be cropped. Read more about placement in the documentation of function fit_into_array().

Basic example of YAML configuration model:

- name: conversion_with_qe_map
  func: pyxel.models.charge_generation.conversion_with_qe_map
  enabled: true
  arguments:
    filename: data/qe_map.npy

Note

You can find an example of this model used in this Jupyter Notebook H2RG pipeline from Pyxel Data.

pyxel.models.charge_generation.conversion_with_qe_map(detector, filename, position=(0, 0), align=None, seed=None, binomial_sampling=True)[source]#

Generate charge from incident photon via photoelectric effect, simple model with custom QE map.

Parameters:
  • detector (Detector) – Pyxel Detector object.

  • filename (str or Path) – File path.

  • position (tuple) – Indices of starting row and column, used when fitting QE map to detector.

  • align (Literal) – Keyword to align the QE map to detector. Can be any from: (“center”, “top_left”, “top_right”, “bottom_left”, “bottom_right”)

  • seed (int, optional)

  • binomial_sampling (bool) – Binomial sampling. Default is True.

Notes

For more information, you can find an example here: H2RG pipeline.

Warning

Model assumes shot noise model was applied to photon array when using binomial sampling.

Apply QE curve#

Note

This model operates multi-wavelength photons.

PhotonCharge

With this model you can create and add charge to Detector via photoelectric effect by converting photons in charge. Loading QE vs wavelength values from a file to apply the QE to the photon array. Accepted file formats are .npy, .fits, .txt, .data and .csv. The column containing wavelength information should be in nanometers. After the photoconversion from photon to charge, applying the QE values to the photon array takes places and finally integrating along the wavelength dimension to get a 2D charge array as output.

Basic example of YAML configuration model:

- name: load_qe_curve
  func: pyxel.models.charge_generation.apply_qe_curve
  enabled: false
  arguments:
    filename: "qe_curve.csv"
    wavelength_col_name: "corrected lambda / nm"
    qe_col_name: "QE"
pyxel.models.charge_generation.apply_qe_curve(detector, filename, wavelength_col_name, qe_col_name)[source]#

Apply wavelength dependent QE with loading a file.

Parameters:
  • detector (Detector) – Pyxel Detector object.

  • filename (str or Path) – CSV File path.

  • wavelength_col_name (str) – Column name of wavelength in loaded file.

  • qe_col_name (str) – Column name of quantum efficiency in loaded file.

Conversion with 3D QE map#

Note

This model operates multi-wavelength photons.

PhotonCharge

With this model you can create and add charge to Detector via photoelectric effect by converting photons in charge. Loading QE values from a file to apply the QE to the photon array. Loading a 3D QE map from a file containing one QE map in the size of the detector per wavelength to apply the QE to the photon array. The file format must be netCDF, so ending with .nc to be able to read in. The file loaded will be interpreted as xarray.DataArray and should have the “wavelength” as coordinate, such that the wavelength resolution of the QE map data can be interpolated to match to the resolution of the wavelength used in the photon array. After that the photoconversion from photon to charge, applying the QE values to the photon array takes places and finally integrating along the wavelength dimension to get a 2D charge array as output.

Basic example of YAML configuration model:

- name: conversion_with_3d_qe_map
  func: pyxel.models.charge_generation.conversion_with_3d_qe_map
  enabled: true
  arguments:
    filename: "qe_map.nc"
pyxel.models.charge_generation.conversion_with_3d_qe_map(detector, filename)[source]#

Generate charge from incident photon via photoelectric effect.

Model converts photon to charge using custom QE map for different wavelengths.

Parameters:
  • detector (Detector) – Pyxel Detector object.

  • filename (str or Path) – File path.

Exponential absorption law#

Note

This model can operate with monochromatic and multi-wavelength photons.

Note

This model is valid only for silicon detectors at the moment.

PhotonCharge

With this model you can create and add charge to Detector via photoelectric effect by estimating the QE of your detector through an exponential absorption law ([2]).

The QE calculation depends on the detector type (Back-Illuminated or Front-Illuminated) and incorporates various parameters. The equations are as follows:

  1. Back-Illuminated (BI) Detector

The QE for a Back-Illuminated detector is given by:

\[QE = CCE \cdot (1 - R) \cdot \left(1 - e^{- \frac{x_{\text{epi}}}{\alpha}} \right)\]

With: \(CCE\): charge collection efficiency (fraction, between 0 and 1); \(R\): reflectivity (fraction, between 0 and 1); \(x_{\text{epi}}\): thickness of the epitaxial layer (in cm); \(\alpha\): absorptivity (in \(\text{cm}^{-1}\)).

  1. Front-Illuminated (FI) Detector

For a Front-Illuminated detector, the QE calculation includes the additional effect of the poly layer. The formula is:

\[QE = CCE \cdot (1 - R) \cdot e^{- \frac{x_{\text{poly}}}{\alpha}} \cdot \left(1 - e^{- \frac{x_{\text{epi}}}{\alpha}} \right)\]

With: \(x_{\text{poly}}\): thickness of the poly layer (in cm), and the other terms as defined above.

Temperature-Dependent Absorptivity Correction

The absorptivity, \(\alpha\), is corrected for temperature changes using the following formula:

\[\alpha' = \alpha \cdot e^{c \cdot \Delta T}\]

With: \(\alpha'\): adjusted absorptivity at the new temperature; \(\alpha\): absorptivity at the reference temperature (300 K); \(c\): temperature correction coefficient (in \(1/\text{K}\)). \(\Delta T\): temperature difference from the reference temperature (in \(\text{K}\)).

The embedded conversion coefficient \(c\) is wavelength and temperature-specific.

User-specified coefficients

In case you want to operate with your own conversion coefficients, you can add an additional column to the .csv file, with a \(c\) value for each working wavelength.

Warning

The reference model has been validated at 300 K; the conversion equation may be inaccurate for coefficients measured at different temperatures.

Once the QE is computed for the required wavelength(s) and at the desired temperature, it is applied to the detector photon array and a charge array is generated for the next steps of the pipeline.

Basic example of YAML configuration model for monochromatic pipeline:

- name: exponential_qe
  func: pyxel.models.charge_generation.exponential_qe
  enabled: true
  arguments:
      filename: qe_cleaned_data.csv
      x_epi: 0.0002
      detector_type: "BI" #or "FI"
      default_wavelength: 750.0 #in nm
      x_poly: 0.0001 #only in case of "FI" detector

Basic example of YAML configuration model for multi-wavelength pipeline:

- name: exponential_qe
  func: pyxel.models.charge_generation.exponential_qe
  enabled: true
  arguments:
      filename: qe_cleaned_data.csv
      x_epi: 0.0002
      detector_type: "BI" #or "FI"
      default_wavelength: 'multi'
      x_poly: 0.0001 #only in case of "FI" detector
pyxel.models.charge_generation.exponential_qe(detector, filename, x_epi, detector_type, default_wavelength=None, x_poly=0.0, cce=1.0, name=None)[source]#

Apply QE with temperature correction for absorptivity using a provided or backup coefficient c.

Parameters:
  • detector (Detector) – Pyxel Detector object.

  • filename (str or Path) – Path to the CSV file containing reflectivity, absorptivity, and c.

  • x_epi (float) – Thickness of the epitaxial layer in cm.

  • detector_type (str) – Type of detector (“BI” for Back-Illuminated, “FI” for Front-Illuminated).

  • default_wavelength (str or float) – Wavelength in nm for 2D photon arrays, or ‘multi’ for multiple wavelengths (no default value).

  • x_poly (float) – Thickness of the poly layer in cm.

  • cce (float) – Charge Collection Efficiency (default: 1.0).

  • name (str, optional) – Name to use for the result.

Load charge#

ChargeCharge

With this model you can add charge to Detector by loading charge values from a file. Accepted file formats are .npy, .fits, .txt, .data, .jpg, .jpeg, .bmp, .png and .tiff. Use argument position to set the offset from (0,0) pixel and set where the input charge is placed onto detector. You can set preset positions with argument align. Values outside of detector shape will be cropped. Read more about placement in the documentation of function fit_into_array(). Use argument time_scale to set the time scale of the input charge, default is 1 second.

Basic example of YAML configuration model:

- name: load_charge
  func: pyxel.models.charge_generation.load_charge
  enabled: true
  arguments:
    filename: data/charge.npy
    position: [0,0]
pyxel.models.charge_generation.load_charge(detector, filename, position=(0, 0), align=None, time_scale=1.0)[source]#

Load charge from txt file for detector, mostly for but not limited to CCDs.

Parameters:
  • detector (Detector) – Pyxel Detector object.

  • filename (str or Path) – File path.

  • position (tuple) – Indices of starting row and column, used when fitting charge to detector.

  • align (Literal) – Keyword to align the charge to detector. Can be any from: (“center”, “top_left”, “top_right”, “bottom_left”, “bottom_right”)

  • time_scale (float) – Time scale of the input charge, default is 1 second. 0.001 would be ms.

Notes

The Pixel coordinate conventions in Pyxel define the coordinate of pixel (0, 0) at the center of the leftmost bottom pyxel.

For example, when using the parameter align = "bottom_left", the image generated will be aligned to the lower leftmost pixel (0, 0) like this:

            ┌─────────────┐
   top    4  ...........           3  ...........      Size: 5 x 10
Y         2  xxxxx......      x: charge injected
          1  xxxxx......    bottom 0  xxxxx......             └─────────────┘
              0 2 4 6 8 9
             left    right
                   X

Charge injection#

Note

This model is specific for the CCD detector.

ChargeCharge

With this model you can inject arbitrary charge block into rows of a CCD detector. Charge will be injected uniformly from row number block_start to row number block_end.

Example of YAML configuration model:

- name: charge_blocks
  func: pyxel.models.charge_generation.charge_blocks
  enabled: true
  arguments:
    charge_level: 100
    block_start: 10
    block_end: 50
pyxel.models.charge_generation.charge_blocks(detector, charge_level, block_start=0, block_end=None)[source]#

Inject a block of charge into the CCD detector.

Parameters:
  • detector (Detector) – Pyxel Detector object.

  • charge_level (int) – Value of charges.

  • block_start (int) – Starting row of the injected charge.

  • block_end (int) – Ending row for the injected charge.

Charge deposition model#

ChargeCharge

With this model it is possible to simulate the deposition of charge in the detector by ionized particles using user-provided stopping power curves. It is possible to simulate mono-energetic beams (with a certain spread in energy) or provide an energy distribution (e.g., representative of the radiation environment). Stopping power curves for protons in silicon and for protons in MCT alloy are provided. Similarly, the proton energy distribution at L2 with and without 11-mm aluminium shielding is provided within Pyxel. This model is not as realistic as CosmiX but it is faster and easier to apply to a wide range of material and particles. In particular due to its simplistic nature, it fails at reproducing the deposition of only a small amount of charge.

Example of the configuration file:

- name: charge_deposition
  func: pyxel.models.charge_generation.charge_deposition
  enabled: true
  arguments:
    flux: 100
    step_size: 1.
    energy_mean: 1.
    energy_spread: .1
    energy_spectrum: data/proton_L2_solarMax_NoShielding.txt
    energy_spectrum_sampling: log
    ehpair_creation: 3.6
    material_density: 2.33
    particle_direction: isotropic
    stopping_power_curve: data/protons-in-silicon_stopping-power.csv
pyxel.models.charge_generation.charge_deposition(detector, flux, step_size=1.0, energy_mean=1.0, energy_spread=0.1, energy_spectrum=None, energy_spectrum_sampling='log', ehpair_creation=3.65, material_density=2.329, particle_direction='isotropic', stopping_power_curve=None, seed=None)[source]#

Simulate charge deposition by ionizing particles using a stopping power curve.

Parameters:
  • detector (Detector) – the detector

  • flux (float) – the flux of incoming particles in particle/s

  • step_size (float) – the size of the considered unitary step in unit length along which energy is deposited

  • energy_mean (float) – the mean energy of the incoming ionizing particles in MeV

  • energy_spread (float) – the spread in energy of the incoming ionizing particles in MeV

  • energy_spectrum (str) – the location of the file describing the energy spectrum of incident particles if no spectrum is provided energies are randomly drawn from a normal distribution with mean and spread defined above note that the energy spectrum is assumed to be a txt file with two columns [energy, flux] with the energy in MeV

  • energy_spectrum_sampling (str) – “log” or None: the energy spectrum is sampled in log space “linear” : the energy spectrum is sampled in linear space

  • ehpair_creation (float) – the energy required to generate a electron-hole pair in eV by default the Si value at room temperature is parsed i.e. 3.65 eV

  • material_density (float) – the material density in g/cm3 by default he Si value at room temperature is parsed i.e. 2.3290 g/cm3

  • particle_direction (str) – “isotropic” : particles are coming from all directions (outside of the sensor) “orthogonal” : particles are coming from the top of the sensor (thickness = 0) and orthogonal to its surface

  • stopping_power_curve (str) – the location of the file describing the total massive stopping power energetic loss per mass of material and per unit path length versus particle energy note that the the stopping power curve is assumed to be a csv file with two columns [energy, stopping power] energy in MeV, stopping power in MeV cm2/g

  • seed (int, optional)

Charge deposition model in MCT#

ChargeCharge

This model is the same as charge deposition model but is specific to MCT material. It computes the e-h pair creation (assuming it is 3 times the bandgap) and the alloy density based on the detector temperature and cut-off wavelength.

Example of the configuration file:

- name: charge_deposition
  func: pyxel.models.charge_generation.charge_deposition_in_mct
  enabled: true
  arguments:
    flux: 100
    step_size: 1.
    energy_mean: 1.
    energy_spread: .1
    energy_spectrum: data/proton_L2_solarMax_NoShielding.txt
    energy_spectrum_sampling: log
    cutoff_wavelength: 2.5
    particle_direction: isotropic
    stopping_power_curve: data/mct-stopping-power.csv
pyxel.models.charge_generation.charge_deposition_in_mct(detector, flux, step_size=1.0, energy_mean=1.0, energy_spread=0.1, energy_spectrum=None, energy_spectrum_sampling='log', cutoff_wavelength=2.5, particle_direction='isotropic', stopping_power_curve=None, seed=None)[source]#

Simulate charge deposition by ionizing particles using a stopping power curve.

Parameters:
  • detector (Detector) – the detector

  • flux (float) – the flux of incoming particles in particle/s

  • step_size (float) – the size of the considered unitary step in unit length along which energy is deposited

  • energy_mean (float) – the mean energy of the incoming ionizing particles in MeV

  • energy_spread (float) – the spread in energy of the incoming ionizing particles in MeV

  • energy_spectrum (str) – the location of the file describing the energy spectrum of incident particles if no spectrum is provided energies are randomly drawn from a normal distribution with mean and spread defined above note that the energy spectrum is assumed to be a txt file with two columns [energy, flux] with the energy in MeV

  • energy_spectrum_sampling (str. Default: 'log') – “log” : the energy spectrum is sampled in log space “linear” : the energy spectrum is sampled in linear space

  • cutoff_wavelength (float) – the longest wavelength in micrometer at which the QE reaches 50% of its maximum, used to compute the bandgap energy, and the corresponding fraction of cadmium

  • particle_direction (str. Default: 'isotropic') – “isotropic” : particles are coming from all directions (outside of the sensor) “orthogonal” : particles are coming from the top of the sensor (thickness = 0) and orthogonal to its surface

  • stopping_power_curve (str) – the location of the file describing the total massive stopping power energetic loss per mass of material and per unit path length versus particle energy note that the the stopping power curve is assumed to be a csv file with two columns [energy, stopping power] energy in MeV, stopping power in MeV cm2/g

  • seed (int, optional)

CosmiX cosmic ray model#

ChargeCharge

A cosmic ray event simulator was the first model added to Pyxel. Initially it was a simple, semi-analytical model in Fortran using the stopping power curve of protons to optimize the on-board source detection algorithm of the Gaia telescope to discriminate between stars and cosmic rays. Then it was reimplemented in Python as TARS (Tools for Astronomical Radiation Simulations) and later as CosmiX.

With this model you can add the effect of cosmic rays to the Charge data structure. See the documentation below for descriptions of parameters. CosmiX model is described in detail in [11].

  • Developed by: David Lucsanyi, ESA

Poppy

CosmiX cosmix ray model#

Example of the configuration file using default running_mode: stepsize with step size files, all with the same incident energy of 100 MeV and for 5 different thicknesses of 40 µm, 50 µm, 60 µm, 70 µm and 100 µm.

- name: cosmix
  func: pyxel.models.charge_generation.cosmix
  enabled: true
  arguments:
    simulation_mode: cosmic_ray
    running_mode: "stepsize"
    particle_type: proton
    initial_energy: 100.          # MeV
    particles_per_second: 100
    incident_angles:
    starting_position:
    spectrum_file: 'data/proton_L2_solarMax_11mm_Shielding.txt'
    seed: 4321

Note

You can find examples of this model in these Jupyter Notebooks from Pyxel Data:

Another example of the configuration file using default running_mode: stepsize with defined step size files.

- name: cosmix
  func: pyxel.models.charge_generation.cosmix
  enabled: true
  arguments:
    simulation_mode: cosmic_ray
    running_mode: "stepsize"
    particle_type: proton
    initial_energy: 100.          # MeV
    particles_per_second: 100
    incident_angles:
    starting_position:
    spectrum_file: 'data/proton_L2_solarMax_11mm_Shielding.txt'
    seed: 4321
    stepsize:
      - type: proton
        energy:    100.0  # MeV
        thickness: 40.0   # um
        filename:  pyxel/models/charge_generation/cosmix/data/stepsize_proton_100MeV_40um_Si_10k.ascii
      - type: proton
        energy:    100.0  # MeV
        thickness: 50.0   # um
        filename:  pyxel/models/charge_generation/cosmix/data/stepsize_proton_100MeV_50um_Si_10k.ascii
      - type: proton
        energy:    100.0  # MeV
        thickness: 60.0   # um
        filename:  pyxel/models/charge_generation/cosmix/data/stepsize_proton_100MeV_60um_Si_10k.ascii
      - type: proton
        energy:    100.0  # MeV
        thickness: 70.0   # um
        filename:  pyxel/models/charge_generation/cosmix/data/stepsize_proton_100MeV_70um_Si_10k.ascii
      - type: proton
        energy:    100.0  # MeV
        thickness: 100.0   # um
        filename:  pyxel/models/charge_generation/cosmix/data/stepsize_proton_100MeV_100um_Si_10k.ascii
pyxel.models.charge_generation.cosmix(detector, simulation_mode, running_mode, particle_type, particles_per_second, spectrum_file, initial_energy='random', incident_angles=('random', 'random'), starting_position=('random', 'random', 'random'), seed=None, ionization_energy=3.6, progressbar=True, stepsize=None)[source]#

Apply CosmiX model.

Parameters:
  • detector (Detector) – Pyxel Detector object.

  • simulation_mode (literal) – Simulation mode: cosmic_rays, radioactive_decay.

  • running_mode (literal) – Mode: stopping, stepsize, geant4, plotting.

  • particle_type – Type of particle: proton, alpha, ion.

  • particles_per_second (float) – Number of particles per second.

  • spectrum_file (str) – Path to input spectrum

  • initial_energy (int or float or literal) – Kinetic energy of particle, set random for random.

  • incident_angles (tuple of str) – Incident angles: (α, β).

  • starting_position (tuple of str) – Starting position: (x, y, z).

  • seed (int, optional) – Random seed.

  • ionization_energy (float) – Mean ionization energy of the semiconductor lattice.

  • progressbar (bool) – Progressbar.

  • stepsize (optional, list of dict) – Define the different step sizes. Only for running mode ‘stepsize’

Notes

For more information, you can find examples here:

Dark current#

ChargeCharge

With this model you can add a temperature dependent dark current to charge data, stored in the a Detector object. The model follows the description in [12]. The average dark current rate (in \(\mathit{e^-/s/pixel}\)) is:

\(D_R = \frac{D_{FM}P_{S}}{q}\frac{T^\frac{3}{2}e^{-\frac{E_{gap}}{2k_{B}T}}}{T_{room}^\frac{3}{2}e^{-\frac{E_{g,room}}{2k_{B}T_{room}}}}\)

where

\(T\) is temperature (in \(K\)), \(T_{room}\) room temperature (\(\mathit{300 K}\)), \(E_{g}\) band gap (in \(eV\)), \(k_B\) Boltzmann constant, \(D_{FM}\) dark current figure of merit (in \(nA/cm^{2}\)), \(P_S\) pixel area (in \(cm^{2}\)), \(q\) charge of an electron (in \(C\))and \(E_{g, room}\) band gap at room temperature. The entire dark current during exposure is:

\(I_{dark}=\mathcal{P}\big(t_{exp}D_R\big)\bigg(1+\mathcal{lnN}\big(0, \sigma^2_{fpn}\big)\bigg)\),

where \(\sigma_{fpn}=t_{exp} D_R D_N\), \(\mathcal{P}\) Poisson distribution, \(\mathcal{lnN}\) log-normal distribution, \(D_N\) the dark current spatial noise factor and \(t_{exp}\) exposure time (in \(s\)).

To use the model, user has to provide arguments figure_of_merit in \(\mathit{nA/cm^2}\) (\(D_{FM}\)), band_gap in \(\mathit{eV}\), band_gap_room_temperature in \(\mathit{eV}\), spatial_noise_factor (\(D_N\)) and temporal_noise. If temporal_noise is true, shot noise will be included. The spatial_noise_factor is typically between 0.1 and 0.4 for CCD and CMOS sensors [12].

Parameter temperature in \(\mathit{K}\) is taken from detector Environment. If arguments band_gap and band_gap_room_temperature are not provided, the model will use the Varshni empirical formula (see [13]) with parameters for Silicon by default:

\(E_{gap}(T) = E_{gap}(0) - \frac{\alpha T^2}{T+\beta}\).

For Silicon, material constants are \(E_{gap}(0)=1.1577\mathit{[eV]}\), \(\alpha=7.021\times10^{-4}\mathit{[eV/K]}\), and \(\beta=1108\mathit{[K]}\).

Example of the configuration file:

- name: dark_current
  func: pyxel.models.charge_generation.dark_current
  enabled: true
  arguments:
    figure_of_merit: 1.  # nA/cm^2
    band_gap: 1.2  # eV, optional
    band_gap_room_temperature: 1.2  # eV, optional
    spatial_noise_factor: 0.1
    temporal_noise: false

Note

You can find an example of this model used in this Jupyter Notebook Dark current versus temperature for silicon detectors from Pyxel Data.

pyxel.models.charge_generation.dark_current(detector, figure_of_merit, spatial_noise_factor=None, band_gap=None, band_gap_room_temperature=None, seed=None, temporal_noise=True)[source]#

Add dark current to the detector charge.

Based on: Konnik, Mikhail V. and James S. Welsh. “High-level numerical simulations of noise in CCD and CMOS photosensors: review and tutorial.” ArXiv abs/1412.4031 (2014): n. pag.

Parameters:
  • detector (Detector) – Pyxel detector object.

  • figure_of_merit (float) – Dark current figure of merit. Unit: nA/cm^2

  • spatial_noise_factor (float) – Dark current fixed pattern noise factor.

  • band_gap (float, optional) – Semiconductor band_gap. If none, the one for silicon is used. Unit: eV

  • band_gap_room_temperature (float, optional) – Semiconductor band gap at 300K. If none, the one for silicon is used. Unit: eV

  • seed (int, optional) – Random seed.

  • temporal_noise (bool, optional) – Shot noise.

Notes

For more information, you can find an example here: Dark current versus temperature for silicon detectors.

Dark current rule07#

Note

This model is specific for the MCT detector.

ChargeCharge

With this model you can add dark current to Charge following the model described in [14]. This model is only valid for MCT hybridised array (MCT). If temporal_noise is true, shot noise will be included. The model has one extra argument: cut-off wavelength, and also takes some values from Detector object, to be precise: temperature, pixel size (assuming it is square), and time step since last read-out. Please make sure the detector Environment, Geometry and Characteristics are properly set in the YAML configuration file.

Example of the configuration file:

- name: dark_current
  func: pyxel.models.charge_generation.dark_current_rule07
  enabled: true
  arguments:
    cutoff_wavelength: 2.5
    spatial_noise_factor: 0.1
    temporal_noise: true

Note

You can find an example of this model used in this Jupyter Notebook Dark current vs temperature for MCT detectors from Pyxel Data.

pyxel.models.charge_generation.dark_current_rule07(detector, cutoff_wavelength=2.5, spatial_noise_factor=None, seed=None, temporal_noise=True)[source]#

Generate charge from dark current process based on the Rule07 model.

The relationship between dark current and the temperature can often be modeled by exponential function, which follows a rule of thumb known as “Rule 7”. This rule states that the dark current approximately doubles for every 7°K increase in temperature.

Based on Rule07 paper by W.E. Tennant Journal of Electronic Materials volume 37, pages1406–1410 (2008).

Parameters:
  • detector (Detector) – The detector object (must be either a CCD or CMOS detector).

  • cutoff_wavelength (float, optional. Default: 2.5) – The wavelength cut-off for the detector (default is 2.5 µm). Unit: µm

  • spatial_noise_factor (float, optional) – Factor for introducing spatial noise (fixed pattern noise) in the dark current.

  • seed (int, optional) – Seed for random number generation to allow reproducibility of noise patterns.

  • temporal_noise (bool, optional) – If True, shot noise (Poisson noise) is added to the dark current.

Raises:
  • TypeError – If the detector is not an instance of either CCD or CMOS.

  • ValueError – If the cut-off wavelength is not within the range of 1.7 to 15.0 µm.

Notes

This function simulates dark current generation in a detector, considering both spatial noise (fixed pattern noise) and temporal noise (shot noise).

For more information, you can find an example here: Dark current vs temperature for MCT detectors.

Simple dark current#

ChargeCharge

With this model you can add dark current to a Detector object.

Example of the configuration file:

- name: simple_dark_current
  func: pyxel.models.charge_generation.simple_dark_current
  enabled: true
  arguments:
    dark_rate: 10.0
pyxel.models.charge_generation.simple_dark_current(detector, dark_rate, seed=None)[source]#

Simulate dark current in a detector.

Parameters:
  • detector (Detector) – Any detector object.

  • dark_rate (float) – Dark current, in e⁻/pixel/second, which is the way manufacturers typically report it.

  • seed (int, optional)

APD gain#

Note

This model is specific to the APD detector.

ChargeCharge

With this model you can apply APD gain to the a APD object. Model simply multiplies the values of charge with the avalanche gain, which should be specified in the detector characteristics.

Example of the configuration file:

- name: apd_gain
  func: pyxel.models.charge_generation.apd_gain
  enabled: true

Note

You can find an example of this model used in this Jupyter Notebook Saphira from Pyxel Data.

pyxel.models.charge_generation.apd_gain(detector)[source]#

Apply APD gain.

Parameters:

detector (APD) – Pyxel APD detector object.

Notes

For more information, you can find an example here: Saphira.

Dark current Saphira#

Note

This model is specific to the APD detector.

ChargeCharge

With this empirical model you can add dark current to a APD object. The model is an approximation the dark current vs. gain vs. temp plot in [15], Fig. 3. We can split it into three linear ‘regimes’: 1) low-gain, low dark current; 2) nominal; and 3) trap-assisted tunneling. The model ignores the first one for now since this only applies at gains less than ~2. All the necessary arguments are provided through the detector characteristics. The model works best for temperature less than 100 and avalanche gain more than 2.

Example of the configuration file:

- name: dark_current_saphira
  func: pyxel.models.charge_generation.dark_current_saphira
  enabled: true

Note

Dark current calculated with this model already takes into account the avalanche gain.

Note

You can find an example of this model used in this Jupyter Notebook Saphira from Pyxel Data.

pyxel.models.charge_generation.dark_current_saphira(detector, seed=None)[source]#

Simulate dark current in a Saphira APD detector.

Reference: I. M. Baker et al., Linear-mode avalanche photodiode arrays in HgCdTe at Leonardo, UK: the current status, in Image Sensing Technologies: Materials, Devices, Systems, and Applications VI, 2019, vol. 10980, no. May, p. 20.

Parameters:
  • detector (APD) – An APD detector object.

  • seed (int, optional)

Notes

For more information, you can find an example here: Saphira.

Radiation induced Dark Current#

ChargeCharge

This model adds dark current induced by radiation. A more detailed description of the models can be found in [16] and [17].

Example of configuration file:

- name: radiation_induced_dark_current
  func: pyxel.models.charge_generation.radiation_induced_dark_current
  enabled: true
  arguments:
    depletion_volume: 64 # µm3
    annealing_time: 0.1 # weeks
    displacement_dose:  50  # TeV/g
    shot_noise: false

Note

You can find an example of this model used in this Jupyter Notebook Simulating radiation induced dark current from Pyxel Data.

pyxel.models.charge_generation.radiation_induced_dark_current(detector, depletion_volume, annealing_time, displacement_dose, shot_noise, seed=None)[source]#

Model to add dark current induced by radiation to the detector charge.

The radiation induced dark current model description can be found in [16] and [17].

Parameters:
  • detector (Detector) – Pyxel detector object.

  • depletion_volume (float) – Depletion volume parameter. Unit: µm3.

  • annealing_time (float) – Annealing time. Unit: weeks

  • displacement_dose (float) – Displacement dose parameter. Unit: TeV/g

  • shot_noise (bool) – True to enable shot noise.

  • seed (int, optional)

Notes

For more information, you can find an example here: Simulating radiation induced dark current.