autoemxsp.core.XSp_quantifier module
X-ray Spectrum Quantification Module
Created on Thu Jun 27 14:23:22 2024
@author: Andrea Giunto
This module provides classes and functions for matrix correction and quantitative analysis of X-ray spectra, implementing the peak-to-background (P/B) method and related ZAF corrections for scanning electron microscopy (SEM) energy-dispersive X-ray spectroscopy (EDS).
Class Structure and Interactions
The main classes are: - XSp_Quantifier
Performs quantification of X-ray spectra using matrix corrections. Interfaces with spectral fitting routines from XSp_Fitter to extract peak/background intensities and applies the full suite of correction factors from Quant_Corrections to obtain quantitative elemental concentrations. Additionally provides functions to plot and print the results.
Quant_Corrections Provides methods for calculating matrix correction factors (Z, A, R) for the P/B method. This class is instanced by XSp_Quantifier to calculate the matrix correction factors.
Typical Usage
- Initialize the quantifier object:
quantifier = XSp_Quantifier(…) # See class docs for initialization
- Quantify a spectrum:
- quant_result = quantifier.quantify_spectrum()
quant_result contains quantified atomic and weight fractions, analytical error, and spectral fitting metrics (reduced chi-square, and R-squared).
- Print quantification results:
quantifier.print_quant_result(quant_result)
- Plot spectrum:
quantifier.plot_quantified_spectrum()
Customization & Calibration
Detector calibration, physical constants, and mass absorption coefficients are handled via the calibs module and supporting utility functions. Users should calibrate the values in the calibs module corresponding to their microscope and EDS settings prior usage.
Dependencies
numpy, scipy, sympy, pandas lmfit, pymatgen.core.Element (from supporting libraries) XSp_Fitter (required for spectral fitting; see separate module) calibs, lib modules
How the classes interact
XSp_Quantifier is the main user-facing class for quantification. It uses Quant_Corrections to compute all correction factors and matrix effects, and relies on the XSp_Fitter module to extract peak/background intensities from measured spectra. Quant_Corrections provides the core physical models for all matrix corrections, and is initialized with the relevant sample and measurement parameters. The XSp_Fitter module must be installed and available for spectral fitting and deconvolution.
- class autoemxsp.core.XSp_quantifier.XSp_Quantifier(spectrum_vals, spectrum_lims, microscope_ID, meas_type, meas_mode, det_ch_offset, det_ch_width, beam_e, emergence_angle, energy_vals=None, background_vals=None, els_sample=None, els_substrate=None, els_w_fr=None, is_particle=False, sp_collection_time=None, max_undetectable_w_fr=0.1, fit_tol=0.0001, standards_dict=None, verbose=False, fitting_verbose=False)[source]
Bases:
objectClass for quantitative analysis of EDS spectra.
- spectrum_vals
The measured EDS spectrum (counts per channel).
- Type:
numpy.ndarray
- energy_vals
Energy scale corresponding to spectrum_vals.
- Type:
numpy.ndarray
- spectrum_lims
Start and end indices for the spectrum region to analyze.
- Type:
tuple of int
- fit_background
Whether to fit a background model. True, if background_vals are not provided.
- Type:
bool
- background_vals
Fitted or provided background spectrum counts.
- Type:
numpy.ndarray
- els_sample
All elements present in the sample, including undetectable elements.
- Type:
list of str
- els_to_quantify
Elements to be quantified (excluding undetectable).
- Type:
list of str
- els_substrate
Elements present in the substrate or sample holder appearing in spectra (excluding undetectable).
- Type:
list of str
- els_w_fr
Dictionary of fixed elemental mass fractions, by element symbol, e.g. {‘Si’:0.33, ‘O’:0.67}.
- Type:
dict or None, optional
- max_undetectable_w_fr
Maximum allowed total mass fraction for undetectable elements during fitting.
- Type:
float
- force_total_w_fr
Whether total mass fraction is normalized 1 during fitting. False, if undetectable elements are present in the sample. True otherwise.
- Type:
bool
- is_particle
Whether the sample is a particle (affects fitting and quantification corrections).
- Type:
bool
- sp_collection_time
Live time of spectrum acquisition (in seconds).
- Type:
float or None
- fit_tol
Tolerance for spectrum fitting convergence.
- Type:
float
- bad_quant_flag
Flag indicating quantification issues, or None if successful.
- Type:
int or None
- microscope_ID
Microscope identifier for calibration.
- Type:
str
- meas_type
Measurement type (e.g., ‘EDS’).
- Type:
str
- meas_mode
Measurement mode, defining detector calibrations and (optionally) beam current (e.g., ‘point’).
- Type:
str
- det_ch_offset
Detector channel energy offset (keV).
- Type:
float
- det_ch_width
Detector channel width (keV).
- Type:
float
- beam_energy
Electron beam energy (keV).
- Type:
float
- emergence_angle
Detector emergence (take-off) angle (degrees).
- Type:
float
- verbose
Print information during quantification.
- Type:
bool
- fitting_verbose
If True, print detailed information during each fitting step.
- Type:
bool
- Class attributes
- ----------------
- xray_quant_ref_lines
List of X-ray lines used as reference
- Type:
list(str)
- xray_quant_ref_lines = ['Ka1', 'La1', 'Ma1', 'Mz1']
- initialize_and_fit_spectrum(params: Parameters | None = None, print_results: bool | None = False) None[source]
Perform a complete fit of the spectrum provided to the instance of XSp_Quantifier. This method is intended for single-iteration fits where sample composition is known and constrained (e.g., EDS standards).
The fit results and all relevant fitted parameters are stored as instance attributes.
- Parameters:
params (Optional[lmfit.Parameters], optional) – Parameters object to pass to the spectrum fitter. If None, default parameters are used.
print_result (bool, optional) – If True, prints the fitting results.
- Returns:
None
Nested Calls
————
Calls the following methods to further store values extracted from the fit –
self._initialize_spectrum_fitter(): initializes the spectrum fitter
self._fit_spectrum(): fits the spectrum and stores the fit results
- get_starting_K_val() float | None[source]
Estimate the initial background scaling factor K by fitting the high-energy portion of the spectrum.
This method fits only the region of the spectrum above a threshold energy to avoid regions heavily affected by absorption. It is intended to provide an optimal starting value for K, especially for particle spectra, to prevent the algorithm from compensating for background intensity via particle geometry parameters, which are only intended for background shape fitting.
- Returns:
K_val – Estimated initial value for the background scaling factor K, or None if a suitable value could not be determined.
- Return type:
Optional[float]
Notes
Uses XSp_Fitter with particle geometry disabled for this quick fit.
Prints diagnostic information if self.verbose is True.
- quantify_spectrum(force_single_iteration=False, interrupt_fits_bad_spectra=True, print_result=True)[source]
Quantifies a spectrum, using iterative fitting to refine elemental fractions. At each iteration, elemental fractions are quantified and enforced in the next iteration. Algorithm converges when quantified elemental fractions all change by less than 0.01%.
Instead, if elemental fractions are provided (i.e., fitting of experimental standards), background values are provided, or if force_single_iteration = True, only a single iteration is performed.
- Parameters:
force_single_iteration (bool, optional) – If True, only a single fit iteration is performed, even if some elemental fractions are undefined.
interrupt_fits_bad_spectra (bool, optional) – If True, fitting is interrupted early if the spectrum is deemed unreliable.
print_result (bool, optional) – If True, prints the quantification results.
- Returns:
quant_result (dict or None) – Dictionary with quantification results, or None if fit failed.
min_bckgrnd_ref_lines (float) – Minimum background counts across reference peaks, or 0 if fit failed.
bad_quant_flag (int or None) – Flag indicating quantification issues, or None if successful.
- print_quant_result(quant_result: dict, Z_sample: dict = None) None[source]
Print the quantification results, including fit quality metrics and elemental composition.
- Parameters:
quant_result (dict) –
- Dictionary containing quantification results with keys:
cnst.REDCHI_SQ_KEY: Reduced Chi-square of the fit.
cnst.R_SQ_KEY: R-squared of the fit.
cnst.COMP_W_FR_KEY: Weight fractions (as decimals).
cnst.COMP_AT_FR_KEY: Atomic fractions (as decimals).
cnst.AN_ER_KEY: Analytical error (as decimal).
Z_sample (dict, optional) –
- Dictionary containing mean atomic numbers (optional), e.g.:
’Statham2016’: Mean atomic number (Z̅) calculated according to Statham (2016).
’mass-averaged’: Mean atomic number (Z̅) weighted by composition.
- Return type:
None
References
Statham P, Penman C, Duncumb P. IOP Conf Ser Mater Sci Eng. 2016;109(1):0–10.
- plot_quantified_spectrum(annotate_peaks: str = 'all', plot_bckgrnd_cnts_ref_peaks: bool = True, plot_initial_guess: bool = False, plot_title: str | None = None, peaks_to_zoom: str | List[str] | None = None) None[source]
Plot the quantified spectrum. The background counts under reference peaks are highlighted for spectra used for quantification.
- Parameters:
annotate_peaks (str, optional) – Which peaks to annotate. Options: ‘all’, ‘most’, ‘main’, ‘none’. Default is ‘all’.
plot_bckgrnd_cnts_ref_peaks (bool, optional) – If True, plot vertical lines that illustrate value of background counts used for quantification. This is shown underneath the corresponding reference characteristic peaks for each element.
plot_initial_guess (bool, optional) – If True, plot the initial guess as well. Default is False.
plot_title (str, optional) – Title printed at the top of the plot. Default is None.
peaks_to_zoom (str or list of str, optional) – Peak label (e.g. ‘Si_Ka1’) or list of labels to zoom in on. If provided, creates a new figure for each.
- Return type:
None.
- class autoemxsp.core.XSp_quantifier.Quant_Corrections(elements: Sequence[str], beam_energy: float, emergence_angle: float, meas_mode: str, energies: Sequence[float] | ndarray | None = None, verbose: bool = False)[source]
Bases:
objectImplements matrix correction factors for quantitative X-ray microanalysis using the peak-to-background (P/B) method.
This class provides methods for calculating Z (atomic number), A (absorption), and R (backscattering) correction factors, as well as mass absorption coefficients for a given set of elements and measurement conditions. It is designed for both standard-based and standardless quantification workflows in electron probe microanalysis (EPMA) or energy-dispersive X-ray spectroscopy (EDS).
References
- Love, V.D. Scott, “Evaluation of a new correction procedure for quantitative electron probe microanalysis”,
Phys. D. Appl. Phys. 11 (1978) 1369–1376. https://doi.org/10.1088/0022-3727/11/10/002
- P.J. Statham, “A ZAF PROCEDURE FOR MICROPROBE ANALYSIS BASED ON MEASUREMENT OF PEAK-TO-BACKGROUND RATIOS”,
in: D.E. Newbury (Ed.), Fourteenth Annu. Conf. Microbeam Anal. Soc., San Francisco Press, 1979: pp. 247–253.
- Essani, E. Brackx, E. Excoffier, “A method for the correction of size effects in microparticles using a peak-to-background approach
in electron-probe microanalysis”, Spectrochim. Acta B 169 (2020) 105880. https://doi.org/10.1016/j.sab.2020.105880
- elements
List of element symbols included in the quantification (excluding undetectable elements).
- Type:
list of str
- energies
X-ray energies (keV) for each element/line.
- Type:
np.ndarray
- emergence_angle
Detector emergence angle (degrees).
- Type:
float
- beam_energy
Incident electron beam energy (keV).
- Type:
float
- meas_mode
Detector mode (for calibration parameters).
- Type:
str
- Z_els
Atomic numbers for each element.
- Type:
np.ndarray
- W_els
Atomic weights for each element.
- Type:
np.ndarray
- els_nu
Backscattering coefficients for each element.
- Type:
np.ndarray
- mass_abs_coeffs_lines
Mass absorption coefficients for each element at each characteristic energy.
- Type:
list of list of float
- verbose
If True, enables verbose output for debugging.
- Type:
bool
- get_ZAF_mult_f_pb(weight_fractions: ndarray, el_lines_energies_d: Dict[str, float] | None = None) Tuple[ndarray, Dict[str, float]][source]
Calculate the ZAF multiplicative correction factors for the measured sample P/B ratio.
- This method accounts for:
Differences in average Z between the sample and the employed standard, which affect continuum intensity.
Second-order corrections for backscattering and absorption due to differential mean generation path between characteristic and continuum X-rays.
Fluorescence and particle-size corrections are ignored.
- Parameters:
weight_fractions (np.ndarray) – Estimated weight fractions of the elements to quantify.
el_lines_energies_d (dict[str, float], optional) – Dictionary mapping peak labels (e.g., ‘Fe_Ka’) to energies (keV). If None, uses all elements and energies in the class.
- Returns:
ZAF_pb_corrections (np.ndarray) – ZAF correction factors to multiply with the measured sample P/B ratio.
Z_sample (dict) – Dictionary with various sample mean atomic numbers.
References
- [1] Statham P, Penman C, Duncumb P. Improved spectrum simulation for validating SEM-EDS analysis.
IOP Conf. Ser. Mater. Sci. Eng. 109, 0 (2016).
- [2] Markowicz AA, Van Grieken RE. Composition dependence of bremsstrahlung background in electron-probe
x-ray microanalysis. Anal. Chem. 56, 2049 (1984).
Potential improvements
Include fluorescence corrections for large particles Include particle size corrections, from:
- [1] J. L. Lábár and S. Török, A peak‐to‐background method for electron‐probe x‐ray microanalysis
applied to individual small particles, X-Ray Spectrom. 21, 183 (1992).
- [2] M. Essani, E. Brackx, and E. Excoffier, A method for the correction of size effects in microparticles
using a peak-to-background approach in electron-probe microanalysis, Spectrochim. Acta - Part B At. Spectrosc. 169, 105880 (2020).