Spectrum
Working with Spectrum objects
The Spectrum type represents a single spectrum with associated properties. The channel data is indexed like a Vector. The property data is indexed using Symbol objects. A set of spectrum properties are defined in the library and the user can create additional ones.
In addition to the Vector-like properties, the Spectrum associates energy bins with each channel. The energy bins are continuous, non-overlapping and monotonic. Functions like energy(...) maps channel index to energies and channel(...) to make energies to channel index.
NeXLSpectrum.Spectrum — TypeSpectrum{T<:Real} <: AbstractVector{T}
Spectrum(energy::EnergyScale, data::Vector{<:Real}, props::Dict{Symbol,Any})Construct a structure to hold spectrum data (energy scale, counts and metadata).
See NeXLSpectrum.EnergyScale or NeXLSpectrum.LinearEnergyScale
Example:
julia> spec = Spectrum(LinearEnergyScale(0.0,10.0),
collect(1:1024),
Dict{Symbol,Any}(:BeamEnergy=>10.0e3, :LiveTime=>30.0))
** 1024.0
*********
*****************
************************
*******************************
**************************************
*********************************************
****************************************************
***********************************************************
******************************************************************
*************************************************************************
******************************************************************************** 10.0 keV
Spectrum[3062][10.0 keV, Unknown, 525000.0 counts]Spectra are usually loaded from disk using:
s1 = loadspectrum(joinpath(path,"ADM-6005a_1.msa")) # Auto-detects the file type (supports ISO/EMSA, Bruker, ASPEX, ??? files)
* 128860.0
*
*
*
*
*
* *
* *
* *
* ** **
* ** ***** ** **
**************************************************************************************************** 20.0 keV
Spectrum{Float64}[ADM-6005a_1, -484.20818 + 5.01716⋅ch eV, 4096 ch, 20.0 keV, Unknown, 6.81e6 counts]Spectrum implements indexing using various different mechanisms. If spec is a Spectrum then
spec[123] # will return the number of counts in channel 123
spec[123:222] # will return a Vector of counts from channel 123:222
spec[134.] # will return the number of counts in the channel at energy 134.0 eV
spec[134.0:270.0] # will return a Vector of counts for channels with energies between 134.0 eV and 270.0 eV
spec[n"Ca L3-M5"] # Counts in the channel at the Ca L3-M5 characteristic X-ray energy
spec[:Comment] # will return the meta-data item named :CommentBasic spectrum math is supported using the operators *, /, ÷, +, and -. If s1, s2 and s3 are Spectrum objects then:
2*s1 # A Spectrum containing twice as many counts in each channel
2.0*s1 # A Spectrum containing twice as many counts in each channel
s1/2.0 # A spectrum containing half as many counts in each channel
s1÷2 # A spectrum containing half as many counts in each channel (Only works for Spectrum{<:Integer})
s1 + s2 # A Spectrum containing channel-by-channel sum of counts in s1 and s1 and common properties.
s1 - s2 # The channel-by-channel differenceSpectrum metadata is identified by a Symbol. These Symbols are used within NeXLSpectrum. You can create others to associate other data items with a Spectrum.
:BeamEnergy # In eV
:Elevation # In radians
:TakeOffAngle # In radians (Detector position)
:Azimuthal # In radians (Detector position)
:WorkingDistance # In cm (not mm!!!!)
:LiveTime # In seconds
:RealTime # In seconds
:DeadFraction # Fractional
:ProbeCurrent # In nano-amps
:Name # A `String`
:Owner # A `String`
:Sample # Description of the sample
:StagePosition # A Dict{Symbol,Real} with entries :X, :Y, :Z, :R, :T, B: in cm and degrees
:Comment # A `String`
:Composition # A `Material` (known composition, not measured)
:Elements # A collection of elements in the material like `[ n"Fe", n"Ca", n"Si" ]`
:Detector # A Detector like a BasicEDS or another EDSDetector
:Filename # Source filename
:Coating # A `Film` or `Film[]` (eg. 10 nm of C|Au etc.)
:AcquisitionTime # Date and time of acquisition (`DateTime` struct)
:Signature # Dict{Element,Real} with the "particle signature"
:SolidAngle # Detector solid angle is steradians (area/dist²)Spectrum Image items:
:ImageRotation # Rotation of the primary scan direction from the :X axis towards the :Y axis in radiansLess common items:
:ImageMag # Magnification (assuming a 3.5" image) of the first image
:ImageZoom # Additional zoom for second image in a two image TIFF
:Operator # Analyst in ASPEX TIFF files
:Image1, :Image2 ... # Images associated with the spectrum
:BrukerThroughtput # Nominal throughtput setting on a Bruker detector
:DetectorSerialNumber # EDS detector serial number
:DetectorModel # Vendor model name
:DetectorThickness # Thickness of detector active area
:DeadLayerThickness # Thickness of Si dead layer on the entrance surface of the detector
:Window # Window construction details
:DetectorSolidAngle # Collection solid angle of the X-ray detector
:ChamberPressure # Vacuum presure in the sample chamber
:ChamberAtmosphere # Nominally the composition of the residual gas in the chamberXRF related items:
:BeamEnergy # Accelarating voltage within X-ray tube (eV)
:XRFTubeAnode # Element from which the X-ray tube is constructed
:ProbeCurrent # Electron current in the X-ray tube
:XRFTubeIncidentAngle # Incident angle of electron beam in tube
:XRFTubeTakeOffAngle # Take-off angle from tube
:XRFExcitationAngle # Angle of incidence of the X-ray beam on the sample
:XRFDetectionAngle # Angle of the detector relative to the sample
:XRFExcitationPathLength # Distance from X-ray source to sample
:XRFDetectionPathLength # Distance from the sample to the X-ray detector
:XRFSampleTilt # Additional tilt of the sample
:XRFTubeWindow # Construction of the tube windowNot all spectra will define all properties. Algorithms can define the NeXLCore.minproperties(ty::Type) method to specify which properties are required by an algorithm of ty::Type. Then hasminrequired and requiredbutmissing methods will determine whether a Spectrum or Dict{Symbol,Any} is suitable for an algorithm.
Acccessing Properties
Most properties are accessed using spec[:Symbol] notation. Some combined or special properties have special methods.
NeXLSpectrum.dose — Functiondose(spec::Spectrum, def=missing)
dose(props::Dict{Symbol, Any}, def=missing)The probe dose in nA⋅s and is equals spec[:LiveTime]⋅spec[:ProbeCurrent].
dose(hss::HyperSpectrum) # Average dose per pixel
dose(hss::HyperSpectrum, idx...) # Dose for the `idx` pixelReturns the product of the probe current and the live-time on a per pixel basis.
NeXLCore.elms — Methodelms(spec::Spectrum, withcoating = false)::Set{Element}Returns a set of the elements associated with this spectrum. withcoating determines whether the coating elements are also added.
Displaying Spectrum Data
NeXLUncertainties.asa — FunctionNeXLUncertainties.asa(::Type{DataFrame}, spec::Spectrum; properties::Bool = false)Converts the spectrum energy and counts data into a DataFrame.
NeXLUncertainties.asa(::Type{DataFrame}, spec::AbstractArray{Spectrum})::DataFrameReturns a DataFrame that summarizes the list of spectra.
NeXLUncertainties.asa(::Type{DataFrame}, spec::AbstractDict{Element, Spectrum})::DataFrameReturns a DataFrame that summarizes a dictionary of standard spectra.
NeXLUncertainties.asa(::Type{DataFrame}, ffrs::AbstractVector{<:FitResult}; charOnly = true, withUnc = false, format = :normal # :pivot or :long)Return generic FitResult as a DataFrame.
Format:
- :normal - One row per spectrum, one column per k-ratio
- :pivot - One row per ROI, one column per spectrum (optional: column for 1σ uncertainty on k-ratio)
- :long - One row per spectrum, feature, measured k-ratio
NeXLUncertainties.asa(::Type{DataFrame}, fr::FitResult; withUnc = false)Summarize the ROI and k-ratio data within a FitResult structure as a DataFrame.
NeXLUncertainties.asa(
::Type{DataFrame},
ffr::FilterFitResult;
charOnly::Bool=true,
material=nothing,
mc = XPP, fc=ReedFluorescence,
columns = ( :counts, ) # Selected from ( :roi, :peakback, :counts, :dose )
)::DataFrameTabulate details about each region-of-interest in the 'FilterFitResult' in a 'DataFrame'.
- If charOnly then only display characteristic X-ray data (not escapes etc.)
- If
materialis a Material then the computed k-ratio (KCalc) will also be tabulated along with kmeas/kcalc (KoKCalc). - columns - Select optional column outputs (see below)
Columns:
- Spectrum : UnknownLabel - Identifies the fit spectrum
- Feature : Label - Identifies the fit feature (Vector{CharXRay} or other)
- Reference: String - Name of reference against which :Spectrum was fit over :Feature
- K : Float64 - The multiplicative fit constant
- dK : Float64 - The 1σ uncertainty in :K
- :Start : Int - Start index for fit channels (:roi ∈ columns)
- :Stop : Int - Stop index for fit channels (:roi ∈ columns)
- :Counts : Float64 - Total counts in characteristic peak (peak-back) (:peakback || :counts ∈ columns)
- :Back : Float64 - Total counts in background under the characteristic peak (:peakback ∈ columns)
- :PtoB : Float64 - Peak-to-Background assuming 10 eV/channel (:peakback ∈ columns)
- :KCalc : Float64 - Computed k-ratio assuming a composition. (Requires
materialargument to be specified.) - :KoKcalc : Float64 - Ratio of measured/computed k-ratio. (Requires
materialargument to be specified.) - :LiveTime : Float64 - Acquisiton live time (s) (:dose ∈ columns)
- :ProbeCurrent : Float64 - Probe current (nA) (:dose ∈ columns)
- :DeadPct : Float64 - Dead time in ProbeCurrent (:dose ∈ columns)
- :RefCountsPernAs : Float64 - Estimated counts in :Reference in :Feature per unit dose. (:counts ∈ columns)
- :CountsPernAs : Float64 - Estimated counts in :Spectrum in :Feature per unit dose. (:counts ∈ columns)
NeXLUncertainties.asa(::Type{DataFrame}, ffp::FilterFitPacket)Summarize the FilteredReference structs within a FilterFitPacket as a DataFrame.
Base.keys — FunctionBase.keys(spec::Spectrum)Return the defined properties as a set of Symbols.
Base.haskey — FunctionBase.haskey(spec::Spectrum, sym::Symbol)Is the specified key defined?
Gadfly.plot — Functionplot(
specs::Union{Spectrum...,AbstractVector{Spectrum{<:Real}}};
klms=[],
edges=[],
escapes=[],
coincidences=[],
autoklms = false,
xmin=0.0,
xmax=missing,
norm=:None,
yscale=1.05,
ytransform = identity,
style=NeXLSpectrumStyle,
palette=NeXLPalette
)::PlotRequired:
specs::AbstractVector{Spectrum};Named:
klms = Union{Element, CharXRay}[ ]
edges = Union{Element, AtomicSubShell}[ ]
escapes = Union{Element, CharXRay}[ ],
coincidences = CharXRay[ ],
autoklms = false # Add KLMs based on elements in spectra
xmin = 0.0 # Min energy (eV)
xmax = missing # Max energy (eV) (defaults to max(:BeamEnergy))
norm = NoScaling() | ScaleDoseWidth() | ScaleDose() | ScaleSum() | ScaleROISum() | ScalePeak() | (<: SpectrumScaling)()
yscale = 1.05 # Fraction of max intensity for ymax over [max(lld,xmin):xmax]
ytransform = identity | log10 | sqrt | ??? # How to transform the counts data before plotting
style=NeXLSpectrumStyle (or another Gadfly.style)
palette = NeXLPalette | Colorant[ ... ] # Colors for spectra...
customlayers = Gadfly.Layer[] # Allows additional plot layers to be addedPlot a multiple spectra on a single plot using Gadfly.
Gadfly.plot(
ffr::FilterFitResult,
roi::Union{Nothing,AbstractUnitRange{<:Integer}} = nothing;
palette = NeXLPalette,
style = NeXLSpectrumStyle,
xmax::Union{AbstractFloat, Nothing} = nothing,
comp::Union{Material, Nothing} = nothing,
det::Union{EDSDetector, Nothing} = nothing,
resp::Union{AbstractArray{<:AbstractFloat,2},Nothing} = nothing,
yscale = 1.0
)Plot the sample spectrum, the residual and fit regions-of-interests and the associated k-ratios.
Gadfly.plot(fr::FilteredReference; palette = NeXLPalette))Plot a filtered reference spectrum.
plot(ff::TopHatFilter, fr::FilteredReference)Plot a color map showing the filter data relevant to filtering the specified FilteredReference.
Gadfly.plot(vq::VectorQuant, chs::UnitRange)Plots the "vectors" used to quantify various elements/regions-of-interest over the range of channels specified.
Gadfly.plot(deteff::Union{DetectorEfficiency,AbstractVector{DetectorEfficiency}}, emax = 20.0e3)Plots the detector efficiency function assuming the detector is perpendicular to the incident X-rays.
Gadfly.plot(ffp::FilterFitPacket; kwargs...)Plots the reference spectra which were used to construct a FilterFitPacket.
plot(wind::Union{AbstractWindow, AbstractArray{<:AbstractWindow}}; xmax=20.0e3, angle=π/2, style=NeXLSpectrumStyle)Plot the window transmission function.
Gadfly.plot(
dfr::DirectFitResult,
roi::Union{Nothing,AbstractUnitRange{<:Integer}} = nothing;
palette = NeXLPalette,
style = NeXLSpectrumStyle,
xmax::Union{AbstractFloat, Nothing} = nothing,
comp::Union{Material, Nothing} = nothing,
det::Union{EDSDetector, Nothing} = nothing,
resp::Union{AbstractArray{<:AbstractFloat,2},Nothing} = nothing,
yscale = 1.0
)Plot the sample spectrum, the residual and fit regions-of-interests and the associated k-ratios.
Energy Scale Functions
These functions handle mapping channel index to energies and vice versa.
NeXLCore.energy — FunctionNeXLCore.energy(ch::Integer, sc::EnergyScale)
NeXLCore.energy(ch::Int, spec::Spectrum)
NeXLCore.energy(ch::Int, det::EDSDetector)Returns the energy (in eV) for the low energy side of the bin representing the ch-th channel.
Example:
les = LinearEnergyScale(3.0, 10.1)
energy(101,lsc) == 10.1*101 + 3.0
energy(101,lsc) - energy(100,lsc) == 10.1
pes = PolyEnergyScale([ 3.0, 10.1, 0.001])
energy(101,pes) == 3.0 + 10.0*101 + 0.001*101^2NeXLCore.energy(ch::Int, spec::Spectrum)::Float64What energy is associated with ch in spec?
NeXLSpectrum.channel — Functionchannel(::Type{Float64}, eV::AbstractFloat, sc::EnergyScale)::Float64Returns the fractional-channel index of the specified energy X-ray (in eV).
channel(eV::AbstractFloat, sc::EnergyScale)::Int
channel(eV::Float64, spec::Spectrum)::Int
channel(eV::Float64, det::EDSDetector)::IntReturns the integer index of the channel for the specified energy X-ray (in eV).
NeXLSpectrum.rangeofenergies — Functionrangeofenergies(ch::Integer, spec::Spectrum)Returns the low and high energy extremes for the channels ch.
NeXLSpectrum.channelwidth — Functionchannelwidth(ch::Int, spec::Spectrum)::Float64
channelwidth(ch::Int, spec::HyperSpectrum)::Float64Returns the width of the ch channel in eV.
channelwidth(spec::Spectrum)::Float64Returns the mean channel width.
NeXLSpectrum.energyscale — Functionenergyscale(es::EnergyScale, channels)
energyscale(det::EDSDetector)Computes the energy associated with a range of channel indexes and returns it as an Array.
energyscale(spec::Spectrum)
energyscale(spec::Spectrum, chs::AbstractRange{<:Integer})Returns an array with the bin-by-bin energies
Extracting the Counts Data
The channel data in a Spectrum object spec can be get/set using standard Array indexing techniques.
> spec[123]
> spec[123:345]
> spec[1234.0] # Channel containing energy 1234.0 eV
> spec[123] = 99Alternatively, the counts(...) method can be used.
NeXLSpectrum.counts — Functioncounts(spec::Spectrum, ::Type{T}, applyLLD=false)::Vector{T} where {T<:Number}
counts(spec::Spectrum, channels::AbstractUnitRange{<:Integer}, ::Type{T}, applyLLD=false)::Vector{T} where {T<:Number}Creates a copy of the spectrum counts data as the specified Number type. If the spectrum has a :Detector property then the detector's lld (low-level discriminator) and applyLLD=true then the lld is applied to the result by setting all channels less-than-or-equal to det.lld to zero. This method throws an error if the counts data cannot be converted without loss to the type T.
counts(hss::HyperSpectrum{T, N, NP})::Array{T,NP}Creates type-friendly view of the counts data array. Use of this function helps to avoid performance penalties associated with boxing/unboxing the counts data.
counts(hss::HyperSpectrum{T, N, NP}, ci::CartesianIndex)::Vector{T}Access the counts data associated with the pixel ci.
counts(hss::HyperSpectrum{T, N, NP}, ci::CartesianIndex, ch::Int)::TAccess the counts data at the pixel represented by ci and the channel represented by ch.
counts(...) can apply the :LLD property to zero counts below a specified channel.
NeXLSpectrum.lld — Functionlld(det::EDSDetector)Low level detection limit in channels. Channels at or below this value will be zeroed when the lld is applied.
lld(spec::Spectrum)Gets the low-level discriminator associated with this spectrum if there is one.
Defining Detectors
Build a detector to match the data in a Spectrum.
NeXLSpectrum.simpleEDS — FunctionsimpleEDS(chCount::Integer, width::Float64, offset::Float64, fwhmatmnka::Float64, lld::Int = channel(150.0 eV))Construct simple model of an EDS detector.
simpleEDS(spec::Spectrum, fwhmatmnka::Float64)Build a EDSDetector object for this spectrum with the specified FWHM at Mn Kα.
NeXLSpectrum.matching — Functionmatching(spec::Spectrum, resMnKa::Float64, lld::Int=1, minByFam::Dict{Shell,Element} = Dict())::BasicEDSBuild an EDSDetector to match the channel count and energy scale in this spectrum.
matching(spec::HyperSpectrum, resMnKa::Float64, lld::Int=1)::BasicEDSBuild an EDSDetector to match the channel count and energy scale in this spectrum.
Statistical Sub-Sampling
NeXLSpectrum.subdivide — Functionsubdivide(spec::Spectrum, n::Int)::Vector{Spectrum}Splits the event data in one spectrum into n spectra by assigning each event to a pseudo-random choice of one of the n result spectra. Produces n spectra that act as though the original spectrum was collected in n time intervals of LiveTime/n. This is quite slow because it needs to call rand() for each count in the spectrum (not just each channel).
NeXLSpectrum.subsample — Functionsubsample(spec::Spectrum, frac::Float64)Subsample the counts data in a spectrum according to a statistically valid algorithm. Returns spec if frac>=1.0.
Processing Channel Data
NeXLSpectrum.kratio — Functionkratio(unk::Spectrum, std::Spectrum, back1::AbstractUnitRange{<:Integer}, peak::AbstractUnitRange{<:Integer}, back2::AbstractUnitRange{<:Integer})::UncertainValue
kratio(unk::Spectrum, std::Spectrum, back1::StepRangeLen{Float64}, peak::StepRangeLen{Float64}, back2::StepRangeLen{Float64})::UncertainValueThe k-ratio of unk relative to std corrected for dose. Requires that unk and std have the properties :LiveTime and :ProbeCurrent defined.
kratio(ffr::FitResult, cxr::CharXRay)::KRatioExtract the k-ratio associated with the specified CharXRay (zero if not measured)
NeXLSpectrum.integrate — Functionintegrate(spec::Spectrum, channels::AbstractUnitRange{<:Integer})
integrate(spec::Spectrum, energyRange::StepRangeLen{Float64})Sums all the counts in the specified channels. No background correction. The energyRange version includes the fractional contributions from partial channels when the energies don't perfectly map to channels.
integrate(spec::Spectrum, back1::AbstractUnitRange{<:Integer}, peak::AbstractUnitRange{<:Integer}, back2::AbstractUnitRange{<:Integer})::UncertainValue
integrate(spec::Spectrum, back1::StepRangeLen{Float64}, peak::StepRangeLen{Float64}, back2::StepRangeLen{Float64})::UncertainValueSums all the counts in the specified channels with background correction using the background intervals.
integrate(spec::Spectrum)Total integral of all counts from the LLD to the beam energy
NeXLSpectrum.estimatebackground — Functionestimatebackground(data::AbstractArray{Float64}, channel::Int, width::Int=5, order::Int=2)Returns the tangent to the a quadratic fit to the counts data centered at channel with width
NeXLSpectrum.modelBackground — FunctionmodelBackground(spec::Spectrum, chs::AbstractUnitRange{<:Integer}, ash::AtomicSubShell)spec: A spectrum containing a peak centered on chs chs: A range of channels containing a peak ash: The edge (as an AtomicSubShell)
A simple model for modeling the background under a characteristic x-ray peak. The model fits a line to low and high energy background regions around firsr(chs) and last(chs). If the low energy line extended out to the edge energy is larger than the high energy line at the same energy, then a negative going edge is fit between the two. Otherwise a line is fit between the low energy side and the high energy side. This model only works when there are no peak interference over the range chs.
modelBackground(spec::Spectrum, chs::AbstractUnitRange{<:Integer})
modelBackground(spec::Spectrum, chs::AbstractUnitRange{<:Integer}, ash::AtomicSubShell)spec: A spectrum containing a peak centered on chs chs: A range of channels containing a peak ash: The largest edge within the range of channels chs associated with the characteristic peak
A simple model for modeling the background under a characteristic x-ray peak. The model fits a line between the low and high energy background regions around first(chs) and last(chs). This model only works when there are no peak interference over the range chs.
NeXLSpectrum.extractcharacteristic — Functionextractcharacteristic(spec::Spectrum, lowBack::AbstractUnitRange{<:Integer}, highBack::AbstractUnitRange{<:Integer})::Vector{Float64}Extract the characteristic intensity for the peak located within chs with an edge at ash.
NeXLSpectrum.peak — Functionpeak(spec::Spectrum, chs::AbstractUnitRange{<:Integer}, ash::AtomicSubShell)::Float64Estimates the peak intensity for the characteristic X-ray in the specified range of channels.
NeXLSpectrum.background — Functionbackground(spec::Spectrum, chs::AbstractUnitRange{<:Integer}, ash::AtomicSubShell)::Float64Estimates the background intensity for the characteristic X-ray in the specified range of channels.
NeXLSpectrum.peaktobackground — Functionpeaktobackground(spec::Spectrum, chs::AbstractUnitRange{<:Integer}, ash::AtomicSubShell)::Float64Estimates the peak-to-background ratio for the characteristic X-ray intensity in the specified range of channels which encompass the specified AtomicSubShell.
peaktobackground(ffr::FilterFitResult, backwidth::Float64=10.0)::Float64The peak-to-background ratio as determined from the raw and residual spectra integrated over the fit region-of-interest and scaled to backwidth eV of continuum (nominally 10 eV).
NeXLSpectrum.estkratio — Functionestkratio(unk::Spectrum, std::Spectrum, chs::AbstractUnitRange{<:Integer})Estimates the k-ratio from niave models of peak and background intensity. Only works if the peak is not interfered.
Miscellaneous Functions
NeXLSpectrum.details — Functiondetails(io, spec::Spectrum)
details(spec::Spectrum)::StringOutputs a description of the data in the spectrum.
NeXLSpectrum.commonproperties — Functioncommonproperties(specs::AbstractArray{Spectrum})
commonproperties(props1::Dict{Symbol,Any}, props2::Dict{Symbol,Any})Return the properties that are held in common by all the spectra.
NeXLSpectrum.maxspectrum — Functionmaxspectrum(specs::AbstractArray{Spectrum{T}})::Spectrum{T} where T<:RealCompute the max-pixel spectrum for the specified spectra.
Base.findmax — FunctionBase.findmax(spec::Spectrum, chs::AbstractRange{<:Integer})
Base.findmax(spec::Spectrum)Returns the (maximum intensity, channel index) over the specified range of channels