mean_charge_capture#

optika.sensors.mean_charge_capture(width_diffusion, width_pixel)[source]#

A function to compute the mean charge capture [Stern et al., 2004], the fraction of charge from each photon event retained in the central pixel.

Parameters:
Return type:

AbstractScalar

Examples

Plot the mean charge capture as a function of wavelength and energy for the sensor parameters in Heymes et al. [2020].

import matplotlib.pyplot as plt
import astropy.units as u
import astropy.visualization
import named_arrays as na
import optika

# Define a grid of wavelengths
wavelength = na.geomspace(1, 10000, axis="wavelength", num=1001) * u.AA

# Convert the grid to energies as well
energy = wavelength.to(u.eV, equivalencies=u.spectral())

# Load the optical properties of silicon
si = optika.chemicals.Chemical("Si")

# Retrieve the absorption coefficient of silicon
# for the given wavelenghts.
absorption = si.absorption(wavelength)

# Compute the charge diffusion
width_diffusion = optika.sensors.charge_diffusion(
    absorption=absorption,
    thickness_substrate=14 * u.um,
    thickness_depletion=2.4 * u.um,
)

# Compute the mean charge capture
mcc = optika.sensors.mean_charge_capture(
    width_diffusion=width_diffusion,
    width_pixel=16 * u.um,
)

# Plot the mean charge capture as a function
# of wavelength and energy
with astropy.visualization.quantity_support():
    fig, ax = plt.subplots()
    ax2 = ax.twiny()
    ax2.invert_xaxis()
    na.plt.plot(
        wavelength,
        mcc,
        ax=ax,
    )
    na.plt.plot(
        energy,
        mcc,
        ax=ax2,
        linestyle="None",
    )
    ax.set_xscale("log")
    ax2.set_xscale("log")
    ax.set_xlabel(f"wavelength ({ax.get_xlabel()})")
    ax2.set_xlabel(f"energy ({ax2.get_xlabel()})")
    ax.set_ylabel(f"mean charge capture")
../_images/optika.sensors.mean_charge_capture_0_0.png

Notes

Naively, the mean charge capture (MCC) is the integral of the charge diffusion kernel over the extent of a pixel. However, since a photon can strike anywhere within the central pixel, the charge diffusion kernel should be convolved with a rectangle function the width of a pixel before integrating. So, our definition for the MCC is

\[P_\text{MCC} = \left\{ \frac{1}{d} \int_{-d/2}^{d/2} \left[ K(x') * \Pi \left( \frac{x'}{d} \right) \right](x) \, dx \right\}^2,\]

where \(K(x)\) is the charge diffusion kernel, \(\Pi(x)\) is the rectangle function, and \(d\) is the width of a pixel. If we assume that the charge diffusion kernel is a Gaussian with standard deviation \(\sigma\),

\[K(x) = \frac{1}{\sqrt{2\pi} \sigma} \exp \left( -\frac{x^2}{2 \sigma^2} \right),\]

then we can analytically solve for the MCC,

\[\begin{split}P_\text{MCC} &= \left\{ \frac{1}{2d} \int_{-d/2}^{d/2} \left[ \text{erf} \left( \frac{d - 2x}{2 \sqrt{2} \sigma} \right) + \text{erf} \left( \frac{d + 2x}{2 \sqrt{2} \sigma} \right) \right] dx \right\}^2 \\ &= \left\{ \sqrt{\frac{2}{\pi}} \frac{\sigma}{d} \left[ \exp \left( -\frac{d^2}{2 \sigma^2} \right) - 1 \right] + \text{erf} \left( \frac{d}{\sqrt{2} \sigma} \right) \right\}^2,\end{split}\]

where \(\text{erf}(x)\) is the error function.