MatrixCorrection - Part of the NeXL X-ray Microanalysis Library

Installing NeXLMatrixCorrection.jl

NeXLMatrixCorrection is available throught the standard Julia registry and can be installed using the standard package manager.

julia> ]add NeXLMatrixCorrection

or

julia> using Pkg
julia> Pkg.add("NeXLMatrixCorrection")

The matrix correction package in the NeXL microanalysis library for Julia. NeXLMatrixCorrection depends upon NeXLUncertainties and NeXLCore which are automatically installed when NeXLMatrixCorrection is installed.

Currently NeXLMatrixCorrection implements the XPP, CitZAF and Riveros matrix correction and Reed fluorescence correction algorithms for bulk and coated samples. The library is designed to make it easy to add additional algorithms.

Primarily the algorithms in NeXLMatrixCorrection are designed to take a Vector{NeXLCore.KRatio} and return a NeXLCore.Material. Since they are intended for both WDS and EDS, the k-ratio can represent one or more characteristic X-ray lines from a single element. K-ratios compare a measured intensity with the intensity from a reference (standard) material. Typically, these two materials are measured at the same beam energy but multiple beam energy measurements are also supported.

The primary method is

function quantify(
    iter::Iteration,  # The Iteration object providing algorithmic details
    label::Label,     # A label for the unknown
    measured::Vector{KRatio}, # A complete set of k-ratios (one per element)
    maxIter::Int = 100, # The maximum number of iterations to try before
)::IterationResult

which is wrapped as

quantify(iter::Iteration, sampleName::String, measured::Vector{KRatio})

and

quantify(ffr::FilterFitResult,...)

to simplify usage.

The KRatio structure is defined in NeXLCore.

KRatio(
    lines::Vector{CharXRay},
    unkProps::Dict{Symbol,<:Any},
    stdProps::Dict{Symbol,<:Any},
    standard::Material,
    kratio::AbstractFloat,
)

Usually it is sufficient to define the unkProps and stdProps corresponding to the :BeamEnergy, the :TakeOffAngle which are, of course, in eV and radians.

KRatio(
  characteristic(n"F",ktransitions), # Builds a vector containing all the K shell characteristic x-rays for F
  # [ n"F K-L3" ], # An alternative with only one transition
  Dict(:BeamEnergy=>15.0e3, :TakeOffAngle=>deg2rad(40.0)),
  Dict(:BeamEnergy=>15.0, :TakeOffAngle=>deg2rad(40.0)),
  mat"CaF2",  # The standard material
  0.324 # The k-ratio
)

Structure and Method Documentation

NeXLMatrixCorrection.ConvergenceTestType

The ConvergenceTest abstract type represents mechanisms to decide when the iteration has converged.

converged(ct::ConvergenceTest, meas::Vector{KRatio}, computed::Dict{Element,Float64})::Bool
source
NeXLMatrixCorrection.IterationType

Collects the information necessary to define the iteration process including the MatrixCorrection and FLuorescenceCorrection algorithms, the iteration UpdateRule, the ConvergenceTest, and an UnmeasuredElementRule.

source
NeXLMatrixCorrection.IterationResultType

IterationResult contains the results of the iteration process including a Label identifying the source of the k-ratios, the resulting Material, the initial and final k-ratios, whether the iteration converged and the number of steps. The results can be output using asa(DataFrame, ir::IterationResult).

source
NeXLMatrixCorrection.KRatioOptimizerType
KRatioOptimizer

Defines an optimizeks(kro::KRatioOptimizer, krs::Vector{KRatio})::Vector{KRatio} method which takes a vector of k-ratios which may have redundant data (more than one KRatio per element) and trims it down to a vector of k-ratios with one KRatio per element.

source
NeXLMatrixCorrection.MatrixCorrectionType

MatrixCorrection is an abstract type for computing ϕ(ρz)-type matrix correction algorithms. A sub-class MCA <: MatrixCorrection should implement

# Integral of the area under the ϕ(ρz)-curve
F(mc::MCA)
# Integral of the transmitted area under the ϕ(ρz)-curve
Fχ(mc::MCA, xray::CharXRay, θtoa::Real)
# Area under 0 to τ under the transmitted ϕ(ρz)-curve
Fχp(mc::MCA, xray::CharXRay, θtoa::Real, τ::Real)
# The ϕ(ρz)-curve
ϕ(mc::MCA, ρz)
# A factory method from MCA
matrixcorrection(::Type{MCA}, mat::Material, ashell::AtomicSubShell, e0)::MCA

The algorithm should precompute as much as possible based on the Material, AtomicSubShell and beam energy. The class MCA should define member variables subshell::AtomicSubShell, material::Material and E0::AbstractFloat to store the input values. See XPP for an example.

From these methods, other methods like Z(...), A(...) are implemented.

source
NeXLMatrixCorrection.MultiZAFType

The MultiZAF structure holds the information necessary to perform matrix correction on a collection of characteristic X-rays that were measured simultaneously from the same element. Use zafcorrection(...) to construct these rather than the internal constructor.

source
NeXLMatrixCorrection.RecordingUpdateRuleType

The RecordingUpdateRule wraps other UpdateRule instances to provide diagnostics which may be tabulated using asa(DataFrame, rur::RecordingUpdateRule) or plotted using Gadfly's plot(rur::RecordingUpdateRule).

source
NeXLMatrixCorrection.Riveros1993Type

@article{riveros1993review, title={Review of ϕ(ρz) curves in electron probe microanalysis}, author={Riveros, Jose and Castellano, Gustavo}, journal={X-Ray Spectrometry}, volume={22}, number={1}, pages={3–10}, year={1993}, publisher={Wiley Online Library} }

The instruction in Packwood1991 in Electron Probe Quantitation is: "For compounds weight averaging is used for all appropriate variables: Z, Z/A, η, and (Z/A)log(1.166(E0-Ec)/2J)"

source
NeXLMatrixCorrection.SimpleKRatioOptimizerType
SimpleKRatioOptimizer

Implements a simple optimizer based on shell first, overvoltage next and brightness last. Once it picks an optimum set of lines for an element, it will not change.

source
NeXLMatrixCorrection.UnmeasuredElementRuleType

The UnmeasuredElementRule mechanism provides a method to implement rules for adding unmeasured elements to the fitting process. Examples include element-by-stoichiometry or element-by-difference.

source
NeXLMatrixCorrection.WegsteinUpdateRuleType

The WegsteinUpdateRule implements the very effective method of Reed and Mason (S.J.B. Reed and P.K. Mason, Transactions ofthe Second National Conference on Electron Microprobe Analysis, Boston, 1967.) for updating the composition estimate between iteration steps.

source
Base.rangeFunction
range(::Type{XPP}, mat::MaterialLabel, e0, inclDensity=true)
range(ty::Type{<:BetheEnergyLoss}, mat::Material, e0::Float64, inclDensity=true; emin=50.0, mip::Type{<:NeXLMeanIonizationPotential}=Berger1982)

Total trajectory (range) of an electron with initial energy e0 (eV). (Units: inclDensity ? cm : cm/(g/cm³))

source
NeXLCore.JMethod
J(::Type{XPP}, mat::Material)

Mean ionization potential for the specified material in eV. (PAP1991 Eqn 6)

source
NeXLCore.atomicsubshellMethod
NeXLCore.atomicsubshell(mc::MatrixCorrection)

The sub-shell for which this MatrixCorrection has been calculated.

source
NeXLCore.materialMethod
NeXLCore.material(mc::MatrixCorrection) = mc.material

The material for which this MatrixCorrection has been calculated.

source
NeXLCore.transmissionMethod
NeXLCore.transmission(zaf::Coating, xray::CharXRay, toa)

Calculate the transmission fraction for the specified X-ray through the coating in the direction of the detector.

source
NeXLCore.transmissionMethod
NeXLCore.transmission(zaf::NullCoating, xray::CharXRay, toa)

Calculate the transmission fraction for the specified X-ray through no coating (1.0).

source
NeXLMatrixCorrection.AMethod
A(unk::MatrixCorrection, std::MatrixCorrection, xray::CharXRay, χcunk=0.0, tcunk=0.0, χcstd=0.0, tcstd=0.0)

The absorption correction factors.

source
NeXLMatrixCorrection.AMethod
A(unk::MultiZAF, std::MultiZAF, θunk::AbstractFloat, θstd::AbstractFloat)

The A (absorption) correction for unk relative to std.

source
NeXLMatrixCorrection.AMethod
A(unk::ZAFCorrection, std::ZAFCorrection, cxr::CharXRay, θunk::AbstractFloat, θstd::AbstractFloat)

Computes the absorption correction.

source
NeXLMatrixCorrection.FMethod
F(unk::MultiZAF, std::MultiZAF, θunk::AbstractFloat, θstd::AbstractFloat)

The F (fluoresence) correction for unk relative to std.

source
NeXLMatrixCorrection.FMethod

F(reed::ReedFluorescence, secondary::CharXRay, toa::Float64)

Compute the enhancement of the secondary characteristic X-ray due to the primaries specified in reed.

source
NeXLMatrixCorrection.FMethod
F(unk::ZAFCorrection, std::ZAFCorrection, cxr::CharXRay, θunk::AbstractFloat, θstd::AbstractFloat)

Computes the secondary fluorescence correction.

source
NeXLMatrixCorrection.FχMethod
Fχ(mc::MatrixCorrection, xray::CharXRay, θtoa::Real)

Integral of the area under the absorption corrected ϕ(ρz)-curve from ρz = 0 to ∞.

source
NeXLMatrixCorrection.FχpMethod
Fχp(mc::NeXLMatrixCorrection, xray::CharXRay, θtoa::Real, τ::Real)

The partial integral of the absorption corrected ϕ(ρz) curve from ρz = 0 to τ #

source
NeXLMatrixCorrection.RMethod
R(mat::Material, u0)

Backscatter factor as a function of material and overvoltage. Compares well to PAP Figure 23.

source
NeXLMatrixCorrection.ZAMethod
ZA(
  unk::MatrixCorrection,
  std::MatrixCorrection,
  xray::CharXRay,
  θunk::AbstractFloat,
  θstd::AbstractFloat
)

The atomic number (Z) and absorption (A) correction factors.

source
NeXLMatrixCorrection.ZAMethod
ZA(unk::ZAFCorrection, std::ZAFCorrection, cxr::CharXRay, θunk::AbstractFloat, θstd::AbstractFloat)

Computes the combined atomic number and fluorescence correction.

source
NeXLMatrixCorrection.ZAFcMethod
ZAFc(unk::ZAFCorrection, std::ZAFCorrection, cxr::CharXRay, θunk::AbstractFloat, θstd::AbstractFloat)

Computes the combined correction for atomic number, absorption, secondary fluorescence and generation.

source
NeXLMatrixCorrection.coatingMethod
coating(unk::MultiZAF, std::MultiZAF, θunk::AbstractFloat, θstd::AbstractFloat)

The conductive (or other) coating correction factor.

source
NeXLMatrixCorrection.coatingMethod
coating(unk::ZAFCorrection, std::ZAFCorrection, cxr::CharXRay, θunk::AbstractFloat, θstd::AbstractFloat)

Computes the coating correction.

source
NeXLMatrixCorrection.computeMethod
compute(::Type{UnmeasuredElementRule}, inp::Dict{Element,Float64})::Dict{Element,Float64}

A null UnmeasuredElementRule. Just returns the inputs.

source
NeXLMatrixCorrection.computeZAFsMethod
computeZAFs(
    iter::Iteration,
    est::Material,
    stdZafs::Dict{KRatio,MultiZAF}
)::Dict{Element, Float64}

Given an estimate of the composition compute the corresponding k-ratios.

source
NeXLMatrixCorrection.correctkratiosMethod
correctkratios(krs::AbstractVector{KRatio}, coating::Material, θtoa::Real, ρz::Real)::Vector{KRatio}

This function is mainly for pedagogical purposes. It takes a KRatio[], a coating Material on the unknown, and a mass-thickness (g/cm²) and creates a new KRatio[] that accounts for the intensity missing due to absorption by the coating. Favor the coating correction built into ZAFCorrection or MultiZAF.

source
NeXLMatrixCorrection.dEdρsMethod
dEdρs(args::Dict{Label,AbstractFloat}, mat::MaterialLabel, Ekev::AbstractFloat, elms)

The function P&P uses to describe the deceleration of an electron in a material. Output units are (keV/cm)/(g/cm^3) = keV cm^2/g. (PAP eqn 5)

source
NeXLMatrixCorrection.detailMethod
detail(::Type{DataFrame}, mzs::AbstractArray{Tuple{MultiZAF, MultiZAF}})::DataFrame

Tabulate the details of a matrix correction relative to the specified unknown and standard in a DataFrame.

source
NeXLMatrixCorrection.detailMethod
detail(::Type{DataFrame}, unk::MultiZAF, std::MultiZAF)::DataFrame

Tabulate each term in the MultiZAF matrix correction in a DataFrame.

source
NeXLMatrixCorrection.familyfactorMethod
familyfactor(shellA::AtomicSubShell, shellB::AtomicSubShell)::Float64

Accounts for the differences in ionization cross section between K , L and M shells

source
NeXLMatrixCorrection.fluorescencecorrectionMethod
fluorescencecorrection(::Type{ReedFluorescence}, comp::Material, primary::Vector{CharXRay}, secondary::AtomicSubShell, e0::Float64)

Construct an instance of a ReedFluorescence correction structure to compute the secondary fluorescence due to a primary characteristic X-ray in the specified material and beam energy.

source
NeXLMatrixCorrection.fluorescencecorrectionMethod
fluorescence(fltype::Type{ReedFluorescence}, comp::Material, secondary::AtomicSubShell, e0::Float64)

Construct an instance of a fltype correction structure to compute the secondary fluorescence in the specified material and beam energy.

source
NeXLMatrixCorrection.gZAFcMethod
gZAFc(unk::MultiZAF, std::MultiZAF, θunk::AbstractFloat, θstd::AbstractFloat)
gZAFc(kr::KRatio, unkComp::Material; mc::Type{<:MatrixCorrection} = XPP, fc::Type{<:FluorescenceCorrection} = ReedFluorescence, cc::Type{<:CoatingCorrection} = Coating)

The combined generation, atomic number, absorption and generation corrections.

source
NeXLMatrixCorrection.generationMethod
generation(unk::MultiZAF, std::MultiZAF)

The generation correction for unk relative to std. Usually, 1.0 unless the standard and unknown were collected at different beam energies.

source
NeXLMatrixCorrection.generationMethod
generation(unk::ZAFCorrection, std::ZAFCorrection, ass::AtomicSubShell)

Computes a correction factor for differences X-ray generation due to differences in beam energy.

source
NeXLMatrixCorrection.kMethod
k(unk::MultiZAF, std::MultiZAF, θunk::AbstractFloat, θstd::AbstractFloat)

The computed k-ratio for the unknown relative to standard.

source
NeXLMatrixCorrection.kMethod
k(unk::MultiZAF, std::MultiZAF, θunk::AbstractFloat, θstd::AbstractFloat)

The computed k-ratio for the unknown relative to standard.

source
NeXLMatrixCorrection.kcoatingMethod
kcoating(ty::Type{<:MatrixCorrection}, subtrate::Material, coating::Material, cxr::CharXRay, e0::Real, toa::Real, τ::Real)

Estimate the k-ratio for a coating of mass-thickness τ (g/cm²) on the specified substrate. The standard for the coating is assumed to be of the same material as the coating.

source
NeXLMatrixCorrection.lenardcoefficientMethod

lenardcoefficient(e0::Float64, ashell::AtomicSubShell)

Computes the Lenard coefficient according to the algorithm of Heinrich. "Heinrich K. F. J. (1967) EPASA 2, paper no. 7"

source
NeXLMatrixCorrection.massthicknessMethod

" massthickness(ty::Type{<:MatrixCorrection}, subtrate::Material, coating::Material, cxr::CharXRay, e0::Real, toa::Real, k::Real)

Estimate the mass-thickness of a ultra-thin layer of a coating material on a substrate from a measured k-ratio k of a characteristic X-ray cxr. Works for k-ratios of the order of 1 %. The standard for the coating is assumed to be of the same material as the coating.

source
NeXLMatrixCorrection.quantifyMethod
quantify(iter::Iteration, label::Label, measured::Vector{KRatio}, maxIter::Int = 100)::IterationResult
quantify(iter::Iteration, name::String, measured::Vector{KRatio})::IterationResult
quantify(ffr::FitResult; strip::AbstractVector{Element} = Element[], mc::Type{<:MatrixCorrection} = XPP, fl::Type{<:FluorescenceCorrection} = ReedFluorescence, cc::Type{<:CoatingCorrection} = Coating)::IterationResult

Perform the iteration procedurer as described in iter using the measured k-ratios to produce the best estimate Material in an IterationResult object. The third form makes it easier to quantify the k-ratios from filter fit spectra.

source
NeXLMatrixCorrection.quantifyMethod

quantify(ffr::FitResult, strip::AbstractVector{Element} = [], mc::Type{<:MatrixCorrection} = XPP, fl::Type{<:FluorescenceCorrection} = ReedFluorescence)::IterationResult

Facilitates quantifying FilterFitResult or BasicFitResult objects from extracting k-ratios from measured spectra.

source
NeXLMatrixCorrection.steps1Function
steps1(sample, elms, shell, all)

steps1 requires as data MassFractionLabel, AtomicWeightLabel, JzLabel, E0Label, mLabel in an UncertainValues

source
NeXLMatrixCorrection.updateMethod
update(
    ::NaiveUpdateRule,
    prevcomp::Material,
    measured::Vector{KRatio},
    zafs::Dict{Element, Float64}
)::Dict{Element,Float64}

Determine the next estimate of the composition that brings the estimate k-ratios closer to measured.

source
NeXLMatrixCorrection.zafcorrectionFunction
zafcorrection(
  mctype::Type{<:MatrixCorrection},
  fctype::Type{<:FluorescenceCorrection},
  cctype::Type{<:CoatingCorrection},
  mat::Material,
  cxrs,
  e0,
  coating=missing
)

Constructs a MultiZAF around the mctype and fctype algorithms for a collection of CharXRay cxrs.

source
NeXLMatrixCorrection.zafcorrectionMethod
zafcorrection(
  mctype::Type{<:MatrixCorrection},
  fctype::Type{<:FluorescenceCorrection},
  cctype::Type{<:CoatingCorrection},
  mat::Material,
  ashell::AtomicSubShell,
  e0::Real,
  coating::Union{Film,AbstractVector{Film},Missing}
)

Constructs an ZAFCorrection object using the mctype correction model with the fluorescence model for the specified parameters.

source
NeXLMatrixCorrection.zafcorrectionMethod
zafcorrection(
  mctype::Type{<:MatrixCorrection},
  fctype::Type{<:FluorescenceCorrection},
  cctype::Type{<:CoatingCorrection},
  unk::Material,
  std::Material,
  cxrs,
  e0;
  unkCoating::Union{Film,AbstractVector{Film},Missing} = missing,
  stdCoating::Union{Film,AbstractVector{Film},Missing} = missing,
)

Constructs a tuple of MultiZAF around the mctype and fctype correction algorithms for the unknown and standard for a collection of CharXRay cxrs.

source
NeXLMatrixCorrection.zafcorrectionMethod
zafcorrection(
   mctype::Type{<:MatrixCorrection},
   fctype::Type{<:FluorescenceCorrection},
   cctype::Type{<:CoatingCorrection},
   unk::Material,
   std::Material,
   ashell::AtomicSubShell,
   e0::Real;
   unkCoating::Union{Film,AbstractVector{Film},Missing} = missing,
   stdCoating::Union{Film,AbstractVector{Film},Missing} = missing,
)

Creates a matched pair of ZAFCorrection objects using the matrix correction algorithm for the specified unknown and standard.

source
NeXLMatrixCorrection.ϕabsMethod
ϕabs(mc::MatrixCorection, ρz, xray::CharXRay, θtoa::AbstractFloat)

Computes the absorbed ϕ(ρz) curve according to the XPP algorithm.

source
NeXLUncertainties.asaMethod
asa(::Type{DataFrame}, mzs::AbstractArray{Tuple{MultiZAF,MultiZAF}}, θunk::AbstractFloat, θstd::AbstractFloat)::DataFrame

Tabulate a matrix correction relative to a specified Dict of unknowns and standards in a DataFrame.

source
NeXLUncertainties.asaMethod
asa(::Type{DataFrame}, unk::MultiZAF, std::MultiZAF, θunk::AbstractFloat, θstd::AbstractFloat)::DataFrame

Tabulate a matrix correction relative to the specified unknown and standard in a DataFrame.

source
NeXLUncertainties.asaMethod
NeXLUncertainties.asa(::Type{DataFrame}, rur::RecordingUpdateRule)::DataFrame

Tabulate the iteration steps in a DataFrame.

source
NeXLUncertainties.asaMethod
NeXLUncertainties.asa(::Type{DataFrame}, unk::ZAFCorrection, std::ZAFCorrection, trans::AbstractVector{Transition},
θunk::AbstractFloat, θstd::AbstractFloat)::DataFrame

Tabulate a matrix correction relative to the specified unknown and standard for the iterable of Transition, trans.

source