e2v_ccd97#
- optika.sensors.materials.e2v_ccd97(temperature=<Quantity 300. K>)[source]#
A model of the light-sensitive material of an e2v CCD97 sensor based on measurements by Moody et al. [2017] and Heymes et al. [2020].
This is a measurement of e2v’s “enhanced” process, which has a narrower partial charge collection region than e2v’s “standard” process.
This model uses
e2v_ccd64_thick()to represent the depletion region.- Parameters:
temperature (Quantity | AbstractScalar) – The temperature of the light-sensitive silicon substrate.
- Return type:
Examples
Plot the measured e2v CCD97 quantum efficiency vs the fitted quantum efficiency calculated using the method of Stern et al. [1994].
import matplotlib.pyplot as plt import astropy.units as u import astropy.visualization import named_arrays as na import optika # Create a new instance of the e2v CCD97 light-sensitive material material = optika.sensors.materials.e2v_ccd97() # Store the wavelengths at which the QE was measured wavelength_measured = material.eqe_measured.inputs # Store the QE measurements eqe_measured = material.eqe_measured.outputs # Define a grid of wavelengths with which to evaluate the fitted QE wavelength_fit = na.geomspace(5, 10000, axis="wavelength", num=1001) * u.AA # Evaluate the fitted QE using the given wavelengths eqe_fit = material.quantum_efficiency_effective( rays=optika.rays.RayVectorArray( wavelength=wavelength_fit, direction=na.Cartesian3dVectorArray(0, 0, 1), ), normal=na.Cartesian3dVectorArray(0, 0, -1), ) # Plot the measured QE vs the fitted QE with astropy.visualization.quantity_support(): fig, ax = plt.subplots(constrained_layout=True) na.plt.scatter( wavelength_measured, eqe_measured, label="measured", ) na.plt.plot( wavelength_fit, eqe_fit, label="fit", ) ax.set_xscale("log") ax.set_xlabel(f"wavelength ({wavelength_fit.unit:latex_inline})") ax.set_ylabel("quantum efficiency") ax.legend()
The thickness of the oxide layer found by the fit is
material.thickness_oxide
\[46.109497 \; \mathrm{\mathring{A}}\]The thickness of the implant layer found by the fit is
material.thickness_implant
\[1534.8864 \; \mathrm{\mathring{A}}\]The thickness of the substrate is
material.thickness_substrate
\[14 \; \mathrm{\mu m}\]The differential charge collection efficiency at the backsurface found by the fit is
material.cce_backsurface
\[0.33363706 \; \mathrm{}\]And the roughness of the substrate found by the fit is
material.roughness_substrate
\[4.7503048 \; \mathrm{nm}\]Now plot the effective quantum efficiency of the fit to this data vs. the fit to the data in
optika.sensors.materials.e2v_ccd203()material_ccd_aia = optika.sensors.materials.e2v_ccd203() eqe_fit_aia = material_ccd_aia.quantum_efficiency_effective( rays=optika.rays.RayVectorArray( wavelength=wavelength_fit, direction=na.Cartesian3dVectorArray(0, 0, 1), ), normal=na.Cartesian3dVectorArray(0, 0, -1), ) with astropy.visualization.quantity_support(): fig, ax = plt.subplots(constrained_layout=True) na.plt.scatter( wavelength_measured, eqe_measured, label="Heymes 2020 measured", ) na.plt.plot( wavelength_fit, eqe_fit, label="Heymes 2020 fit", ) na.plt.scatter( material_ccd_aia.eqe_measured.inputs, material_ccd_aia.eqe_measured.outputs, label="Boerner 2012 measured", ) na.plt.plot( wavelength_fit, eqe_fit_aia, label="Boerner 2012 fit", ) ax.set_xscale("log") ax.set_xlabel(f"wavelength ({wavelength_fit.unit:latex_inline})") ax.set_ylabel("quantum efficiency") ax.legend()
Plot the width of the charge diffusion kernel for this sensor as a function of wavelength.
# Compute the width of the charge diffusion kernel # for each wavelength. width = material.width_charge_diffusion( rays=optika.rays.RayVectorArray( wavelength=wavelength_fit, direction=na.Cartesian3dVectorArray(0, 0, 1), ), normal=na.Cartesian3dVectorArray(0, 0, -1), ) # Plot the results with astropy.visualization.quantity_support(): fig, ax = plt.subplots() na.plt.plot( wavelength_fit, width, ax=ax, ) ax.set_xscale("log") ax.set_xlabel(f"wavelength ({ax.get_xlabel()})") ax.set_ylabel(f"width ({ax.get_ylabel()})")
References to
optika.sensors.materials.e2v_ccd97